libMesh
miscellaneous_ex16.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 // <h1>Miscellaneous Example 16 - Static Condensation</h1>
19 // \author Alexander D. Lindsay
20 // \date 2024
21 //
22 // This example is equivalent to the third introductory example except that
23 // we add and use a StaticCondensation object. The StaticCondensation class
24 // forward eliminates internal degrees of freedom and then performs a global
25 // solve only with the trace degrees of freedom. We then perform backward
26 // substitution to recover the solution for the internal degrees of freedom.
27 // We verify that this solution matches the solution when the global solve
28 // is performed on both interior and trace degrees of freedom
29 
30 #include "libmesh/libmesh_config.h"
31 
32 // C++ include files that we need
33 #include <iostream>
34 #include <algorithm>
35 #include <math.h>
36 
37 // Basic include files needed for the mesh functionality.
38 #include "libmesh/libmesh.h"
39 #include "libmesh/mesh.h"
40 #include "libmesh/mesh_generation.h"
41 #include "libmesh/linear_implicit_system.h"
42 #include "libmesh/equation_systems.h"
43 
44 // Define the Finite Element object.
45 #include "libmesh/fe.h"
46 
47 // Define Gauss quadrature rules.
48 #include "libmesh/quadrature_gauss.h"
49 
50 // Define useful datatypes for finite element
51 // matrix and vector components.
52 #include "libmesh/sparse_matrix.h"
53 #include "libmesh/numeric_vector.h"
54 #include "libmesh/dense_matrix.h"
55 #include "libmesh/dense_vector.h"
56 #include "libmesh/elem.h"
57 #include "libmesh/enum_solver_package.h"
58 #include "libmesh/static_condensation.h"
59 
60 // Define the DofMap, which handles degree of freedom
61 // indexing.
62 #include "libmesh/dof_map.h"
63 
64 // For mesh refinement
65 #include "libmesh/mesh_refinement.h"
66 #include "libmesh/error_vector.h"
67 #include "libmesh/kelly_error_estimator.h"
68 
69 // I/O utilities.
70 #include "libmesh/getpot.h"
71 #include "libmesh/exodusII_io.h"
72 
73 // For the solver for the system with static condensation
74 #include "libmesh/petsc_linear_solver.h"
75 
76 // Bring in everything from the libMesh namespace
77 using namespace libMesh;
78 
79 #if defined(LIBMESH_HAVE_EIGEN_DENSE) && defined(LIBMESH_HAVE_PETSC)
80 // Function prototype. This is the function that will assemble
81 // the linear system for our Poisson problem. Note that the
82 // function will take the EquationSystems object and the
83 // name of the system we are assembling as input. From the
84 // EquationSystems object we have access to the Mesh and
85 // other objects we might need.
86 void assemble_poisson(EquationSystems & es, const std::string & system_name);
87 #endif
88 
89 // Function prototype for the exact solution.
90 Real exact_solution(const Real x, const Real y, const Real z = 0.);
91 
92 int
93 main(int argc, char ** argv)
94 {
95  // Initialize libraries, like in example 2.
96  LibMeshInit init(argc, argv);
97 
98 #if !defined(LIBMESH_HAVE_EIGEN_DENSE)
99  libmesh_example_requires(false, "--enable-eigen");
100 #elif !defined(LIBMESH_HAVE_PETSC)
101  libmesh_example_requires(false, "--enable-petsc");
102 #else
103 
104  // Brief message to the user regarding the program name
105  // and command line arguments.
106  libMesh::out << "Running " << argv[0];
107 
108  for (int i = 1; i < argc; i++)
109  libMesh::out << " " << argv[i];
110  libMesh::out << std::endl << std::endl;
111 
112  // Skip this 2D example if libMesh was compiled as 1D-only.
113  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
114 
115  // Parse the input file.
116  GetPot infile("miscellaneous_ex16.in");
117 
118  // But allow the command line to override it.
119  infile.parse_command_line(argc, argv);
120 
121  const unsigned int grid_size = infile("grid_size", 5);
122 
123  // Create a mesh, with dimension to be overridden later, distributed
124  // across the default MPI communicator.
125  Mesh mesh(init.comm());
126 
127  // Use the MeshTools::Generation mesh generator to create a uniform
128  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
129  // to build a mesh of 15x15 QUAD9 elements. Building QUAD9
130  // elements instead of the default QUAD4's we used in example 2
131  // allow us to use higher-order approximation.
132  MeshTools::Generation::build_square(mesh, grid_size, grid_size, -1., 1., -1., 1., QUAD9);
133 
134  // Print information about the mesh to the screen.
135  // Note that 5x5 QUAD9 elements actually has 11x11 nodes,
136  // so this mesh is significantly larger than the one in example 2.
137  mesh.print_info();
138 
139  // Create an equation systems object.
140  EquationSystems equation_systems(mesh);
141 
142  // Declare the Poisson system and its variables.
143  // The Poisson system is another example of a steady system.
144  auto & sys = equation_systems.add_system<LinearImplicitSystem>("Poisson");
145 
146  // Adds the variable "u" to "Poisson". "u"
147  // will be approximated using second-order approximation.
148  sys.add_variable("u", SECOND);
149 
150  // Give the system a pointer to the matrix assembly
151  // function. This will be called when needed by the
152  // library.
153  sys.attach_assemble_function(assemble_poisson);
154 
155  // Now perform same steps for system with static condensation enabled
156  auto & sc_sys = equation_systems.add_system<LinearImplicitSystem>("SC_Poisson");
157  sc_sys.add_variable("v", SECOND);
158  sc_sys.attach_assemble_function(assemble_poisson);
159  sc_sys.create_static_condensation();
160 
161 #ifdef LIBMESH_ENABLE_AMR
162  // Define the mesh refinement object that takes care of adaptively
163  // refining the mesh.
164  MeshRefinement mesh_refinement(mesh);
165 
166  // These parameters determine the proportion of elements that will
167  // be refined and coarsened. Any element within 30% of the maximum
168  // error on any element will be refined, and any element within 30%
169  // of the minimum error on any element might be coarsened
170  mesh_refinement.refine_fraction() = 0.7;
171  mesh_refinement.coarsen_fraction() = 0.3;
172  // We won't refine any element more than 2 times in total
173  mesh_refinement.max_h_level() = 2;
174 #endif
175 
176  // Initialize the data structures for the equation system.
177  equation_systems.init();
178 
179 #ifdef LIBMESH_ENABLE_AMR
180  // Refinement parameters
181  const unsigned int max_r_steps = 2; // Refine the mesh 2 times
182 
183  for (const auto r_step : make_range(max_r_steps + 1))
184  {
185 #endif
186  // Prints information about the system to the screen.
187  equation_systems.print_info();
188 
189  // Solve
190  sys.solve();
191  auto * sc_solver = dynamic_cast<PetscLinearSolver<Number> *>(sc_sys.get_linear_solver());
192  libmesh_assert(sc_solver);
193  KSP sc_ksp = sc_solver->ksp();
194  LibmeshPetscCall2(sc_solver->comm(), KSPSetType(sc_ksp, KSPPREONLY));
195  LibmeshPetscCall2(sc_solver->comm(), KSPSetInitialGuessNonzero(sc_ksp, PETSC_FALSE));
196  sc_sys.solve();
197 
198  libmesh_error_msg_if(!libMesh::relative_fuzzy_equals(*sys.solution, *sc_sys.solution, 1e-4),
199  "mismatching solution");
200  libMesh::out << "Static condensation reduced problem size to "
201  << sc_sys.get_static_condensation().get_condensed_mat().m() << std::endl << std::endl;
202 
203 #if defined(LIBMESH_HAVE_EXODUS_API) && !defined(LIBMESH_ENABLE_PARMESH)
204  // After solving the system write the solution
205  // to an Exodus-formatted file.
206  ExodusII_IO exii_io(mesh);
207  const std::string file_name =
208 #ifdef LIBMESH_ENABLE_AMR
209  "out_" + std::to_string(r_step) + ".e";
210 #else
211  "out.e";
212 #endif
213  exii_io.write_equation_systems(file_name, equation_systems);
214 #endif
215 
216 #ifdef LIBMESH_ENABLE_AMR
217  // We need to ensure that the mesh is not refined on the last iteration
218  // of this loop, since we do not want to refine the mesh unless we are
219  // going to solve the equation system for that refined mesh.
220  if (r_step != max_r_steps)
221  {
222  // Error estimation objects, see Adaptivity Example 2 for details
223  ErrorVector error;
224  KellyErrorEstimator error_estimator;
225  error_estimator.use_unweighted_quadrature_rules = true;
226 
227  // Compute the error for each active element
228  error_estimator.estimate_error(sys, error);
229 
230  // Output error estimate magnitude
231  libMesh::out << "Error estimate\nl2 norm = " << error.l2_norm()
232  << "\nmaximum = " << error.maximum() << std::endl << std::endl;
233 
234  // Flag elements to be refined and coarsened
235  mesh_refinement.flag_elements_by_error_fraction(error);
236 
237  // Perform refinement and coarsening
238  mesh_refinement.refine_and_coarsen_elements();
239 
240  // Reinitialize the equation_systems object for the newly refined
241  // mesh. One of the steps in this is project the solution onto the
242  // new mesh
243  equation_systems.reinit();
244  }
245  }
246 #endif // LIBMESH_ENABLE_AMR
247 #endif // defined(LIBMESH_HAVE_EIGEN_DENSE) && defined(LIBMESH_HAVE_PETSC)
248 
249  // All done.
250  return 0;
251 }
252 
253 #if defined(LIBMESH_HAVE_EIGEN_DENSE) && defined(LIBMESH_HAVE_PETSC)
254 
255 // We now define the matrix assembly function for the
256 // Poisson system. We need to first compute element
257 // matrices and right-hand sides, and then take into
258 // account the boundary conditions, which will be handled
259 // via a penalty method.
260 void
261 assemble_poisson(EquationSystems & es, const std::string & system_name)
262 {
263  // Get a constant reference to the mesh object.
264  const MeshBase & mesh = es.get_mesh();
265 
266  // The dimension that we are running
267  const unsigned int dim = mesh.mesh_dimension();
268 
269  // Get a reference to the LinearImplicitSystem we are solving
270  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>(system_name);
271 
272  // Get a pointer to the StaticCondensation class if it exists
273  StaticCondensation * sc = nullptr;
274  if (system.has_static_condensation())
275  sc = &system.get_static_condensation();
276 
277  // A reference to the DofMap object for this system. The DofMap
278  // object handles the index translation from node and element numbers
279  // to degree of freedom numbers. We will talk more about the DofMap
280  // in future examples.
281  const DofMap & dof_map = system.get_dof_map();
282 
283  // Get a constant reference to the Finite Element type
284  // for the first (and only) variable in the system.
285  FEType fe_type = dof_map.variable_type(0);
286 
287  // Build a Finite Element object of the specified type. Since the
288  // FEBase::build() member dynamically creates memory we will
289  // store the object as a std::unique_ptr<FEBase>. This can be thought
290  // of as a pointer that will clean up after itself. Introduction Example 4
291  // describes some advantages of std::unique_ptr's in the context of
292  // quadrature rules.
293  std::unique_ptr<FEBase> fe(FEBase::build(dim, fe_type));
294 
295  // A 5th order Gauss quadrature rule for numerical integration.
296  QGauss qrule(dim, FIFTH);
297 
298  // Tell the finite element object to use our quadrature rule.
299  fe->attach_quadrature_rule(&qrule);
300 
301  // Declare a special finite element object for
302  // boundary integration.
303  std::unique_ptr<FEBase> fe_face(FEBase::build(dim, fe_type));
304 
305  // Boundary integration requires one quadrature rule,
306  // with dimensionality one less than the dimensionality
307  // of the element.
308  QGauss qface(dim - 1, FIFTH);
309 
310  // Tell the finite element object to use our
311  // quadrature rule.
312  fe_face->attach_quadrature_rule(&qface);
313 
314  // Here we define some references to cell-specific data that
315  // will be used to assemble the linear system.
316  //
317  // The element Jacobian * quadrature weight at each integration point.
318  const std::vector<Real> & JxW = fe->get_JxW();
319 
320  // The physical XY locations of the quadrature points on the element.
321  // These might be useful for evaluating spatially varying material
322  // properties at the quadrature points.
323  const std::vector<Point> & q_point = fe->get_xyz();
324 
325  // The element shape functions evaluated at the quadrature points.
326  const std::vector<std::vector<Real>> & phi = fe->get_phi();
327 
328  // The element shape function gradients evaluated at the quadrature
329  // points.
330  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
331 
332  // Define data structures to contain the element matrix
333  // and right-hand-side vector contribution. Following
334  // basic finite element terminology we will denote these
335  // "Ke" and "Fe". These datatypes are templated on
336  // Number, which allows the same code to work for real
337  // or complex numbers.
340 
341  // This vector will hold the degree of freedom indices for
342  // the element. These define where in the global system
343  // the element degrees of freedom get mapped.
344  std::vector<dof_id_type> dof_indices;
345 
346  // The global system matrix
347  SparseMatrix<Number> & matrix = system.get_system_matrix();
348 
349  // Now we will loop over all the elements in the mesh.
350  // We will compute the element matrix and right-hand-side
351  // contribution.
352  //
353  // Element ranges are a nice way to iterate through all the
354  // elements, or all the elements that have some property. The
355  // range will iterate from the first to the last element on
356  // the local processor.
357  // It is smart to make this one const so that we don't accidentally
358  // mess it up! In case users later modify this program to include
359  // refinement, we will be safe and will only consider the active
360  // elements; hence we use a variant of the
361  // active_local_element_ptr_range.
362  for (const auto & elem : mesh.active_local_element_ptr_range())
363  {
364  // Get the degree of freedom indices for the
365  // current element. These define where in the global
366  // matrix and right-hand-side this element will
367  // contribute to.
368  dof_map.dof_indices(elem, dof_indices);
369 
370  // Cache the number of degrees of freedom on this element, for
371  // use as a loop bound later. We use cast_int to explicitly
372  // convert from size() (which may be 64-bit) to unsigned int
373  // (which may be 32-bit but which is definitely enough to count
374  // *local* degrees of freedom.
375  const unsigned int n_dofs = cast_int<unsigned int>(dof_indices.size());
376 
377  // Compute the element-specific data for the current
378  // element. This involves computing the location of the
379  // quadrature points (q_point) and the shape functions
380  // (phi, dphi) for the current element.
381  fe->reinit(elem);
382 
383  // With one variable, we should have the same number of degrees
384  // of freedom as shape functions.
385  libmesh_assert_equal_to(n_dofs, phi.size());
386 
387  // Zero the element matrix and right-hand side before
388  // summing them. We use the resize member here because
389  // the number of degrees of freedom might have changed from
390  // the last element. Note that this will be the case if the
391  // element type is different (i.e. the last element was a
392  // triangle, now we are on a quadrilateral).
393 
394  // The DenseMatrix::resize() and the DenseVector::resize()
395  // members will automatically zero out the matrix and vector.
396  Ke.resize(n_dofs, n_dofs);
397 
398  Fe.resize(n_dofs);
399 
400  // Now loop over the quadrature points. This handles
401  // the numeric integration.
402  for (unsigned int qp = 0; qp < qrule.n_points(); qp++)
403  {
404 
405  // Now we will build the element matrix. This involves
406  // a double loop to integrate the test functions (i) against
407  // the trial functions (j).
408  for (unsigned int i = 0; i != n_dofs; i++)
409  for (unsigned int j = 0; j != n_dofs; j++)
410  {
411  Ke(i, j) += JxW[qp] * (dphi[i][qp] * dphi[j][qp]);
412  }
413 
414  // This is the end of the matrix summation loop
415  // Now we build the element right-hand-side contribution.
416  // This involves a single loop in which we integrate the
417  // "forcing function" in the PDE against the test functions.
418  {
419  const Real x = q_point[qp](0);
420  const Real y = q_point[qp](1);
421  const Real eps = 1.e-3;
422 
423  // "fxy" is the forcing function for the Poisson equation.
424  // In this case we set fxy to be a finite difference
425  // Laplacian approximation to the (known) exact solution.
426  //
427  // We will use the second-order accurate FD Laplacian
428  // approximation, which in 2D is
429  //
430  // u_xx + u_yy = (u(i,j-1) + u(i,j+1) +
431  // u(i-1,j) + u(i+1,j) +
432  // -4*u(i,j))/h^2
433  //
434  // Since the value of the forcing function depends only
435  // on the location of the quadrature point (q_point[qp])
436  // we will compute it here, outside of the i-loop
437  const Real fxy =
438  -(exact_solution(x, y - eps) + exact_solution(x, y + eps) + exact_solution(x - eps, y) +
439  exact_solution(x + eps, y) - 4. * exact_solution(x, y)) /
440  eps / eps;
441 
442  for (unsigned int i = 0; i != n_dofs; i++)
443  Fe(i) += JxW[qp] * fxy * phi[i][qp];
444  }
445  }
446 
447  // We have now reached the end of the RHS summation,
448  // and the end of quadrature point loop, so
449  // the interior element integration has
450  // been completed. However, we have not yet addressed
451  // boundary conditions. For this example we will only
452  // consider simple Dirichlet boundary conditions.
453  //
454  // There are several ways Dirichlet boundary conditions
455  // can be imposed. A simple approach, which works for
456  // interpolary bases like the standard Lagrange polynomials,
457  // is to assign function values to the
458  // degrees of freedom living on the domain boundary. This
459  // works well for interpolary bases, but is more difficult
460  // when non-interpolary (e.g Legendre or Hierarchic) bases
461  // are used.
462  //
463  // Dirichlet boundary conditions can also be imposed with a
464  // "penalty" method. In this case essentially the L2 projection
465  // of the boundary values are added to the matrix. The
466  // projection is multiplied by some large factor so that, in
467  // floating point arithmetic, the existing (smaller) entries
468  // in the matrix and right-hand-side are effectively ignored.
469  //
470  // This amounts to adding a term of the form (in latex notation)
471  //
472  // \frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta
473  // \Omega} u \phi_i
474  //
475  // where
476  //
477  // \frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
478  {
479 
480  // The following loop is over the sides of the element.
481  // If the element has no neighbor on a side then that
482  // side MUST live on a boundary of the domain.
483  for (auto side : elem->side_index_range())
484  if (elem->neighbor_ptr(side) == nullptr)
485  {
486  // The value of the shape functions at the quadrature
487  // points.
488  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
489 
490  // The Jacobian * Quadrature Weight at the quadrature
491  // points on the face.
492  const std::vector<Real> & JxW_face = fe_face->get_JxW();
493 
494  // The XYZ locations (in physical space) of the
495  // quadrature points on the face. This is where
496  // we will interpolate the boundary value function.
497  const std::vector<Point> & qface_point = fe_face->get_xyz();
498 
499  // Compute the shape function values on the element
500  // face.
501  fe_face->reinit(elem, side);
502 
503  // Some shape functions will be 0 on the face, but for
504  // ease of indexing and generality of code we loop over
505  // them anyway
506  libmesh_assert_equal_to(n_dofs, phi_face.size());
507 
508  // Loop over the face quadrature points for integration.
509  for (unsigned int qp = 0; qp < qface.n_points(); qp++)
510  {
511  // The location on the boundary of the current
512  // face quadrature point.
513  const Real xf = qface_point[qp](0);
514  const Real yf = qface_point[qp](1);
515 
516  // The penalty value. \frac{1}{\epsilon}
517  // in the discussion above.
518  const Real penalty = 1.e10;
519 
520  // The boundary value.
521  const Real value = exact_solution(xf, yf);
522 
523  // Matrix contribution of the L2 projection.
524  for (unsigned int i = 0; i != n_dofs; i++)
525  for (unsigned int j = 0; j != n_dofs; j++)
526  Ke(i, j) += JxW_face[qp] * penalty * phi_face[i][qp] * phi_face[j][qp];
527 
528  // Right-hand-side contribution of the L2
529  // projection.
530  for (unsigned int i = 0; i != n_dofs; i++)
531  Fe(i) += JxW_face[qp] * penalty * value * phi_face[i][qp];
532  }
533  }
534  }
535 
536  // We have now finished the quadrature point loop,
537  // and have therefore applied all the boundary conditions.
538 
539  // If this assembly program were to be used on an adaptive mesh,
540  // we would have to apply any hanging node constraint equations
541  dof_map.constrain_element_matrix_and_vector(Ke, Fe, dof_indices);
542 
543  if (sc)
544  sc->set_current_elem(*elem);
545 
546  // The element matrix and right-hand-side are now built
547  // for this element. Add them to the global matrix and
548  // right-hand-side vector. The SparseMatrix::add_matrix()
549  // and NumericVector::add_vector() members do this for us.
550  matrix.add_matrix(Ke, dof_indices);
551  system.rhs->add_vector(Fe, dof_indices);
552  }
553 
554  matrix.close();
555 }
556 
557 #endif // defined(LIBMESH_HAVE_EIGEN_DENSE) && defined(LIBMESH_HAVE_PETSC)
virtual T maximum() const
Definition: statistics.C:62
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
This is the EquationSystems class.
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
void constrain_element_matrix_and_vector(DenseMatrix< Number > &matrix, DenseVector< Number > &rhs, std::vector< dof_id_type > &elem_dofs, bool asymmetric_constraint_rows=true) const
Constrains the element matrix and vector.
Definition: dof_map.h:2330
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:2164
virtual Real l2_norm() const
Definition: statistics.C:37
unsigned int dim
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
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]...
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:2220
Real & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
bool has_static_condensation() const
Definition: system.C:2871
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.
const T_sys & get_system(std::string_view name) const
unsigned int & max_h_level()
The max_h_level is the greatest refinement level an element should reach.
This is the MeshBase class.
Definition: mesh_base.h:75
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
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.
int main(int argc, char **argv)
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
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
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.
libmesh_assert(ctx)
StaticCondensation & get_static_condensation()
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
unsigned int n_points() const
Definition: quadrature.h:131
void assemble_poisson(EquationSystems &es, const std::string &system_name)
This class implements the Kelly error indicator which is based on the flux jumps between elements...
virtual void close()=0
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
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
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:140
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
void flag_elements_by_error_fraction(const ErrorVector &error_per_cell, const Real refine_fraction=0.3, const Real coarsen_fraction=0.0, const unsigned int max_level=libMesh::invalid_uint)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell...
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class&#39;s jump error estimate formula to estimate the error on each cell...
bool relative_fuzzy_equals(const T &var1, const T2 &var2, const Real tol=TOLERANCE *TOLERANCE)
Function to check whether two variables are equal within a relative tolerance.
Definition: fuzzy_equals.h:78
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