libMesh
Functions
introduction_ex3.C File Reference

Go to the source code of this file.

Functions

void assemble_poisson (EquationSystems &es, const std::string &system_name)
 
Real exact_solution (const Real x, const Real y, const Real z=0.)
 This is the exact solution that we are trying to obtain. More...
 
int main (int argc, char **argv)
 
void assemble_poisson (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_poisson() [1/2]

void assemble_poisson ( EquationSystems es,
const std::string &  system_name 
)

Referenced by main().

◆ assemble_poisson() [2/2]

void assemble_poisson ( EquationSystems es,
const std::string &  libmesh_dbg_varsystem_name 
)

Definition at line 179 of file introduction_ex3.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, exact_solution(), libMesh::FIFTH, libMesh::System::get_dof_map(), libMesh::EquationSystems::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::ImplicitSystem::get_system_matrix(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::QBase::n_points(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and value.

181 {
182 
183  // It is a good idea to make sure we are assembling
184  // the proper system.
185  libmesh_assert_equal_to (system_name, "Poisson");
186 
187  // Get a constant reference to the mesh object.
188  const MeshBase & mesh = es.get_mesh();
189 
190  // The dimension that we are running
191  const unsigned int dim = mesh.mesh_dimension();
192 
193  // Get a reference to the LinearImplicitSystem we are solving
194  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Poisson");
195 
196  // A reference to the DofMap object for this system. The DofMap
197  // object handles the index translation from node and element numbers
198  // to degree of freedom numbers. We will talk more about the DofMap
199  // in future examples.
200  const DofMap & dof_map = system.get_dof_map();
201 
202  // Get a constant reference to the Finite Element type
203  // for the first (and only) variable in the system.
204  FEType fe_type = dof_map.variable_type(0);
205 
206  // Build a Finite Element object of the specified type. Since the
207  // FEBase::build() member dynamically creates memory we will
208  // store the object as a std::unique_ptr<FEBase>. This can be thought
209  // of as a pointer that will clean up after itself. Introduction Example 4
210  // describes some advantages of std::unique_ptr's in the context of
211  // quadrature rules.
212  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
213 
214  // A 5th order Gauss quadrature rule for numerical integration.
215  QGauss qrule (dim, FIFTH);
216 
217  // Tell the finite element object to use our quadrature rule.
218  fe->attach_quadrature_rule (&qrule);
219 
220  // Declare a special finite element object for
221  // boundary integration.
222  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
223 
224  // Boundary integration requires one quadrature rule,
225  // with dimensionality one less than the dimensionality
226  // of the element.
227  QGauss qface(dim-1, FIFTH);
228 
229  // Tell the finite element object to use our
230  // quadrature rule.
231  fe_face->attach_quadrature_rule (&qface);
232 
233  // Here we define some references to cell-specific data that
234  // will be used to assemble the linear system.
235  //
236  // The element Jacobian * quadrature weight at each integration point.
237  const std::vector<Real> & JxW = fe->get_JxW();
238 
239  // The physical XY locations of the quadrature points on the element.
240  // These might be useful for evaluating spatially varying material
241  // properties at the quadrature points.
242  const std::vector<Point> & q_point = fe->get_xyz();
243 
244  // The element shape functions evaluated at the quadrature points.
245  const std::vector<std::vector<Real>> & phi = fe->get_phi();
246 
247  // The element shape function gradients evaluated at the quadrature
248  // points.
249  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
250 
251  // Define data structures to contain the element matrix
252  // and right-hand-side vector contribution. Following
253  // basic finite element terminology we will denote these
254  // "Ke" and "Fe". These datatypes are templated on
255  // Number, which allows the same code to work for real
256  // or complex numbers.
259 
260  // This vector will hold the degree of freedom indices for
261  // the element. These define where in the global system
262  // the element degrees of freedom get mapped.
263  std::vector<dof_id_type> dof_indices;
264 
265  // The global system matrix
266  SparseMatrix<Number> & matrix = system.get_system_matrix();
267 
268  // Now we will loop over all the elements in the mesh.
269  // We will compute the element matrix and right-hand-side
270  // contribution.
271  //
272  // Element ranges are a nice way to iterate through all the
273  // elements, or all the elements that have some property. The
274  // range will iterate from the first to the last element on
275  // the local processor.
276  // It is smart to make this one const so that we don't accidentally
277  // mess it up! In case users later modify this program to include
278  // refinement, we will be safe and will only consider the active
279  // elements; hence we use a variant of the
280  // active_local_element_ptr_range.
281  for (const auto & elem : mesh.active_local_element_ptr_range())
282  {
283  // Get the degree of freedom indices for the
284  // current element. These define where in the global
285  // matrix and right-hand-side this element will
286  // contribute to.
287  dof_map.dof_indices (elem, dof_indices);
288 
289  // Cache the number of degrees of freedom on this element, for
290  // use as a loop bound later. We use cast_int to explicitly
291  // convert from size() (which may be 64-bit) to unsigned int
292  // (which may be 32-bit but which is definitely enough to count
293  // *local* degrees of freedom.
294  const unsigned int n_dofs =
295  cast_int<unsigned int>(dof_indices.size());
296 
297  // Compute the element-specific data for the current
298  // element. This involves computing the location of the
299  // quadrature points (q_point) and the shape functions
300  // (phi, dphi) for the current element.
301  fe->reinit (elem);
302 
303  // With one variable, we should have the same number of degrees
304  // of freedom as shape functions.
305  libmesh_assert_equal_to (n_dofs, phi.size());
306 
307  // Zero the element matrix and right-hand side before
308  // summing them. We use the resize member here because
309  // the number of degrees of freedom might have changed from
310  // the last element. Note that this will be the case if the
311  // element type is different (i.e. the last element was a
312  // triangle, now we are on a quadrilateral).
313 
314  // The DenseMatrix::resize() and the DenseVector::resize()
315  // members will automatically zero out the matrix and vector.
316  Ke.resize (n_dofs, n_dofs);
317 
318  Fe.resize (n_dofs);
319 
320  // Now loop over the quadrature points. This handles
321  // the numeric integration.
322  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
323  {
324 
325  // Now we will build the element matrix. This involves
326  // a double loop to integrate the test functions (i) against
327  // the trial functions (j).
328  for (unsigned int i=0; i != n_dofs; i++)
329  for (unsigned int j=0; j != n_dofs; j++)
330  {
331  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
332  }
333 
334  // This is the end of the matrix summation loop
335  // Now we build the element right-hand-side contribution.
336  // This involves a single loop in which we integrate the
337  // "forcing function" in the PDE against the test functions.
338  {
339  const Real x = q_point[qp](0);
340  const Real y = q_point[qp](1);
341  const Real eps = 1.e-3;
342 
343 
344  // "fxy" is the forcing function for the Poisson equation.
345  // In this case we set fxy to be a finite difference
346  // Laplacian approximation to the (known) exact solution.
347  //
348  // We will use the second-order accurate FD Laplacian
349  // approximation, which in 2D is
350  //
351  // u_xx + u_yy = (u(i,j-1) + u(i,j+1) +
352  // u(i-1,j) + u(i+1,j) +
353  // -4*u(i,j))/h^2
354  //
355  // Since the value of the forcing function depends only
356  // on the location of the quadrature point (q_point[qp])
357  // we will compute it here, outside of the i-loop
358  const Real fxy = -(exact_solution(x, y-eps) +
359  exact_solution(x, y+eps) +
360  exact_solution(x-eps, y) +
361  exact_solution(x+eps, y) -
362  4.*exact_solution(x, y))/eps/eps;
363 
364  for (unsigned int i=0; i != n_dofs; i++)
365  Fe(i) += JxW[qp]*fxy*phi[i][qp];
366  }
367  }
368 
369  // We have now reached the end of the RHS summation,
370  // and the end of quadrature point loop, so
371  // the interior element integration has
372  // been completed. However, we have not yet addressed
373  // boundary conditions. For this example we will only
374  // consider simple Dirichlet boundary conditions.
375  //
376  // There are several ways Dirichlet boundary conditions
377  // can be imposed. A simple approach, which works for
378  // interpolary bases like the standard Lagrange polynomials,
379  // is to assign function values to the
380  // degrees of freedom living on the domain boundary. This
381  // works well for interpolary bases, but is more difficult
382  // when non-interpolary (e.g Legendre or Hierarchic) bases
383  // are used.
384  //
385  // Dirichlet boundary conditions can also be imposed with a
386  // "penalty" method. In this case essentially the L2 projection
387  // of the boundary values are added to the matrix. The
388  // projection is multiplied by some large factor so that, in
389  // floating point arithmetic, the existing (smaller) entries
390  // in the matrix and right-hand-side are effectively ignored.
391  //
392  // This amounts to adding a term of the form (in latex notation)
393  //
394  // \frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
395  //
396  // where
397  //
398  // \frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
399  {
400 
401  // The following loop is over the sides of the element.
402  // If the element has no neighbor on a side then that
403  // side MUST live on a boundary of the domain.
404  for (auto side : elem->side_index_range())
405  if (elem->neighbor_ptr(side) == nullptr)
406  {
407  // The value of the shape functions at the quadrature
408  // points.
409  const std::vector<std::vector<Real>> & phi_face = 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 = fe_face->get_JxW();
414 
415  // The XYZ locations (in physical space) of the
416  // quadrature points on the face. This is where
417  // we will interpolate the boundary value function.
418  const std::vector<Point> & qface_point = fe_face->get_xyz();
419 
420  // Compute the shape function values on the element
421  // face.
422  fe_face->reinit(elem, side);
423 
424  // Some shape functions will be 0 on the face, but for
425  // ease of indexing and generality of code we loop over
426  // them anyway
427  libmesh_assert_equal_to (n_dofs, phi_face.size());
428 
429  // Loop over the face quadrature points for integration.
430  for (unsigned int qp=0; qp<qface.n_points(); qp++)
431  {
432  // The location on the boundary of the current
433  // face quadrature point.
434  const Real xf = qface_point[qp](0);
435  const Real yf = qface_point[qp](1);
436 
437  // The penalty value. \frac{1}{\epsilon}
438  // in the discussion above.
439  const Real penalty = 1.e10;
440 
441  // The boundary value.
442  const Real value = exact_solution(xf, yf);
443 
444  // Matrix contribution of the L2 projection.
445  for (unsigned int i=0; i != n_dofs; i++)
446  for (unsigned int j=0; j != n_dofs; j++)
447  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
448 
449  // Right-hand-side contribution of the L2
450  // projection.
451  for (unsigned int i=0; i != n_dofs; i++)
452  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
453  }
454  }
455  }
456 
457  // We have now finished the quadrature point loop,
458  // and have therefore applied all the boundary conditions.
459 
460  // If this assembly program were to be used on an adaptive mesh,
461  // we would have to apply any hanging node constraint equations
462  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
463 
464  // The element matrix and right-hand-side are now built
465  // for this element. Add them to the global matrix and
466  // right-hand-side vector. The SparseMatrix::add_matrix()
467  // and NumericVector::add_vector() members do this for us.
468  matrix.add_matrix (Ke, dof_indices);
469  system.rhs->add_vector (Fe, dof_indices);
470  }
471 
472  // All done!
473 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
unsigned int dim
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:374
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]...
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:74
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
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.
Real exact_solution(const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const MeshBase & get_mesh() const
static const bool value
Definition: xdr_io.C:54
void resize(const unsigned int new_m, const unsigned int new_n)
Resizes the matrix to the specified size and calls zero().
Definition: dense_matrix.h:895
This class implements specific orders of Gauss quadrature.
unsigned int mesh_dimension() const
Definition: mesh_base.C:324
const DofMap & get_dof_map() const
Definition: system.h:2293
const SparseMatrix< Number > & get_system_matrix() const

◆ exact_solution()

Real exact_solution ( const Real  x,
const Real  y,
const Real  t 
)

This is the exact solution that we are trying to obtain.

We will solve

  • (u_xx + u_yy) = f

and take a finite difference approximation using this function to get f. This is the well-known "method of manufactured solutions".

Definition at line 43 of file exact_solution.C.

Referenced by assemble_poisson().

46 {
47  static const Real pi = acos(-1.);
48 
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real pi
.
Definition: libmesh.h:274

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 80 of file introduction_ex3.C.

References libMesh::EquationSystems::add_system(), assemble_poisson(), libMesh::MeshTools::Generation::build_square(), libMesh::default_solver_package(), libMesh::EquationSystems::get_system(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD9, libMesh::SECOND, and libMesh::MeshOutput< MT >::write_equation_systems().

81 {
82  // Initialize libraries, like in example 2.
83  LibMeshInit init (argc, argv);
84 
85  // This example requires a linear solver package.
86  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
87  "--enable-petsc, --enable-trilinos, or --enable-eigen");
88 
89  // Brief message to the user regarding the program name
90  // and command line arguments.
91  libMesh::out << "Running " << argv[0];
92 
93  for (int i=1; i<argc; i++)
94  libMesh::out << " " << argv[i];
95 
96  libMesh::out << std::endl << std::endl;
97 
98  // Skip this 2D example if libMesh was compiled as 1D-only.
99  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
100 
101  // Create a mesh, with dimension to be overridden later, distributed
102  // across the default MPI communicator.
103  Mesh mesh(init.comm());
104 
105  // Use the MeshTools::Generation mesh generator to create a uniform
106  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
107  // to build a mesh of 15x15 QUAD9 elements. Building QUAD9
108  // elements instead of the default QUAD4's we used in example 2
109  // allow us to use higher-order approximation.
111  15, 15,
112  -1., 1.,
113  -1., 1.,
114  QUAD9);
115 
116  // Print information about the mesh to the screen.
117  // Note that 5x5 QUAD9 elements actually has 11x11 nodes,
118  // so this mesh is significantly larger than the one in example 2.
119  mesh.print_info();
120 
121  // Create an equation systems object.
122  EquationSystems equation_systems (mesh);
123 
124  // Declare the Poisson system and its variables.
125  // The Poisson system is another example of a steady system.
126  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
127 
128  // Adds the variable "u" to "Poisson". "u"
129  // will be approximated using second-order approximation.
130  equation_systems.get_system("Poisson").add_variable("u", SECOND);
131 
132  // Give the system a pointer to the matrix assembly
133  // function. This will be called when needed by the
134  // library.
135  equation_systems.get_system("Poisson").attach_assemble_function (assemble_poisson);
136 
137  // Initialize the data structures for the equation system.
138  equation_systems.init();
139 
140  // Prints information about the system to the screen.
141  equation_systems.print_info();
142 
143  // Solve the system "Poisson". Note that calling this
144  // member will assemble the linear system and invoke
145  // the default numerical solver. With PETSc the solver can be
146  // controlled from the command line. For example,
147  // you can invoke conjugate gradient with:
148  //
149  // ./introduction_ex3 -ksp_type cg
150  //
151  // You can also get a nice X-window that monitors the solver
152  // convergence with:
153  //
154  // ./introduction-ex3 -ksp_xmonitor
155  //
156  // if you linked against the appropriate X libraries when you
157  // built PETSc.
158  equation_systems.get_system("Poisson").solve();
159 
160 #if defined(LIBMESH_HAVE_VTK) && !defined(LIBMESH_ENABLE_PARMESH)
161 
162  // After solving the system write the solution
163  // to a VTK-formatted plot file.
164  VTKIO (mesh).write_equation_systems ("out.pvtu", equation_systems);
165 
166 #endif // #ifdef LIBMESH_HAVE_VTK
167 
168  // All done.
169  return 0;
170 }
This is the EquationSystems class.
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
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
MeshBase & mesh
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
This class implements reading and writing meshes in the VTK format.
Definition: vtk_io.h:60
SolverPackage default_solver_package()
Definition: libmesh.C:1050
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:1489
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
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:1305
void assemble_poisson(EquationSystems &es, const std::string &system_name)
OStreamProxy out
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50