libMesh
introduction_ex5.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 #include "libmesh/enum_quadrature_type.h"
70 
71 // Bring in everything from the libMesh namespace
72 using namespace libMesh;
73 
74 
75 
76 // Function prototype, as before.
78  const std::string & system_name);
79 
80 // Exact solution function prototype, as before.
81 Real exact_solution (const Real x,
82  const Real y,
83  const Real z = 0.);
84 
85 // Define a wrapper for exact_solution that will be needed below
87  const Point & p,
88  const Real)
89 {
90  output(0) = exact_solution(p(0), p(1), p(2));
91 }
92 
93 
94 // The quadrature type the user requests.
96 
97 
98 
99 // Begin the main program.
100 int main (int argc, char ** argv)
101 {
102  // Initialize libMesh and any dependent libraries, like in example 2.
103  LibMeshInit init (argc, argv);
104 
105  // This example requires a linear solver package.
106  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
107  "--enable-petsc, --enable-trilinos, or --enable-eigen");
108 
109  // Check for proper usage. The quadrature rule
110  // must be given at run time.
111  if (argc < 3)
112  {
113  libmesh_error_msg("Usage: " << argv[0] << " -q <rule>\n" \
114  << " where <rule> is one of QGAUSS, QSIMPSON, or QTRAP.");
115  }
116 
117 
118  // Tell the user what we are doing.
119  else
120  {
121  libMesh::out << "Running " << argv[0];
122 
123  for (int i=1; i<argc; i++)
124  libMesh::out << " " << argv[i];
125 
126  libMesh::out << std::endl << std::endl;
127  }
128 
129 
130  // Set the quadrature rule type that the user wants from argv[2]
131  quad_type = static_cast<QuadratureType>(std::atoi(argv[2]));
132 
133  // Skip this 3D example if libMesh was compiled as 1D-only.
134  libmesh_example_requires(3 <= LIBMESH_DIM, "3D support");
135 
136  // We use Dirichlet boundary conditions here
137 #ifndef LIBMESH_ENABLE_DIRICHLET
138  libmesh_example_requires(false, "--enable-dirichlet");
139 #endif
140 
141  // The following is identical to example 4, and therefore
142  // not commented. Differences are mentioned when present.
143  Mesh mesh(init.comm());
144 
145  // We will use a linear approximation space in this example,
146  // hence 8-noded hexahedral elements are sufficient. This
147  // is different than example 4 where we used 27-noded
148  // hexahedral elements to support a second-order approximation
149  // space.
151  16, 16, 16,
152  -1., 1.,
153  -1., 1.,
154  -1., 1.,
155  HEX8);
156 
157  mesh.print_info();
158 
159  EquationSystems equation_systems (mesh);
160 
161  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
162 
163  unsigned int u_var = equation_systems.get_system("Poisson").add_variable("u", FIRST);
164 
165  equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
166 
167  // Construct a Dirichlet boundary condition object
168 
169  // Indicate which boundary IDs we impose the BC on
170  // We either build a line, a square or a cube, and
171  // here we indicate the boundaries IDs in each case
172  std::set<boundary_id_type> boundary_ids;
173  // the dim==1 mesh has two boundaries with IDs 0 and 1
174  boundary_ids.insert(0);
175  boundary_ids.insert(1);
176  boundary_ids.insert(2);
177  boundary_ids.insert(3);
178  boundary_ids.insert(4);
179  boundary_ids.insert(5);
180 
181  // Create a vector storing the variable numbers which the BC applies to
182  std::vector<unsigned int> variables(1);
183  variables[0] = u_var;
184 
185  // Create an AnalyticFunction object that we use to project the BC
186  // This function just calls the function exact_solution via exact_solution_wrapper
187  AnalyticFunction<> exact_solution_object(exact_solution_wrapper);
188 
189 #ifdef LIBMESH_ENABLE_DIRICHLET
190  // In general, when reusing a system-indexed exact solution, we want
191  // to use the default system-ordering constructor for
192  // DirichletBoundary, so we demonstrate that here. In this case,
193  // though, we have only one variable, so system- and local-
194  // orderings are the same.
195  DirichletBoundary dirichlet_bc
196  (boundary_ids, variables, exact_solution_object);
197 
198  // We must add the Dirichlet boundary condition _before_
199  // we call equation_systems.init()
200  equation_systems.get_system("Poisson").get_dof_map().add_dirichlet_boundary(dirichlet_bc);
201 #endif
202 
203  equation_systems.init();
204 
205  equation_systems.print_info();
206 
207  equation_systems.get_system("Poisson").solve();
208 
209  // "Personalize" the output, with the
210  // number of the quadrature rule appended.
211  std::ostringstream f_name;
212  f_name << "out_" << quad_type << ".e";
213 
214 #ifdef LIBMESH_HAVE_EXODUS_API
215  ExodusII_IO(mesh).write_equation_systems (f_name.str(),
216  equation_systems);
217 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
218 
219  // All done.
220  return 0;
221 }
222 
223 
224 
225 
227  const std::string & libmesh_dbg_var(system_name))
228 {
229  libmesh_assert_equal_to (system_name, "Poisson");
230 
231  const MeshBase & mesh = es.get_mesh();
232 
233  const unsigned int dim = mesh.mesh_dimension();
234 
235  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
236 
237  const DofMap & dof_map = system.get_dof_map();
238 
239  FEType fe_type = dof_map.variable_type(0);
240 
241  // Build a Finite Element object of the specified type. Since the
242  // FEBase::build() member dynamically creates memory we will
243  // store the object as a std::unique_ptr<FEBase>. Below, the
244  // functionality of std::unique_ptr's is described more detailed in
245  // the context of building quadrature rules.
246  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
247 
248  // Now this deviates from example 4. we create a
249  // 5th order quadrature rule of user-specified type
250  // for numerical integration. Note that not all
251  // quadrature rules support this order.
252  std::unique_ptr<QBase> qrule(QBase::build(quad_type, dim, THIRD));
253 
254  // Tell the finite element object to use our
255  // quadrature rule. Note that a std::unique_ptr<QBase> returns
256  // a QBase* pointer to the object it handles with get().
257  // However, using get(), the std::unique_ptr<QBase> qrule is
258  // still in charge of this pointer. I.e., when qrule goes
259  // out of scope, it will safely delete the QBase object it
260  // points to. This behavior may be overridden using
261  // std::unique_ptr<Xyz>::release(), but is currently not
262  // recommended.
263  fe->attach_quadrature_rule (qrule.get());
264 
265  // Declare a special finite element object for
266  // boundary integration.
267  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
268 
269  // As already seen in example 3, boundary integration
270  // requires a quadrature rule. Here, however,
271  // we use the more convenient way of building this
272  // rule at run-time using quad_type. Note that one
273  // could also have initialized the face quadrature rules
274  // with the type directly determined from qrule, namely
275  // through:
276  // \verbatim
277  // std::unique_ptr<QBase> qface (QBase::build(qrule->type(),
278  // dim-1,
279  // THIRD));
280  // \endverbatim
281  // And again: using the std::unique_ptr<QBase> relaxes
282  // the need to delete the object afterward,
283  // they clean up themselves.
284  std::unique_ptr<QBase> qface (QBase::build(quad_type,
285  dim-1,
286  THIRD));
287 
288  // Tell the finite element object to use our
289  // quadrature rule. Note that a std::unique_ptr<QBase> returns
290  // a QBase* pointer to the object it handles with get().
291  // However, using get(), the std::unique_ptr<QBase> qface is
292  // still in charge of this pointer. I.e., when qface goes
293  // out of scope, it will safely delete the QBase object it
294  // points to. This behavior may be overridden using
295  // std::unique_ptr<Xyz>::release(), but is not recommended.
296  fe_face->attach_quadrature_rule (qface.get());
297 
298  // This is again identical to example 4, and not commented.
299  const std::vector<Real> & JxW = fe->get_JxW();
300 
301  const std::vector<Point> & q_point = fe->get_xyz();
302 
303  const std::vector<std::vector<Real>> & phi = fe->get_phi();
304 
305  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
306 
309  std::vector<dof_id_type> dof_indices;
310 
311  // Now we will loop over all the elements in the mesh.
312  // See example 3 for details.
313  for (const auto & elem : mesh.active_local_element_ptr_range())
314  {
315  dof_map.dof_indices (elem, dof_indices);
316 
317  const unsigned int n_dofs =
318  cast_int<unsigned int>(dof_indices.size());
319 
320  fe->reinit (elem);
321 
322  libmesh_assert_equal_to (n_dofs, phi.size());
323 
324  Ke.resize (n_dofs, n_dofs);
325 
326  Fe.resize (n_dofs);
327 
328  // Now loop over the quadrature points. This handles
329  // the numeric integration. Note the slightly different
330  // access to the QBase members!
331  for (unsigned int qp=0; qp<qrule->n_points(); qp++)
332  {
333  // Add the matrix contribution
334  for (unsigned int i=0; i != n_dofs; i++)
335  for (unsigned int j=0; j != n_dofs; j++)
336  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
337 
338  // fxy is the forcing function for the Poisson equation.
339  // In this case we set fxy to be a finite difference
340  // Laplacian approximation to the (known) exact solution.
341  //
342  // We will use the second-order accurate FD Laplacian
343  // approximation, which in 2D on a structured grid is
344  //
345  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
346  // u(i,j-1) + u(i,j+1) +
347  // -4*u(i,j))/h^2
348  //
349  // Since the value of the forcing function depends only
350  // on the location of the quadrature point (q_point[qp])
351  // we will compute it here, outside of the i-loop
352  const Real x = q_point[qp](0);
353  const Real y = q_point[qp](1);
354  const Real z = q_point[qp](2);
355  const Real eps = 1.e-3;
356 
357  const Real uxx = (exact_solution(x-eps, y, z) +
358  exact_solution(x+eps, y, z) +
359  -2.*exact_solution(x, y, z))/eps/eps;
360 
361  const Real uyy = (exact_solution(x, y-eps, z) +
362  exact_solution(x, y+eps, z) +
363  -2.*exact_solution(x, y, z))/eps/eps;
364 
365  const Real uzz = (exact_solution(x, y, z-eps) +
366  exact_solution(x, y, z+eps) +
367  -2.*exact_solution(x, y, z))/eps/eps;
368 
369  const Real fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
370 
371 
372  // Add the RHS contribution
373  for (unsigned int i=0; i != n_dofs; i++)
374  Fe(i) += JxW[qp]*fxy*phi[i][qp];
375  }
376 
377  // If this assembly program were to be used on an adaptive mesh,
378  // we would have to apply any hanging node constraint equations
379  // Call heterogenously_constrain_element_matrix_and_vector to impose
380  // non-homogeneous Dirichlet BCs
381  dof_map.heterogenously_constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
382 
383  // The element matrix and right-hand-side are now built
384  // for this element. Add them to the global matrix and
385  // right-hand-side vector. The SparseMatrix::add_matrix()
386  // and NumericVector::add_vector() members do this for us.
387  system.matrix->add_matrix (Ke, dof_indices);
388  system.rhs->add_vector (Fe, dof_indices);
389 
390  } // end of element loop
391 }
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
exact_solution_wrapper
void exact_solution_wrapper(DenseVector< Number > &output, const Point &p, const Real)
Definition: introduction_ex5.C:86
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::EquationSystems::get_mesh
const MeshBase & get_mesh() const
Definition: equation_systems.h:637
libMesh::ExplicitSystem::rhs
NumericVector< Number > * rhs
The system matrix.
Definition: explicit_system.h:114
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
assemble_poisson
void assemble_poisson(EquationSystems &es, const std::string &system_name)
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::QuadratureType
QuadratureType
Defines an enum for currently available quadrature rules.
Definition: enum_quadrature_type.h:33
libMesh::MeshTools::Generation::build_cube
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.
Definition: mesh_generation.C:298
libMesh::DenseMatrix< Number >
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::QBase::build
static std::unique_ptr< QBase > build(const std::string &name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
Definition: quadrature_build.C:41
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::NumericVector::add_vector
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].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::INVALID_Q_RULE
Definition: enum_quadrature_type.h:48
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::FEGenericBase::build
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libMesh::AnalyticFunction
Wraps a function pointer into a FunctionBase object.
Definition: analytic_function.h:48
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
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
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_matrix
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.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
quad_type
QuadratureType quad_type
Definition: introduction_ex5.C:95
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::DirichletBoundary
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Definition: dirichlet_boundaries.h:88
main
int main(int argc, char **argv)
Definition: introduction_ex5.C:100
libMesh::THIRD
Definition: enum_order.h:44
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::LinearImplicitSystem
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
Definition: linear_implicit_system.h:55
libMesh::out
OStreamProxy out
libMesh::FIRST
Definition: enum_order.h:42
exact_solution
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
Definition: exact_solution.C:43
libMesh::DenseVector< Number >