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