Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
vector_fe_ex10.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 // <h1>Vector Finite Elements Example 10 - Raviart-Thomas elements (grad-div)</h1>
20 // \author Nuno Nobre
21 // \date 2025
22 //
23 // This example uses Raviart-Thomas elements to solve a model grad-div problem
24 // in H(div) in both 2d and 3d: -\nabla (\nabla \cdot \vec{u}) + \vec{u} = f.
25 // Note that, unlike the problem in Vector Finite Elements Example 6, this is
26 // _not_ an elliptic problem.
27 
28 // Basic utilities.
29 #include "libmesh/string_to_enum.h"
30 
31 // The solver packages supported by libMesh.
32 #include "libmesh/enum_solver_package.h"
33 
34 // The mesh object and mesh generation and modification utilities.
35 #include "libmesh/mesh.h"
36 #include "libmesh/mesh_generation.h"
37 #include "libmesh/mesh_modification.h"
38 
39 // Matrix and vector types.
40 #include "libmesh/dense_matrix.h"
41 #include "libmesh/sparse_matrix.h"
42 #include "libmesh/dense_vector.h"
43 #include "libmesh/numeric_vector.h"
44 
45 // The finite element object and the geometric element type.
46 #include "libmesh/fe.h"
47 #include "libmesh/elem.h"
48 
49 // Gauss quadrature rules.
50 #include "libmesh/quadrature_gauss.h"
51 
52 // The dof map, which handles degree of freedom indexing.
53 #include "libmesh/dof_map.h"
54 
55 // The system of equations.
56 #include "libmesh/equation_systems.h"
57 #include "libmesh/linear_implicit_system.h"
58 
59 // The exact solution and error computation.
60 #include "libmesh/exact_solution.h"
61 #include "libmesh/enum_norm_type.h"
62 #include "solution_function.h"
63 
64 // I/O utilities.
65 #include "libmesh/getpot.h"
66 #include "libmesh/exodusII_io.h"
67 
68 
69 // Bring in everything from the libMesh namespace.
70 using namespace libMesh;
71 
72 // Function prototype. This is the function that will assemble
73 // the linear system for our grad-div problem. Note that the
74 // function will take the EquationSystems object and the
75 // name of the system we are assembling as input. From the
76 // EquationSystems object we have access to the Mesh and
77 // other objects we might need.
79  const std::string & system_name);
80 
81 int main (int argc, char ** argv)
82 {
83  // Initialize libMesh.
84  LibMeshInit init (argc, argv);
85 
86  // This example requires a linear solver package.
87  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
88  "--enable-petsc, --enable-trilinos, or --enable-eigen");
89 
90  // Parse the input file.
91  GetPot infile("vector_fe_ex10.in");
92 
93  // But allow the command line to override it.
94  infile.parse_command_line(argc, argv);
95 
96  // Read in parameters from the command line and the input file.
97  const unsigned int dimension = infile("dim", 2);
98  const unsigned int grid_size = infile("grid_size", 15);
99 
100  // Skip higher-dimensional examples on a lower-dimensional libMesh build.
101  libmesh_example_requires(dimension <= LIBMESH_DIM, dimension << "D support");
102 
103  // Create a mesh, with dimension to be overridden later, distributed
104  // across the default MPI communicator.
105  Mesh mesh(init.comm());
106 
107  // Use the MeshTools::Generation mesh generator to create a uniform
108  // grid on the cube [-1,1]^D. To accomodate Raviart-Thomas elements, we must
109  // use TRI6/7 or QUAD8/9 elements in 2d, or TET14 or HEX27 in 3d.
110  const std::string elem_str = infile("element_type", std::string("TRI6"));
111 
112  libmesh_error_msg_if((dimension == 2 && elem_str != "TRI6" && elem_str != "TRI7" && elem_str != "QUAD8" && elem_str != "QUAD9") ||
113  (dimension == 3 && elem_str != "TET14" && elem_str != "HEX27"),
114  "You selected " << elem_str <<
115  " but this example must be run with TRI6, TRI7, QUAD8, or QUAD9 in 2d" <<
116  " or with TET14, or HEX27 in 3d.");
117 
118  if (dimension == 2)
120  grid_size,
121  grid_size,
122  -1., 1.,
123  -1., 1.,
124  Utility::string_to_enum<ElemType>(elem_str));
125  else if (dimension == 3)
127  grid_size,
128  grid_size,
129  grid_size,
130  -1., 1.,
131  -1., 1.,
132  -1., 1.,
133  Utility::string_to_enum<ElemType>(elem_str));
134 
135  // Make sure the code is robust against nodal reorderings.
137 
138  // Print information about the mesh to the screen.
139  mesh.print_info();
140 
141  // Create an equation systems object.
142  EquationSystems equation_systems (mesh);
143 
144  // Declare the system "GradDiv" and its variable.
145  LinearImplicitSystem & system = equation_systems.add_system<LinearImplicitSystem>("GradDiv");
146 
147  // Set the FE approximation order for the vector field variable.
148  const Order vector_order = static_cast<Order>(infile("order", 1u));
149 
150  libmesh_error_msg_if(vector_order < FIRST || vector_order > ((dimension == 3) ? FIRST : FIFTH),
151  "You selected: " << vector_order <<
152  " but this example must be run with either 1 <= order <= 5 in 2d"
153  " or with order 1 in 3d.");
154 
155  // Adds the variable "u" to "GradDiv". "u" will be our vector field.
156  system.add_variable("u", vector_order, RAVIART_THOMAS);
157 
158  // Give the system a pointer to the matrix assembly
159  // function. This will be called when needed by the library.
161 
162  // Initialize the data structures for the equation system.
163  equation_systems.init();
164 
165  // Prints information about the system to the screen.
166  equation_systems.print_info();
167 
168  // Solve the system "GradDiv". Note that calling this
169  // member will assemble the linear system and invoke
170  // the default numerical solver.
171  system.solve();
172 
173  ExactSolution exact_sol(equation_systems);
174 
175  SolutionFunction soln_func;
176  SolutionGradient soln_grad;
177 
178  // Build FunctionBase* containers to attach to the ExactSolution object.
179  std::vector<FunctionBase<Number> *> sols(1, &soln_func);
180  std::vector<FunctionBase<Gradient> *> grads(1, &soln_grad);
181 
182  exact_sol.attach_exact_values(sols);
183  exact_sol.attach_exact_derivs(grads);
184 
185  // Use higher quadrature order for more accurate error results.
186  int extra_error_quadrature = infile("extra_error_quadrature", 2);
187  exact_sol.extra_quadrature_order(extra_error_quadrature);
188 
189  // Compute the error.
190  exact_sol.compute_error("GradDiv", "u");
191 
192  // Print out the error values.
193  libMesh::out << "~~ Vector field (u) ~~"
194  << std::endl;
195  libMesh::out << "L2 error is: "
196  << exact_sol.l2_error("GradDiv", "u")
197  << std::endl;
198  libMesh::out << "HDiv semi-norm error is: "
199  << exact_sol.error_norm("GradDiv", "u", HDIV_SEMINORM)
200  << std::endl;
201  libMesh::out << "HDiv error is: "
202  << exact_sol.hdiv_error("GradDiv", "u")
203  << std::endl;
204 
205 #ifdef LIBMESH_HAVE_EXODUS_API
206 
207  // We write the file in the ExodusII format.
208  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
209 
210 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
211 
212  // All done.
213  return 0;
214 }
215 
216 
217 
218 // We now define the matrix assembly function for the
219 // grad-div system. We need to first compute element
220 // matrices and right-hand sides, and then take into
221 // account the boundary conditions, which will be handled
222 // via a penalty method.
224  const std::string & libmesh_dbg_var(system_name))
225 {
226 
227  // It is a good idea to make sure we are assembling
228  // the proper system.
229  libmesh_assert_equal_to (system_name, "GradDiv");
230 
231  // Get a constant reference to the mesh object.
232  const MeshBase & mesh = es.get_mesh();
233 
234  // The dimension that we are running.
235  const unsigned int dim = mesh.mesh_dimension();
236 
237  // Get a reference to the LinearImplicitSystem we are solving.
238  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("GradDiv");
239 
240  // A reference to the DofMap object for this system. The DofMap
241  // object handles the index translation from node and element numbers
242  // to degree of freedom numbers.
243  const DofMap & dof_map = system.get_dof_map();
244 
245  // Get a constant reference to the Finite Element type
246  // for the variable in the system.
247  FEType vector_fe_type = dof_map.variable_type(system.variable_number("u"));
248 
249  // Build the Finite Element object. Since the
250  // FEBase::build() member dynamically creates memory we will
251  // store the object as a std::unique_ptr<FEBase>. This can be thought
252  // of as a pointer that will clean up after itself. Introduction Example 4
253  // describes some advantages of std::unique_ptr's in the context of
254  // quadrature rules.
255  std::unique_ptr<FEVectorBase> vector_fe (FEVectorBase::build(dim, vector_fe_type));
256 
257  // A just-high-enough Gauss quadrature rule for numerical integration.
258  QGauss qrule (dim, vector_fe_type.default_quadrature_order());
259 
260  // Tell the finite element object to use our quadrature rule.
261  vector_fe->attach_quadrature_rule (&qrule);
262 
263  // Declare a special finite element object for boundary integration.
264  std::unique_ptr<FEVectorBase> vector_fe_face (FEVectorBase::build(dim, vector_fe_type));
265 
266  // Boundary integration requires one quadrature rule with dimensionality one
267  // less than the dimensionality of the element.
268  QGauss qface(dim-1, vector_fe_type.default_quadrature_order());
269 
270  // Tell the finite element object to use our quadrature rule.
271  vector_fe_face->attach_quadrature_rule (&qface);
272 
273  // Here we define some references to cell-specific data that
274  // will be used to assemble the linear system.
275  //
276  // The element Jacobian * quadrature weight at each integration point.
277  const std::vector<Real> & JxW = vector_fe->get_JxW();
278 
279  // The physical XY locations of the quadrature points on the element.
280  // These might be useful for evaluating spatially varying material
281  // properties at the quadrature points.
282  const std::vector<Point> & q_point = vector_fe->get_xyz();
283 
284  // The element shape functions evaluated at the quadrature points.
285  const std::vector<std::vector<RealGradient>> & vector_phi = vector_fe->get_phi();
286 
287  // The divergence of the element vector shape functions evaluated at the
288  // quadrature points.
289  const std::vector<std::vector<Real>> & div_vector_phi = vector_fe->get_div_phi();
290 
291  // Define data structures to contain the element matrix
292  // and right-hand-side vector contribution. Following
293  // basic finite element terminology we will denote these
294  // "Ke" and "Fe". These datatypes are templated on
295  // Number, which allows the same code to work for real
296  // or complex numbers.
299 
300  // These vectors will hold the degree of freedom indices for
301  // the element. These define where in the global system
302  // the element degrees of freedom get mapped.
303  std::vector<dof_id_type> dof_indices;
304 
305  // The global system matrix
306  SparseMatrix<Number> & matrix = system.get_system_matrix();
307 
308  // Now we will loop over all the elements in the mesh.
309  // We will compute the element matrix and right-hand-side
310  // contribution.
311  //
312  // Element ranges are a nice way to iterate through all the
313  // elements, or all the elements that have some property. The
314  // range will iterate from the first to the last element on
315  // the local processor.
316  // It is smart to make this one const so that we don't accidentally
317  // mess it up! In case users later modify this program to include
318  // refinement, we will be safe and will only consider the active
319  // elements; hence we use a variant of the
320  // active_local_element_ptr_range.
321  for (const auto & elem : mesh.active_local_element_ptr_range())
322  {
323  // Get the degree of freedom indices for the
324  // current element. These define where in the global
325  // matrix and right-hand-side this element will
326  // contribute to.
327  dof_map.dof_indices (elem, dof_indices);
328 
329  // Cache the total number of degrees of freedom on this element,
330  // for use as array and loop bounds later.
331  // We use cast_int to explicitly convert from size() (which may be
332  // 64-bit) to unsigned int (which may be 32-bit but which is definitely
333  // enough to count *local* degrees of freedom.
334  const unsigned int n_dofs =
335  cast_int<unsigned int>(dof_indices.size());
336 
337  // Compute the element-specific data for the current
338  // element. This involves computing the location of the
339  // quadrature points (q_point) and the shape functions
340  // and their divergences for the current element.
341  vector_fe->reinit (elem);
342 
343  // We should also have the same number of degrees of freedom as
344  // shape functions for our variable.
345  libmesh_assert_equal_to (n_dofs, vector_phi.size());
346 
347  // Zero the element matrix and right-hand side before
348  // summing them. We use the resize member here because
349  // the number of degrees of freedom might have changed from
350  // the last element. Note that this will be the case if the
351  // element type is different (i.e. the last element was a
352  // triangle, now we are on a quadrilateral).
353 
354  // The DenseMatrix::resize() and the DenseVector::resize()
355  // members will automatically zero out the matrix and vector.
356  Ke.resize (n_dofs, n_dofs);
357  Fe.resize (n_dofs);
358 
359  // Now loop over the quadrature points. This handles
360  // the numeric integration.
361  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
362  {
363 
364  // Now we will build the element matrix.
365  // This a double loop to integrate the vector test functions (i)
366  // against the vector trial functions (j) and their divergences
367  for (unsigned int i = 0; i != n_dofs; i++)
368  for (unsigned int j = 0; j != n_dofs; j++)
369  {
370  Ke(i, j) += JxW[qp]*(div_vector_phi[i][qp]*div_vector_phi[j][qp]+
371  vector_phi[i][qp]*vector_phi[j][qp]);
372  }
373 
374  // This is the end of the matrix summation loop
375  // Now we build the element right-hand-side contribution.
376  // This involves a single loop in which we integrate the "forcing
377  // function" in the PDE against the vector test functions (k).
378  {
379  // The location of the current quadrature point.
380  const Real x = q_point[qp](0);
381  const Real y = q_point[qp](1);
382  const Real z = q_point[qp](2);
383 
384  // "f" is the forcing function, given by the well-known
385  // "method of manufactured solutions".
387 
388  // Loop to integrate the vector test functions (k) against the
389  // forcing function.
390  for (unsigned int k = 0; k != n_dofs; k++)
391  {
392  Fe(k) += JxW[qp]*f*vector_phi[k][qp];
393  }
394  }
395  }
396 
397  // We have now reached the end of the quadrature point loop, so
398  // the interior element integration has been completed. However, we have
399  // not yet addressed boundary conditions.
400  {
401 
402  // The following loop is over the sides of the element.
403  // If the element has no neighbor on a side then that
404  // side MUST live on a boundary of the domain.
405  for (auto side : elem->side_index_range())
406  if (elem->neighbor_ptr(side) == nullptr)
407  {
408  // The value of the shape functions at the quadrature points.
409  const std::vector<std::vector<RealGradient>> & vector_phi_face = vector_fe_face->get_phi();
410 
411  // The Jacobian * Quadrature Weight at the quadrature
412  // points on the face.
413  const std::vector<Real> & JxW_face = vector_fe_face->get_JxW();
414 
415  // The XYZ locations (in physical space) of, and the normals at,
416  // the quadrature points on the face. This is where
417  // we will interpolate the boundary value function.
418  const std::vector<Point> & qface_point = vector_fe_face->get_xyz();
419  const std::vector<Point> & normals = vector_fe_face->get_normals();
420 
421  // Compute the vector shape function values on the element face.
422  vector_fe_face->reinit(elem, side);
423 
424  // Some shape functions will be 0 on the face, but for ease of
425  // indexing and generality of code we loop over them anyway.
426  libmesh_assert_equal_to (n_dofs, vector_phi_face.size());
427 
428  // Loop over the face quadrature points for integration.
429  for (unsigned int qp=0; qp<qface.n_points(); qp++)
430  {
431  // The location on the boundary of the current
432  // face quadrature point.
433  const Real xf = qface_point[qp](0);
434  const Real yf = qface_point[qp](1);
435  const Real zf = qface_point[qp](2);
436 
437  // The boundary value for the vector variable.
438  RealGradient vector_value = GradDivExactSolution()(xf, yf, zf);
439 
440  // We use the penalty method to set the flux of the vector
441  // variable at the boundary, i.e. the RT vector boundary dof.
442  const Real penalty = 1.e10;
443 
444  // A double loop to integrate the normal component of the
445  // vector test functions (i) against the normal component of
446  // the vector trial functions (j).
447  for (unsigned int i = 0; i != n_dofs; i++)
448  for (unsigned int j = 0; j != n_dofs; j++)
449  {
450  Ke(i, j) += JxW_face[qp]*penalty*vector_phi_face[i][qp]*
451  normals[qp]*vector_phi_face[j][qp]*normals[qp];
452  }
453 
454  // Loop to integrate the normal component of the vector test
455  // functions (i) against the normal component of the
456  // exact solution for the vector variable.
457  for (unsigned int i = 0; i != n_dofs; i++)
458  {
459  Fe(i) += JxW_face[qp]*penalty*vector_phi_face[i][qp]*normals[qp]*
460  vector_value*normals[qp];
461  }
462  }
463  }
464  }
465 
466  // We have now finished the quadrature point loop,
467  // and have therefore applied all the boundary conditions.
468 
469  // If this assembly program were to be used on an adaptive mesh,
470  // we would have to apply any hanging node constraint equations.
471  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
472 
473  // The element matrix and right-hand-side are now built
474  // for this element. Add them to the global matrix and
475  // right-hand-side vector. The SparseMatrix::add_matrix()
476  // and NumericVector::add_vector() members do this for us.
477  matrix.add_matrix (Ke, dof_indices);
478  system.rhs->add_vector (Fe, dof_indices);
479  }
480 
481  // All done!
482 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
This is the EquationSystems class.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
Real hdiv_error(std::string_view sys_name, std::string_view unknown_name)
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void assemble_graddiv(EquationSystems &es, const std::string &system_name)
void attach_exact_values(const std::vector< FunctionBase< Number > *> &f)
Clone and attach arbitrary functors which compute the exact values of the EquationSystems&#39; solutions ...
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
void extra_quadrature_order(const int extraorder)
Increases or decreases the order of the quadrature rule used for numerical integration.
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
Definition: fe_type.h:371
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.
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.
void attach_exact_derivs(const std::vector< FunctionBase< Gradient > *> &g)
Clone and attach arbitrary functors which compute the exact gradients of the EquationSystems&#39; solutio...
virtual void solve() override
Assembles & solves the linear system A*x=b.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
unsigned int variable_number(std::string_view var) const
Definition: system.C:1587
void permute_elements(MeshBase &mesh)
Randomly permute the nodal ordering of each element (without twisting the element mapping)...
SolverPackage default_solver_package()
Definition: libmesh.C:1111
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:178
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.
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:2013
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:1544
void compute_error(std::string_view sys_name, std::string_view unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
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.
RealGradient forcing(Real x, Real y, Real z)
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:1335
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:2139
Real l2_error(std::string_view sys_name, std::string_view unknown_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
const MeshBase & get_mesh() const
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:354
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
int main(int argc, char **argv)
const DofMap & get_dof_map() const
Definition: system.h:2370
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.
Real error_norm(std::string_view sys_name, std::string_view unknown_name, const FEMNormType &norm)