libMesh
subdomains_ex1.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1>Subdomains Example 1 - Solving on a Subdomain</h1>
21 // \author Tim Kroger
22 // \date 2010
23 //
24 // This example builds on the example 4 by showing what to do in
25 // order to solve an equation only on a subdomain.
26 
27 
28 // C++ include files that we need
29 #include <iostream>
30 #include <algorithm>
31 #include <math.h>
32 
33 // Basic include file needed for the mesh functionality.
34 #include "libmesh/libmesh.h"
35 #include "libmesh/mesh.h"
36 #include "libmesh/mesh_generation.h"
37 #include "libmesh/exodusII_io.h"
38 #include "libmesh/gmv_io.h"
39 #include "libmesh/gnuplot_io.h"
40 #include "libmesh/linear_implicit_system.h"
41 #include "libmesh/equation_systems.h"
42 
43 // Define the Finite Element object.
44 #include "libmesh/fe.h"
45 
46 // Define Gauss quadrature rules.
47 #include "libmesh/quadrature_gauss.h"
48 
49 // Define the DofMap, which handles degree of freedom
50 // indexing.
51 #include "libmesh/dof_map.h"
52 
53 // Define useful datatypes for finite element
54 // matrix and vector components.
55 #include "libmesh/sparse_matrix.h"
56 #include "libmesh/numeric_vector.h"
57 #include "libmesh/dense_matrix.h"
58 #include "libmesh/dense_vector.h"
59 
60 // Define the PerfLog, a performance logging utility.
61 // It is useful for timing events in a code and giving
62 // you an idea where bottlenecks lie.
63 #include "libmesh/perf_log.h"
64 
65 // The definition of a geometric element
66 #include "libmesh/elem.h"
67 
68 #include "libmesh/mesh_refinement.h"
69 
70 // Classes needed for subdomain computation.
71 #include "libmesh/system_subset_by_subdomain.h"
72 
73 #include "libmesh/string_to_enum.h"
74 #include "libmesh/getpot.h"
75 
76 // Bring in everything from the libMesh namespace
77 using namespace libMesh;
78 
79 
80 
81 // Function prototype. This is the function that will assemble
82 // the linear system for our Poisson problem. Note that the
83 // function will take the EquationSystems object and the
84 // name of the system we are assembling as input. From the
85 // EquationSystems object we have access to the Mesh and
86 // other objects we might need.
88  const std::string & system_name);
89 
90 // Exact solution function prototype.
91 Real exact_solution (const Real x,
92  const Real y = 0.,
93  const Real z = 0.);
94 
95 // Begin the main program.
96 int main (int argc, char ** argv)
97 {
98  // Initialize libMesh and any dependent libraries, like in example 2.
99  LibMeshInit init (argc, argv);
100 
101  // Only our PETSc interface currently supports solves restricted to
102  // subdomains
103  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
104 
105  // Skip adaptive examples on a non-adaptive libMesh build
106 #ifndef LIBMESH_ENABLE_AMR
107  libmesh_example_requires(false, "--enable-amr");
108 #else
109 
110  // Declare a performance log for the main program
111  // PerfLog perf_main("Main Program");
112 
113  // Create a GetPot object to parse the command line
114  GetPot command_line (argc, argv);
115 
116  // Check for proper calling arguments.
117  if (argc < 3)
118  {
119  // This handy function will print the file name, line number,
120  // specified message, and then throw an exception.
121  libmesh_error_msg("Usage:\n" << "\t " << argv[0] << " -d 2(3)" << " -n 15");
122  }
123 
124  // Brief message to the user regarding the program name
125  // and command line arguments.
126  else
127  {
128  libMesh::out << "Running " << argv[0];
129 
130  for (int i=1; i<argc; i++)
131  libMesh::out << " " << argv[i];
132 
133  libMesh::out << std::endl << std::endl;
134  }
135 
136 
137  // Read problem dimension from command line. Use int
138  // instead of unsigned since the GetPot overload is ambiguous
139  // otherwise.
140  int dim = 2;
141  if (command_line.search(1, "-d"))
142  dim = command_line.next(dim);
143 
144  // Skip higher-dimensional examples on a lower-dimensional libMesh build
145  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
146 
147  // Create a mesh with user-defined dimension.
148  // Read number of elements from command line
149  int ps = 15;
150  if (command_line.search(1, "-n"))
151  ps = command_line.next(ps);
152 
153  // Read FE order from command line
154  std::string order = "FIRST";
155  if (command_line.search(2, "-Order", "-o"))
156  order = command_line.next(order);
157 
158  // Read FE Family from command line
159  std::string family = "LAGRANGE";
160  if (command_line.search(2, "-FEFamily", "-f"))
161  family = command_line.next(family);
162 
163  // Cannot use discontinuous basis.
164  if ((family == "MONOMIAL") || (family == "XYZ"))
165  libmesh_error_msg("This example requires a C^0 (or higher) FE basis.");
166 
167  // Create a mesh, with dimension to be overridden later, on the
168  // default MPI communicator.
169  Mesh mesh(init.comm());
170 
171  // Use the MeshTools::Generation mesh generator to create a uniform
172  // grid on the square [-1,1]^D. We instruct the mesh generator
173  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
174  // elements in 3D. Building these higher-order elements allows
175  // us to use higher-order approximation, as in example 3.
176 
177  Real halfwidth = dim > 1 ? 1. : 0.;
178  Real halfheight = dim > 2 ? 1. : 0.;
179 
180  if ((family == "LAGRANGE") && (order == "FIRST"))
181  {
182  // No reason to use high-order geometric elements if we are
183  // solving with low-order finite elements.
185  ps,
186  (dim>1) ? ps : 0,
187  (dim>2) ? ps : 0,
188  -1., 1.,
189  -halfwidth, halfwidth,
190  -halfheight, halfheight,
191  (dim==1) ? EDGE2 :
192  ((dim == 2) ? QUAD4 : HEX8));
193  }
194 
195  else
196  {
198  ps,
199  (dim>1) ? ps : 0,
200  (dim>2) ? ps : 0,
201  -1., 1.,
202  -halfwidth, halfwidth,
203  -halfheight, halfheight,
204  (dim==1) ? EDGE3 :
205  ((dim == 2) ? QUAD9 : HEX27));
206  }
207 
208 
209  // To demonstrate solving on a subdomain, we will solve only on the
210  // interior of a circle (ball in 3d) with radius 0.8. So show that
211  // this also works well on locally refined meshes, we refine once
212  // all elements that are located on the boundary of this circle (or
213  // ball).
214  {
215  // A MeshRefinement object is needed to refine meshes.
216  MeshRefinement meshRefinement(mesh);
217 
218  // Loop over all elements.
219  for (auto & elem : mesh.element_ptr_range())
220  if (elem->active())
221  {
222  // Just check whether the current element has at least one
223  // node inside and one node outside the circle.
224  bool node_in = false;
225  bool node_out = false;
226  for (auto & n : elem->node_ref_range())
227  {
228  double d = n.norm();
229  if (d<0.8)
230  node_in = true;
231  else
232  node_out = true;
233  }
234  if (node_in && node_out)
235  elem->set_refinement_flag(Elem::REFINE);
236  else
237  elem->set_refinement_flag(Elem::DO_NOTHING);
238  }
239  else
240  elem->set_refinement_flag(Elem::INACTIVE);
241 
242  // Now actually refine.
243  meshRefinement.refine_elements();
244  }
245 
246  // Print information about the mesh to the screen.
247  mesh.print_info();
248 
249  // Now set the subdomain_id of all elements whose centroid is inside
250  // the circle to 1.
251  for (auto elem : mesh.element_ptr_range())
252  {
253  double d = elem->centroid().norm();
254  if (d < 0.8)
255  elem->subdomain_id() = 1;
256  }
257 
258  // Create an equation systems object.
259  EquationSystems equation_systems (mesh);
260 
261  // Declare the system and its variables.
262  // Create a system named "Poisson"
263  LinearImplicitSystem & system =
264  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
265 
266 
267  // Add the variable "u" to "Poisson". "u"
268  // will be approximated using second-order approximation.
269  system.add_variable("u",
270  Utility::string_to_enum<Order> (order),
271  Utility::string_to_enum<FEFamily>(family));
272 
273  // Give the system a pointer to the matrix assembly
274  // function.
276 
277  // Initialize the data structures for the equation system.
278  equation_systems.init();
279 
280  // Print information about the system to the screen.
281  equation_systems.print_info();
282  mesh.print_info();
283 
284  // Restrict solves to those elements that have subdomain_id set to 1.
285  std::set<subdomain_id_type> id_list;
286  id_list.insert(1);
288  SystemSubsetBySubdomain subset(system, selection);
289  system.restrict_solve_to(&subset, SUBSET_ZERO);
290 
291  // Note that using SUBSET_ZERO will cause all dofs outside the
292  // subdomain to be cleared. This will, however, cause some hanging
293  // nodes outside the subdomain to have inconsistent values.
294 
295  // Solve the system "Poisson", just like example 2.
296  equation_systems.get_system("Poisson").solve();
297 
298  // After solving the system write the solution
299  // to a GMV-formatted plot file.
300  if (dim == 1)
301  {
302  GnuPlotIO plot(mesh, "Subdomains Example 1, 1D", GnuPlotIO::GRID_ON);
303  plot.write_equation_systems("gnuplot_script", equation_systems);
304  }
305  else
306  {
307  GMVIO (mesh).write_equation_systems ((dim == 3) ?
308  "out_3.gmv" : "out_2.gmv", equation_systems);
309 #ifdef LIBMESH_HAVE_EXODUS_API
310  ExodusII_IO (mesh).write_equation_systems ((dim == 3) ?
311  "out_3.e" : "out_2.e", equation_systems);
312 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
313  }
314 
315 #endif // #ifndef LIBMESH_ENABLE_AMR
316 
317  // All done.
318  return 0;
319 }
320 
321 
322 
323 
324 // We now define the matrix assembly function for the
325 // Poisson system. We need to first compute element
326 // matrices and right-hand sides, and then take into
327 // account the boundary conditions, which will be handled
328 // via a penalty method.
330  const std::string & libmesh_dbg_var(system_name))
331 {
332  // It is a good idea to make sure we are assembling
333  // the proper system.
334  libmesh_assert_equal_to (system_name, "Poisson");
335 
336  // Declare a performance log. Give it a descriptive
337  // string to identify what part of the code we are
338  // logging, since there may be many PerfLogs in an
339  // application.
340  PerfLog perf_log ("Matrix Assembly");
341 
342  // Get a constant reference to the mesh object.
343  const MeshBase & mesh = es.get_mesh();
344 
345  // The dimension that we are running
346  const unsigned int dim = mesh.mesh_dimension();
347 
348  // Get a reference to the LinearImplicitSystem we are solving
349  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
350 
351  // A reference to the DofMap object for this system. The DofMap
352  // object handles the index translation from node and element numbers
353  // to degree of freedom numbers. We will talk more about the DofMap
354  // in future examples.
355  const DofMap & dof_map = system.get_dof_map();
356 
357  // Get a constant reference to the Finite Element type
358  // for the first (and only) variable in the system.
359  FEType fe_type = dof_map.variable_type(0);
360 
361  // Build a Finite Element object of the specified type. Since the
362  // FEBase::build() member dynamically creates memory we will
363  // store the object as a UniquePtr<FEBase>. This can be thought
364  // of as a pointer that will clean up after itself.
365  UniquePtr<FEBase> fe (FEBase::build(dim, fe_type));
366 
367  // A 5th order Gauss quadrature rule for numerical integration.
368  QGauss qrule (dim, FIFTH);
369 
370  // Tell the finite element object to use our quadrature rule.
371  fe->attach_quadrature_rule (&qrule);
372 
373  // Declare a special finite element object for
374  // boundary integration.
375  UniquePtr<FEBase> fe_face (FEBase::build(dim, fe_type));
376 
377  // Boundary integration requires one quadrature rule,
378  // with dimensionality one less than the dimensionality
379  // of the element.
380  QGauss qface(dim-1, FIFTH);
381 
382  // Tell the finite element object to use our
383  // quadrature rule.
384  fe_face->attach_quadrature_rule (&qface);
385 
386  // Here we define some references to cell-specific data that
387  // will be used to assemble the linear system.
388  // We begin with the element Jacobian * quadrature weight at each
389  // integration point.
390  const std::vector<Real> & JxW = fe->get_JxW();
391 
392  // The physical XY locations of the quadrature points on the element.
393  // These might be useful for evaluating spatially varying material
394  // properties at the quadrature points.
395  const std::vector<Point> & q_point = fe->get_xyz();
396 
397  // The element shape functions evaluated at the quadrature points.
398  const std::vector<std::vector<Real>> & phi = fe->get_phi();
399 
400  // The element shape function gradients evaluated at the quadrature
401  // points.
402  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
403 
404  // Define data structures to contain the element matrix
405  // and right-hand-side vector contribution. Following
406  // basic finite element terminology we will denote these
407  // "Ke" and "Fe". More detail is in example 3.
410 
411  // This vector will hold the degree of freedom indices for
412  // the element. These define where in the global system
413  // the element degrees of freedom get mapped.
414  std::vector<dof_id_type> dof_indices;
415 
416  // Now we will loop over all the elements in the mesh.
417  // We will compute the element matrix and right-hand-side
418  // contribution. See example 3 for a discussion of the
419  // element iterators.
420  for (const auto & elem : mesh.active_local_element_ptr_range())
421  {
422  // Elements with subdomain_id other than 1 are not in the active
423  // subdomain. We don't assemble anything for them.
424  if (elem->subdomain_id()==1)
425  {
426  // Start logging the shape function initialization.
427  // This is done through a simple function call with
428  // the name of the event to log.
429  perf_log.push("elem init");
430 
431  // Get the degree of freedom indices for the
432  // current element. These define where in the global
433  // matrix and right-hand-side this element will
434  // contribute to.
435  dof_map.dof_indices (elem, dof_indices);
436 
437  // Compute the element-specific data for the current
438  // element. This involves computing the location of the
439  // quadrature points (q_point) and the shape functions
440  // (phi, dphi) for the current element.
441  fe->reinit (elem);
442 
443  // Zero the element matrix and right-hand side before
444  // summing them. We use the resize member here because
445  // the number of degrees of freedom might have changed from
446  // the last element. Note that this will be the case if the
447  // element type is different (i.e. the last element was a
448  // triangle, now we are on a quadrilateral).
449  Ke.resize (dof_indices.size(),
450  dof_indices.size());
451 
452  Fe.resize (dof_indices.size());
453 
454  // Stop logging the shape function initialization.
455  // If you forget to stop logging an event the PerfLog
456  // object will probably catch the error and abort.
457  perf_log.pop("elem init");
458 
459  // Now we will build the element matrix. This involves
460  // a double loop to integrate the test functions (i) against
461  // the trial functions (j).
462  //
463  // We have split the numeric integration into two loops
464  // so that we can log the matrix and right-hand-side
465  // computation separately.
466  //
467  // Now start logging the element matrix computation
468  perf_log.push ("Ke");
469 
470  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
471  for (std::size_t i=0; i<phi.size(); i++)
472  for (std::size_t j=0; j<phi.size(); j++)
473  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
474 
475 
476  // Stop logging the matrix computation
477  perf_log.pop ("Ke");
478 
479  // Now we build the element right-hand-side contribution.
480  // This involves a single loop in which we integrate the
481  // "forcing function" in the PDE against the test functions.
482  //
483  // Start logging the right-hand-side computation
484  perf_log.push ("Fe");
485 
486  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
487  {
488  // fxy is the forcing function for the Poisson equation.
489  // In this case we set fxy to be a finite difference
490  // Laplacian approximation to the (known) exact solution.
491  //
492  // We will use the second-order accurate FD Laplacian
493  // approximation, which in 2D on a structured grid is
494  //
495  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
496  // u(i,j-1) + u(i,j+1) +
497  // -4*u(i,j))/h^2
498  //
499  // Since the value of the forcing function depends only
500  // on the location of the quadrature point (q_point[qp])
501  // we will compute it here, outside of the i-loop
502  const Real x = q_point[qp](0);
503 #if LIBMESH_DIM > 1
504  const Real y = q_point[qp](1);
505 #else
506  const Real y = 0.;
507 #endif
508 #if LIBMESH_DIM > 2
509  const Real z = q_point[qp](2);
510 #else
511  const Real z = 0.;
512 #endif
513  const Real eps = 1.e-3;
514 
515  const Real uxx = (exact_solution(x-eps, y, z) +
516  exact_solution(x+eps, y, z) +
517  -2.*exact_solution(x, y, z))/eps/eps;
518 
519  const Real uyy = (exact_solution(x, y-eps, z) +
520  exact_solution(x, y+eps, z) +
521  -2.*exact_solution(x, y, z))/eps/eps;
522 
523  const Real uzz = (exact_solution(x, y, z-eps) +
524  exact_solution(x, y, z+eps) +
525  -2.*exact_solution(x, y, z))/eps/eps;
526 
527  Real fxy;
528  if (dim==1)
529  {
530  // In 1D, compute the rhs by differentiating the
531  // exact solution twice.
532  const Real pi = libMesh::pi;
533  fxy = (0.25*pi*pi)*sin(.5*pi*x);
534  }
535  else
536  {
537  fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
538  }
539 
540  // Add the RHS contribution
541  for (std::size_t i=0; i<phi.size(); i++)
542  Fe(i) += JxW[qp]*fxy*phi[i][qp];
543  }
544 
545  // Stop logging the right-hand-side computation
546  perf_log.pop ("Fe");
547 
548  // At this point the interior element integration has
549  // been completed. However, we have not yet addressed
550  // boundary conditions. For this example we will only
551  // consider simple Dirichlet boundary conditions imposed
552  // via the penalty method. This is discussed at length in
553  // example 3.
554  {
555 
556  // Start logging the boundary condition computation
557  perf_log.push ("BCs");
558 
559  // The following loops over the sides of the element. If
560  // the element has no neighbor on a side then that side
561  // MUST live on a boundary of the domain. If there is a
562  // neighbor, check that neighbor's subdomain_id; if that
563  // is different from 1, the side is also located on the
564  // boundary.
565  for (auto side : elem->side_index_range())
566  if ((elem->neighbor_ptr(side) == libmesh_nullptr) ||
567  (elem->neighbor_ptr(side)->subdomain_id()!=1))
568  {
569 
570  // The penalty value. \frac{1}{\epsilon}
571  // in the discussion above.
572  const Real penalty = 1.e10;
573 
574  // The value of the shape functions at the quadrature
575  // points.
576  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
577 
578  // The Jacobian * Quadrature Weight at the quadrature
579  // points on the face.
580  const std::vector<Real> & JxW_face = fe_face->get_JxW();
581 
582  // The XYZ locations (in physical space) of the
583  // quadrature points on the face. This is where
584  // we will interpolate the boundary value function.
585  const std::vector<Point> & qface_point = fe_face->get_xyz();
586 
587  // Compute the shape function values on the element
588  // face.
589  fe_face->reinit(elem, side);
590 
591  // Loop over the face quadrature points for integration.
592  for (unsigned int qp=0; qp<qface.n_points(); qp++)
593  {
594  // The location on the boundary of the current
595  // face quadrature point.
596  const Real xf = qface_point[qp](0);
597 #if LIBMESH_DIM > 1
598  const Real yf = qface_point[qp](1);
599 #else
600  const Real yf = 0.;
601 #endif
602 #if LIBMESH_DIM > 2
603  const Real zf = qface_point[qp](2);
604 #else
605  const Real zf = 0.;
606 #endif
607 
608 
609  // The boundary value.
610  const Real value = exact_solution(xf, yf, zf);
611 
612  // Matrix contribution of the L2 projection.
613  for (std::size_t i=0; i<phi_face.size(); i++)
614  for (std::size_t j=0; j<phi_face.size(); j++)
615  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
616 
617  // Right-hand-side contribution of the L2
618  // projection.
619  for (std::size_t i=0; i<phi_face.size(); i++)
620  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
621  }
622  }
623 
624 
625  // Stop logging the boundary condition computation
626  perf_log.pop ("BCs");
627  }
628 
629  // If this assembly program were to be used on an adaptive mesh,
630  // we would have to apply any hanging node constraint equations
631  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
632 
633  // The element matrix and right-hand-side are now built
634  // for this element. Add them to the global matrix and
635  // right-hand-side vector. The SparseMatrix::add_matrix()
636  // and NumericVector::add_vector() members do this for us.
637  // Start logging the insertion of the local (element)
638  // matrix and vector into the global matrix and vector
639  perf_log.push ("matrix insertion");
640 
641  system.matrix->add_matrix (Ke, dof_indices);
642  system.rhs->add_vector (Fe, dof_indices);
643 
644  // Start logging the insertion of the local (element)
645  // matrix and vector into the global matrix and vector
646  perf_log.pop ("matrix insertion");
647  }
648  }
649 
650  // That's it. We don't need to do anything else to the
651  // PerfLog. When it goes out of scope (at this function return)
652  // it will print its log to the screen. Pretty easy, huh?
653 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
This is the EquationSystems class.
void pop(const char *label, const char *header="")
Pop the event label off the stack, resuming any lower event.
Definition: perf_log.C:162
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
Set dofs outside the subset to zero.
bool refine_elements()
Only refines the user-requested elements.
This class provides a specific system class.
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:350
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=libmesh_nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1101
unsigned short int side
Definition: xdr_io.C:49
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
Real exact_solution(const Real x, const Real y=0., const Real z=0.)
This is the exact solution that we are trying to obtain.
NumericVector< Number > * rhs
The system matrix.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:62
The libMesh namespace provides an interface to certain functionality in the library.
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
This class represents a subset of the dofs of a System, selected by the subdomain_id and possible the...
This is the MeshBase class.
Definition: mesh_base.h:68
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
int main(int argc, char **argv)
SolverPackage default_solver_package()
Definition: libmesh.C:995
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:125
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:167
This class implements writing meshes using GNUplot, designed for use only with 1D meshes...
Definition: gnuplot_io.h:43
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
PetscErrorCode Vec x
const DofMap & get_dof_map() const
Definition: system.h:2030
This is the MeshRefinement class.
virtual SimpleRange< element_iterator > element_ptr_range()=0
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
void attach_assemble_function(void fptr(EquationSystems &es, const std::string &name))
Register a user function to use in assembling the system matrix and RHS.
Definition: system.C:1795
void init(triangulateio &t)
Initializes the fields of t to NULL/0 as necessary.
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
void push(const char *label, const char *header="")
Push the event label onto the stack, pausing any active event.
Definition: perf_log.C:132
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) libmesh_override
After calling this method, any solve will be limited to the given subset.
const Parallel::Communicator & comm() const
Definition: libmesh.h:81
PetscErrorCode Vec Mat libmesh_dbg_var(j)
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=libmesh_nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void assemble_poisson(EquationSystems &es, const std::string &system_name)
OStreamProxy out
SparseMatrix< Number > * matrix
The system matrix.
unsigned int mesh_dimension() const
Definition: mesh_base.C:148
static const bool value
Definition: xdr_io.C:108
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:776
This class implements specific orders of Gauss quadrature.
const MeshBase & get_mesh() const
const T_sys & get_system(const std::string &name) const
virtual void init()
Initialize all the systems.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:448
unsigned int n_points() const
Definition: quadrature.h:113
void build_cube(UnstructuredMesh &mesh, const unsigned int nx=0, const unsigned int ny=0, const unsigned int nz=0, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const Real zmin=0., const Real zmax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
Builds a (elements) cube.
const Real pi
.
Definition: libmesh.h:172
static UniquePtr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.