libMesh
Functions
vector_fe_ex1.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 int component, 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 
)

Definition at line 261 of file miscellaneous_ex16.C.

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

Referenced by main().

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 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
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
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: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
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
bool has_static_condensation() const
Definition: system.C:2871
const T_sys & get_system(std::string_view name) const
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.
StaticCondensation & get_static_condensation()
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
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:372
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_system_matrix() const

◆ assemble_poisson() [2/2]

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

Definition at line 236 of file vector_fe_ex1.C.

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

238 {
239 
240  // It is a good idea to make sure we are assembling
241  // the proper system.
242  libmesh_assert_equal_to (system_name, "Poisson");
243 
244  // Get a constant reference to the mesh object.
245  const MeshBase & mesh = es.get_mesh();
246 
247  // The dimension that we are running
248  const unsigned int dim = mesh.mesh_dimension();
249 
250  // Get a reference to the LinearImplicitSystem we are solving
251  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem> ("Poisson");
252 
253  // A reference to the DofMap object for this system. The DofMap
254  // object handles the index translation from node and element numbers
255  // to degree of freedom numbers. We will talk more about the DofMap
256  // in future examples.
257  const DofMap & dof_map = system.get_dof_map();
258 
259  // Get a constant reference to the Finite Element type
260  // for the first (and only) variable in the system.
261  FEType fe_type = dof_map.variable_type(0);
262 
263  // Build a Finite Element object of the specified type.
264  // Note that FEVectorBase is a typedef for the templated FE
265  // class.
266  std::unique_ptr<FEVectorBase> fe (FEVectorBase::build(dim, fe_type));
267 
268  // A 2*p+1 order Gauss quadrature rule for numerical integration.
269  QGauss qrule (dim, fe_type.default_quadrature_order());
270 
271  // Tell the finite element object to use our quadrature rule.
272  fe->attach_quadrature_rule (&qrule);
273 
274  // Declare a special finite element object for
275  // boundary integration.
276  std::unique_ptr<FEVectorBase> fe_face (FEVectorBase::build(dim, fe_type));
277 
278  // Boundary integration requires one quadrature rule,
279  // with dimensionality one less than the dimensionality
280  // of the element.
281  QGauss qface(dim-1, fe_type.default_quadrature_order());
282 
283  // Tell the finite element object to use our
284  // quadrature rule.
285  fe_face->attach_quadrature_rule (&qface);
286 
287  // Here we define some references to cell-specific data that
288  // will be used to assemble the linear system.
289  //
290  // The element Jacobian * quadrature weight at each integration point.
291  const std::vector<Real> & JxW = fe->get_JxW();
292 
293  // The physical XY locations of the quadrature points on the element.
294  // These might be useful for evaluating spatially varying material
295  // properties at the quadrature points.
296  const std::vector<Point> & q_point = fe->get_xyz();
297 
298  // The element shape functions evaluated at the quadrature points.
299  // Notice the shape functions are a vector rather than a scalar.
300  const std::vector<std::vector<RealGradient>> & phi = fe->get_phi();
301 
302  // The element shape function gradients evaluated at the quadrature
303  // points. Notice that the shape function gradients are a tensor.
304  const std::vector<std::vector<RealTensor>> & dphi = fe->get_dphi();
305 
306  // Define data structures to contain the element matrix
307  // and right-hand-side vector contribution. Following
308  // basic finite element terminology we will denote these
309  // "Ke" and "Fe". These datatypes are templated on
310  // Number, which allows the same code to work for real
311  // or complex numbers.
314 
315  // This vector will hold the degree of freedom indices for
316  // the element. These define where in the global system
317  // the element degrees of freedom get mapped.
318  std::vector<dof_id_type> dof_indices;
319 
320  // The global system matrix
321  SparseMatrix<Number> & matrix = system.get_system_matrix();
322 
323  // Now we will loop over all the elements in the mesh.
324  // We will compute the element matrix and right-hand-side
325  // contribution.
326  //
327  // Element iterators are a nice way to iterate through all the
328  // elements, or all the elements that have some property. The
329  // iterator el will iterate from the first to the last element on
330  // the local processor. The iterator end_el tells us when to stop.
331  // It is smart to make this one const so that we don't accidentally
332  // mess it up! In case users later modify this program to include
333  // refinement, we will be safe and will only consider the active
334  // elements; hence we use a variant of the active_elem_iterator.
335  for (const auto & elem : mesh.active_local_element_ptr_range())
336  {
337  // Get the degree of freedom indices for the
338  // current element. These define where in the global
339  // matrix and right-hand-side this element will
340  // contribute to.
341  dof_map.dof_indices (elem, dof_indices);
342 
343  // Compute the element-specific data for the current
344  // element. This involves computing the location of the
345  // quadrature points (q_point) and the shape functions
346  // (phi, dphi) for the current element.
347  fe->reinit (elem);
348 
349  // Zero the element matrix and right-hand side before
350  // summing them. We use the resize member here because
351  // the number of degrees of freedom might have changed from
352  // the last element. Note that this will be the case if the
353  // element type is different (i.e. the last element was a
354  // triangle, now we are on a quadrilateral).
355 
356  // The DenseMatrix::resize() and the DenseVector::resize()
357  // members will automatically zero out the matrix and vector.
358  Ke.resize (dof_indices.size(),
359  dof_indices.size());
360 
361  Fe.resize (dof_indices.size());
362 
363  // We'll use an element-size-dependent h below, so the FDM error
364  // doesn't easily dominate FEM error.
365  const Real eps = 1.e-3 * elem->hmin();
366 
367  // Now loop over the quadrature points. This handles
368  // the numeric integration.
369  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
370  {
371  // Now we will build the element matrix. This involves
372  // a double loop to integrate the test functions (i) against
373  // the trial functions (j).
374  for (std::size_t i=0; i<phi.size(); i++)
375  for (std::size_t j=0; j<phi.size(); j++)
376  Ke(i,j) += JxW[qp] * dphi[i][qp].contract(dphi[j][qp]);
377 
378  // This is the end of the matrix summation loop
379  // Now we build the element right-hand-side contribution.
380  // This involves a single loop in which we integrate the
381  // "forcing function" in the PDE against the test functions.
382  {
383  const Real x = q_point[qp](0);
384  const Real y = q_point[qp](1);
385 
386  // "f" is the forcing function for the Poisson equation.
387  // In this case we set f to be a finite difference
388  // Laplacian approximation to the (known) exact solution.
389  //
390  // We will use the second-order accurate FD Laplacian
391  // approximation, which in 2D is
392  //
393  // u_xx + u_yy = (u(i,j-1) + u(i,j+1) +
394  // u(i-1,j) + u(i+1,j) +
395  // -4*u(i,j))/h^2
396 
397  // Since the value of the forcing function depends only
398  // on the location of the quadrature point (q_point[qp])
399  // we will compute it here, outside of the i-loop
400  const Real fx = -(exact_solution(0, x, y-eps) +
401  exact_solution(0, x, y+eps) +
402  exact_solution(0, x-eps, y) +
403  exact_solution(0, x+eps, y) -
404  4.*exact_solution(0, x, y))/eps/eps;
405 
406  const Real fy = -(exact_solution(1, x, y-eps) +
407  exact_solution(1, x, y+eps) +
408  exact_solution(1, x-eps, y) +
409  exact_solution(1, x+eps, y) -
410  4.*exact_solution(1, x, y))/eps/eps;
411 
412  const RealGradient f(fx, fy);
413 
414  for (std::size_t i=0; i<phi.size(); i++)
415  Fe(i) += JxW[qp]*f*phi[i][qp];
416  }
417  }
418 
419  // We have now reached the end of the RHS summation,
420  // and the end of quadrature point loop, so
421  // the interior element integration has
422  // been completed. However, we have not yet addressed
423  // boundary conditions. For this example we will only
424  // consider simple Dirichlet boundary conditions.
425  //
426  // There are several ways Dirichlet boundary conditions
427  // can be imposed. A simple approach, which works for
428  // interpolary bases like the standard Lagrange polynomials,
429  // is to assign function values to the
430  // degrees of freedom living on the domain boundary. This
431  // works well for interpolary bases, but is more difficult
432  // when non-interpolary (e.g Legendre or Hierarchic) bases
433  // are used.
434  //
435  // Dirichlet boundary conditions can also be imposed with a
436  // "penalty" method. In this case essentially the L2 projection
437  // of the boundary values are added to the matrix. The
438  // projection is multiplied by some large factor so that, in
439  // floating point arithmetic, the existing (smaller) entries
440  // in the matrix and right-hand-side are effectively ignored.
441  //
442  // This amounts to adding a term of the form (in latex notation)
443  //
444  // \frac{1}{\epsilon} \int_{\delta \Omega} \phi_i \phi_j = \frac{1}{\epsilon} \int_{\delta \Omega} u \phi_i
445  //
446  // where
447  //
448  // \frac{1}{\epsilon} is the penalty parameter, defined such that \epsilon << 1
449  {
450  // The following loop is over the sides of the element.
451  // If the element has no neighbor on a side then that
452  // side MUST live on a boundary of the domain.
453  for (auto side : elem->side_index_range())
454  if (elem->neighbor_ptr(side) == nullptr)
455  {
456  // The value of the shape functions at the quadrature
457  // points.
458  const std::vector<std::vector<RealGradient>> & phi_face = fe_face->get_phi();
459 
460  // The Jacobian * Quadrature Weight at the quadrature
461  // points on the face.
462  const std::vector<Real> & JxW_face = fe_face->get_JxW();
463 
464  // The XYZ locations (in physical space) of the
465  // quadrature points on the face. This is where
466  // we will interpolate the boundary value function.
467  const std::vector<Point> & qface_point = fe_face->get_xyz();
468 
469  // Compute the shape function values on the element
470  // face.
471  fe_face->reinit(elem, side);
472 
473  // Loop over the face quadrature points for integration.
474  for (unsigned int qp=0; qp<qface.n_points(); qp++)
475  {
476  // The location on the boundary of the current
477  // face quadrature point.
478  const Real xf = qface_point[qp](0);
479  const Real yf = qface_point[qp](1);
480 
481  // The penalty value. \frac{1}{\epsilon}
482  // in the discussion above.
483  const Real penalty = 1.e10;
484 
485  // The boundary values.
486  const RealGradient f(exact_solution(0, xf, yf),
487  exact_solution(1, xf, yf));
488 
489  // Matrix contribution of the L2 projection.
490  for (std::size_t i=0; i<phi_face.size(); i++)
491  for (std::size_t j=0; j<phi_face.size(); j++)
492  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
493 
494  // Right-hand-side contribution of the L2
495  // projection.
496  for (std::size_t i=0; i<phi_face.size(); i++)
497  Fe(i) += JxW_face[qp]*penalty*f*phi_face[i][qp];
498  }
499  }
500  }
501 
502  // We have now finished the quadrature point loop,
503  // and have therefore applied all the boundary conditions.
504 
505  // If this assembly program were to be used on an adaptive mesh,
506  // we would have to apply any hanging node constraint equations
507  //dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
508 
509  // The element matrix and right-hand-side are now built
510  // for this element. Add them to the global matrix and
511  // right-hand-side vector. The SparseMatrix::add_matrix()
512  // and NumericVector::add_vector() members do this for us.
513  matrix.add_matrix (Ke, dof_indices);
514  system.rhs->add_vector (Fe, dof_indices);
515  }
516 
517  // All done!
518 }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
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: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]...
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
Order default_quadrature_order() const
Definition: fe_type.h:371
const T_sys & get_system(std::string_view name) const
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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:372
Real exact_solution(const int component, const Real x, const Real y, const Real z=0.)
This is the exact solution that we are trying to obtain.
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_system_matrix() const

◆ exact_solution()

Real exact_solution ( const int  component,
const Real  x,
const Real  y,
const Real  z = 0. 
)

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 39 of file exact_solution.C.

Referenced by assemble_poisson().

43 {
44  static const Real pi = acos(-1.);
45 
46  switch (component)
47  {
48  case 0:
49  return cos(.5*pi*x)*sin(.5*pi*y)*cos(.5*pi*z);
50  case 1:
51  return sin(.5*pi*x)*cos(.5*pi*y)*cos(.5*pi*z);
52  case 2:
53  return sin(.5*pi*x)*cos(.5*pi*y)*cos(.5*pi*z)*cos(.5*pi*x*y*z);
54  default:
55  libmesh_error_msg("Invalid component = " << component);
56  }
57 
58  // dummy
59  return 0.0;
60 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real pi
.
Definition: libmesh.h:299

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 87 of file vector_fe_ex1.C.

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_poisson(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_square(), libMesh::System::calculate_norm(), libMesh::command_line_next(), libMesh::default_solver_package(), libMesh::EIGEN_SOLVERS, libMesh::FEInterface::field_type(), libMesh::LinearImplicitSystem::get_linear_solver(), libMesh::GMRES, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::L2, libMesh::libmesh_isnan(), mesh, libMesh::out, libMesh::Utility::pow(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::QUAD9, libMesh::Real, libMesh::LinearSolver< T >::set_solver_type(), libMesh::System::solution, libMesh::LinearImplicitSystem::solve(), libMesh::TYPE_VECTOR, libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::ExodusII_IO::write_equation_systems().

88 {
89  // Initialize libraries.
90  LibMeshInit init (argc, argv);
91 
92  // This example requires a linear solver package.
93  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
94  "--enable-petsc, --enable-trilinos, or --enable-eigen");
95 
96  // Brief message to the user regarding the program name
97  // and command line arguments.
98  libMesh::out << "Running " << argv[0];
99 
100  for (int i=1; i<argc; i++)
101  libMesh::out << " " << argv[i];
102 
103  libMesh::out << std::endl << std::endl;
104 
105  // Skip this 2D example if libMesh was compiled as 1D-only.
106  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
107 
108  // Get the mesh size from the command line.
109  const int nx = libMesh::command_line_next("-nx", 15),
110  ny = libMesh::command_line_next("-ny", 15);
111 
112  // Create a mesh, with dimension to be overridden later, on the
113  // default MPI communicator.
114  Mesh mesh(init.comm());
115 
116  // Use the MeshTools::Generation mesh generator to create a uniform
117  // 2D grid on the square [-1,1]^2. We instruct the mesh generator
118  // to build a mesh of 15x15 QUAD9 elements.
120  nx, ny,
121  -1., 1.,
122  -1., 1.,
123  QUAD9);
124 
125  // Print information about the mesh to the screen.
126  mesh.print_info();
127 
128  // Create an equation systems object.
129  EquationSystems equation_systems (mesh);
130 
131  // Declare the Poisson system and its variables.
132  // The Poisson system is another example of a steady system.
133  LinearImplicitSystem & poisson = equation_systems.add_system<LinearImplicitSystem> ("Poisson");
134 
135  // Read FE order from command line
136  std::string order_str = "SECOND";
137  order_str = libMesh::command_line_next("-o", order_str);
138  order_str = libMesh::command_line_next("-Order", order_str);
139  const Order order = Utility::string_to_enum<Order>(order_str);
140 
141  // Read FE Family from command line
142  std::string family_str = "LAGRANGE_VEC";
143  family_str = libMesh::command_line_next("-f", family_str);
144  family_str = libMesh::command_line_next("-FEFamily", family_str);
145  const FEFamily family = Utility::string_to_enum<FEFamily>(family_str);
146 
147  libmesh_error_msg_if(FEInterface::field_type(family) != TYPE_VECTOR,
148  "FE family " + family_str + " isn't vector-valued");
149 
150  // Adds the variable "u" to "Poisson". "u" will be approximated
151  // using the requested order of approximation and vector element
152  // type. Since the mesh is 2-D, "u" will have two components.
153  poisson.add_variable("u", order, family);
154 
155  // Give the system a pointer to the matrix assembly
156  // function. This will be called when needed by the
157  // library.
159 
160  // Initialize the data structures for the equation system.
161  equation_systems.init();
162 
163  // Prints information about the system to the screen.
164  equation_systems.print_info();
165 
166  // If we're using Eigen, the default BiCGStab solver does not seem
167  // to converge robustly for this system. Let's try some other
168  // settings for them.
171 
172  // Solve the system "Poisson". Note that calling this
173  // member will assemble the linear system and invoke
174  // the default numerical solver. With PETSc the solver can be
175  // controlled from the command line. For example,
176  // you can invoke conjugate gradient with:
177  //
178  // ./vector_fe_ex1 -ksp_type cg
179  //
180  // You can also get a nice X-window that monitors the solver
181  // convergence with:
182  //
183  // ./vector_fe_ex1 -ksp_xmonitor
184  //
185  // if you linked against the appropriate X libraries when you
186  // built PETSc.
187  poisson.solve();
188 
189  const Real l2_norm =
190  poisson.calculate_norm(*poisson.solution, 0, L2);
191 
192  libMesh::out << "L2 norm of solution = " << std::setprecision(17) <<
193  l2_norm << std::endl;
194 
195  libmesh_error_msg_if (libmesh_isnan(l2_norm),
196  "Failed to calculate solution");
197 
198  const Real error_in_norm = std::abs(l2_norm - sqrt(Real(2)));
199 
200  libMesh::out << "error in L2 norm = " << std::setprecision(17) <<
201  error_in_norm << std::endl;
202 
203  // The error in the norm converges faster than the norm of the
204  // error, at least until it gets low enough that floating-point
205  // roundoff (and the penalty method) kill us.
206  const int n = std::min(nx, ny);
207  const int p = static_cast<int>(order);
208  const Real expected_error_bound = 2*std::pow(n, -p*2);
209  libMesh::out << "error bound = " << std::setprecision(17) <<
210  expected_error_bound << std::endl;
211  libMesh::out << "error ratio = " << std::setprecision(17) <<
212  error_in_norm / expected_error_bound << std::endl;
213  libmesh_error_msg_if (error_in_norm > expected_error_bound,
214  "Error exceeds expected bound of " <<
215  expected_error_bound);
216 
217 #ifdef LIBMESH_HAVE_EXODUS_API
218  ExodusII_IO(mesh).write_equation_systems("out.e", equation_systems);
219 #endif
220 
221 #ifdef LIBMESH_HAVE_GMV
222  GMVIO(mesh).write_equation_systems("out.gmv", equation_systems);
223 #endif
224 
225  // All done.
226  return 0;
227 }
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.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
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
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 ...
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.
virtual LinearSolver< Number > * get_linear_solver() const override
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:46
virtual void solve() override
Assembles & solves the linear system A*x=b.
bool libmesh_isnan(T x)
SolverPackage default_solver_package()
Definition: libmesh.C:1117
T pow(const T &x)
Definition: utility.h:328
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
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1593
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
Definition: system.C:1724
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:1357
void set_solver_type(const SolverType st)
Sets the type of solver to use.
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:2161
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
FEFamily
defines an enum for finite element families.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
void assemble_poisson(EquationSystems &es, const std::string &system_name)