libMesh
introduction_ex5.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>Introduction Example 5 - Run-Time Quadrature Rule Selection</h1>
21 // \author Benjamin S. Kirk
22 // \date 2003
23 //
24 // This is the fifth example program. It builds on
25 // the previous two examples, and extends the use
26 // of the std::unique_ptr as a convenient build method to
27 // determine the quadrature rule at run time.
28 
29 
30 // C++ include files that we need
31 #include <iostream>
32 #include <sstream>
33 #include <algorithm>
34 #include <math.h>
35 
36 // Basic include file needed for the mesh functionality.
37 #include "libmesh/libmesh.h"
38 #include "libmesh/mesh.h"
39 #include "libmesh/mesh_generation.h"
40 #include "libmesh/exodusII_io.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 the base quadrature class, with which
48 // specialized quadrature rules will be built.
49 #include "libmesh/quadrature.h"
50 
51 // Define useful datatypes for finite element
52 // matrix and vector components.
53 #include "libmesh/sparse_matrix.h"
54 #include "libmesh/numeric_vector.h"
55 #include "libmesh/dense_matrix.h"
56 #include "libmesh/dense_vector.h"
57 
58 // Define the DofMap, which handles degree of freedom
59 // indexing.
60 #include "libmesh/dof_map.h"
61 
62 // To impose Dirichlet boundary conditions
63 #include "libmesh/dirichlet_boundaries.h"
64 #include "libmesh/analytic_function.h"
65 
66 // The definition of a geometric element
67 #include "libmesh/elem.h"
68 #include "libmesh/enum_solver_package.h"
69 
70 // Reading a quadrature rule from user arguments
71 #include "libmesh/enum_quadrature_type.h"
72 #include "libmesh/string_to_enum.h"
73 
74 // Bring in everything from the libMesh namespace
75 using namespace libMesh;
76 
77 
78 
79 // Function prototype, as before.
81  const std::string & system_name);
82 
83 // Exact solution function prototype, as before.
84 Real exact_solution (const Real x,
85  const Real y,
86  const Real z = 0.);
87 
88 // Define a wrapper for exact_solution that will be needed below
90  const Point & p,
91  const Real)
92 {
93  output(0) = exact_solution(p(0), p(1), p(2));
94 }
95 
96 
97 // The quadrature type the user requests.
99 
100 
101 
102 // Begin the main program.
103 int main (int argc, char ** argv)
104 {
105  // Initialize libMesh and any dependent libraries, like in example 2.
106  LibMeshInit init (argc, argv);
107 
108  // This example requires a linear solver package.
109  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
110  "--enable-petsc, --enable-trilinos, or --enable-eigen");
111 
112  // Check for proper usage. The quadrature rule
113  // must be given at run time.
114  libmesh_error_msg_if(argc < 3,
115  "Usage: " << argv[0] << " -q <rule>\n"
116  " where <rule> is one of QGAUSS, QSIMPSON, or QTRAP.");
117 
118  // Tell the user what we are doing.
119  libMesh::out << "Running " << argv[0];
120 
121  for (int i=1; i<argc; i++)
122  libMesh::out << " " << argv[i];
123 
124  libMesh::out << std::endl << std::endl;
125 
126  // Set the quadrature rule type that the user wants
127  quad_type = Utility::string_to_enum<QuadratureType>
128  (libMesh::command_line_next("-q", std::string("QGAUSS")));
129 
130  // Skip this 3D example if libMesh was compiled as 1D-only.
131  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
132 
133  // We use Dirichlet boundary conditions here
134 #ifndef LIBMESH_ENABLE_DIRICHLET
135  libmesh_example_requires(false, "--enable-dirichlet");
136 #endif
137 
138  // The following is identical to example 4, and therefore
139  // not commented. Differences are mentioned when present.
140  Mesh mesh(init.comm());
141 
142  // We will use a linear approximation space in this example,
143  // hence 8-noded hexahedral elements are sufficient. This
144  // is different than example 4 where we used 27-noded
145  // hexahedral elements to support a second-order approximation
146  // space.
148  16, 16, 16,
149  -1., 1.,
150  -1., 1.,
151  -1., 1.,
152  HEX8);
153 
154  mesh.print_info();
155 
156  EquationSystems equation_systems (mesh);
157 
158  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
159 
160  unsigned int u_var = equation_systems.get_system("Poisson").add_variable("u", FIRST);
161 
162  equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
163 
164  // Construct a Dirichlet boundary condition object
165 
166  // Indicate which boundary IDs we impose the BC on
167  // We either build a line, a square or a cube, and
168  // here we indicate boundaries covering each case
169  std::set<boundary_id_type> boundary_ids {0,1,2,3,4,5};
170 
171  // Create an AnalyticFunction object that we use to project the BC
172  // This function just calls the function exact_solution via exact_solution_wrapper
173  AnalyticFunction<> exact_solution_object(exact_solution_wrapper);
174 
175 #ifdef LIBMESH_ENABLE_DIRICHLET
176  // In general, when reusing a system-indexed exact solution, we want
177  // to use the default system-ordering constructor for
178  // DirichletBoundary, so we demonstrate that here. In this case,
179  // though, we have only one variable, so system- and local-
180  // orderings are the same.
181  DirichletBoundary dirichlet_bc
182  (boundary_ids, {u_var}, exact_solution_object);
183 
184  // We must add the Dirichlet boundary condition _before_
185  // we call equation_systems.init()
186  equation_systems.get_system("Poisson").get_dof_map().add_dirichlet_boundary(dirichlet_bc);
187 #endif
188 
189  equation_systems.init();
190 
191  equation_systems.print_info();
192 
193  equation_systems.get_system("Poisson").solve();
194 
195  // "Personalize" the output, with the
196  // number of the quadrature rule appended.
197  std::ostringstream f_name;
198  f_name << "out_" << quad_type << ".e";
199 
200 #ifdef LIBMESH_HAVE_EXODUS_API
201  ExodusII_IO(mesh).write_equation_systems (f_name.str(),
202  equation_systems);
203 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
204 
205  // All done.
206  return 0;
207 }
208 
209 
210 
211 
213  const std::string & libmesh_dbg_var(system_name))
214 {
215  libmesh_assert_equal_to (system_name, "Poisson");
216 
217  const MeshBase & mesh = es.get_mesh();
218 
219  const unsigned int dim = mesh.mesh_dimension();
220 
221  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
222 
223  const DofMap & dof_map = system.get_dof_map();
224 
225  FEType fe_type = dof_map.variable_type(0);
226 
227  // Build a Finite Element object of the specified type. Since the
228  // FEBase::build() member dynamically creates memory we will
229  // store the object as a std::unique_ptr<FEBase>. Below, the
230  // functionality of std::unique_ptr's is described more detailed in
231  // the context of building quadrature rules.
232  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
233 
234  // Now this deviates from example 4. we create a
235  // 5th order quadrature rule of user-specified type
236  // for numerical integration. Note that not all
237  // quadrature rules support this order.
238  std::unique_ptr<QBase> qrule(QBase::build(quad_type, dim, THIRD));
239 
240  // Tell the finite element object to use our
241  // quadrature rule. Note that a std::unique_ptr<QBase> returns
242  // a QBase* pointer to the object it handles with get().
243  // However, using get(), the std::unique_ptr<QBase> qrule is
244  // still in charge of this pointer. I.e., when qrule goes
245  // out of scope, it will safely delete the QBase object it
246  // points to. This behavior may be overridden using
247  // std::unique_ptr<Xyz>::release(), but is currently not
248  // recommended.
249  fe->attach_quadrature_rule (qrule.get());
250 
251  // Declare a special finite element object for
252  // boundary integration.
253  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
254 
255  // As already seen in example 3, boundary integration
256  // requires a quadrature rule. Here, however,
257  // we use the more convenient way of building this
258  // rule at run-time using quad_type. Note that one
259  // could also have initialized the face quadrature rules
260  // with the type directly determined from qrule, namely
261  // through:
262  // \verbatim
263  // std::unique_ptr<QBase> qface (QBase::build(qrule->type(),
264  // dim-1,
265  // THIRD));
266  // \endverbatim
267  // And again: using the std::unique_ptr<QBase> relaxes
268  // the need to delete the object afterward,
269  // they clean up themselves.
270  std::unique_ptr<QBase> qface (QBase::build(quad_type,
271  dim-1,
272  THIRD));
273 
274  // Tell the finite element object to use our
275  // quadrature rule. Note that a std::unique_ptr<QBase> returns
276  // a QBase* pointer to the object it handles with get().
277  // However, using get(), the std::unique_ptr<QBase> qface is
278  // still in charge of this pointer. I.e., when qface goes
279  // out of scope, it will safely delete the QBase object it
280  // points to. This behavior may be overridden using
281  // std::unique_ptr<Xyz>::release(), but is not recommended.
282  fe_face->attach_quadrature_rule (qface.get());
283 
284  // This is again identical to example 4, and not commented.
285  const std::vector<Real> & JxW = fe->get_JxW();
286 
287  const std::vector<Point> & q_point = fe->get_xyz();
288 
289  const std::vector<std::vector<Real>> & phi = fe->get_phi();
290 
291  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
292 
295  std::vector<dof_id_type> dof_indices;
296 
297  // The global system matrix
298  SparseMatrix<Number> & matrix = system.get_system_matrix();
299 
300  // Now we will loop over all the elements in the mesh.
301  // See example 3 for details.
302  for (const auto & elem : mesh.active_local_element_ptr_range())
303  {
304  dof_map.dof_indices (elem, dof_indices);
305 
306  const unsigned int n_dofs =
307  cast_int<unsigned int>(dof_indices.size());
308 
309  fe->reinit (elem);
310 
311  libmesh_assert_equal_to (n_dofs, phi.size());
312 
313  Ke.resize (n_dofs, n_dofs);
314 
315  Fe.resize (n_dofs);
316 
317  // Now loop over the quadrature points. This handles
318  // the numeric integration. Note the slightly different
319  // access to the QBase members!
320  for (unsigned int qp=0; qp<qrule->n_points(); qp++)
321  {
322  // Add the matrix contribution
323  for (unsigned int i=0; i != n_dofs; i++)
324  for (unsigned int j=0; j != n_dofs; j++)
325  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
326 
327  // fxy is the forcing function for the Poisson equation.
328  // In this case we set fxy to be a finite difference
329  // Laplacian approximation to the (known) exact solution.
330  //
331  // We will use the second-order accurate FD Laplacian
332  // approximation, which in 2D on a structured grid is
333  //
334  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
335  // u(i,j-1) + u(i,j+1) +
336  // -4*u(i,j))/h^2
337  //
338  // Since the value of the forcing function depends only
339  // on the location of the quadrature point (q_point[qp])
340  // we will compute it here, outside of the i-loop
341  const Real x = q_point[qp](0);
342  const Real y = q_point[qp](1);
343  const Real z = q_point[qp](2);
344  const Real eps = 1.e-3;
345 
346  const Real uxx = (exact_solution(x-eps, y, z) +
347  exact_solution(x+eps, y, z) +
348  -2.*exact_solution(x, y, z))/eps/eps;
349 
350  const Real uyy = (exact_solution(x, y-eps, z) +
351  exact_solution(x, y+eps, z) +
352  -2.*exact_solution(x, y, z))/eps/eps;
353 
354  const Real uzz = (exact_solution(x, y, z-eps) +
355  exact_solution(x, y, z+eps) +
356  -2.*exact_solution(x, y, z))/eps/eps;
357 
358  const Real fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
359 
360 
361  // Add the RHS contribution
362  for (unsigned int i=0; i != n_dofs; i++)
363  Fe(i) += JxW[qp]*fxy*phi[i][qp];
364  }
365 
366  // If this assembly program were to be used on an adaptive mesh,
367  // we would have to apply any hanging node constraint equations
368  // Call heterogenously_constrain_element_matrix_and_vector to impose
369  // non-homogeneous Dirichlet BCs
370  dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
371 
372  // The element matrix and right-hand-side are now built
373  // for this element. Add them to the global matrix and
374  // right-hand-side vector. The SparseMatrix::add_matrix()
375  // and NumericVector::add_vector() members do this for us.
376  matrix.add_matrix (Ke, dof_indices);
377  system.rhs->add_vector (Fe, dof_indices);
378 
379  } // end of element loop
380 }
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.
unsigned int dim
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 ...
QuadratureType quad_type
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
NumericVector< Number > * rhs
The system matrix.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
QuadratureType
Defines an enum for currently available quadrature rules.
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
Wraps a function pointer into a FunctionBase object.
This is the MeshBase class.
Definition: mesh_base.h:75
SolverPackage default_solver_package()
Definition: libmesh.C:1117
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
void exact_solution_wrapper(DenseVector< Number > &output, const Point &p, const Real)
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.
int main(int argc, char **argv)
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
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.
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
unsigned int mesh_dimension() const
Definition: mesh_base.C:372
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
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
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
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.
void assemble_poisson(EquationSystems &es, const std::string &system_name)