libMesh
Functions
subdomains_ex2.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=0., 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.

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 273 of file subdomains_ex2.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::as_range(), 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::pi, libMesh::PerfLog::pop(), libMesh::PerfLog::push(), libMesh::Real, libMesh::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), libMesh::ExplicitSystem::rhs, and value.

275 {
276  // It is a good idea to make sure we are assembling
277  // the proper system.
278  libmesh_assert_equal_to (system_name, "Poisson");
279 
280  // Declare a performance log. Give it a descriptive
281  // string to identify what part of the code we are
282  // logging, since there may be many PerfLogs in an
283  // application.
284  PerfLog perf_log ("Matrix Assembly");
285 
286  // Get a constant reference to the mesh object.
287  const MeshBase & mesh = es.get_mesh();
288 
289  // The dimension that we are running
290  const unsigned int dim = mesh.mesh_dimension();
291 
292  // Get a reference to the LinearImplicitSystem we are solving
293  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
294 
295  // A reference to the DofMap object for this system. The DofMap
296  // object handles the index translation from node and element numbers
297  // to degree of freedom numbers. We will talk more about the DofMap
298  // in future examples.
299  const DofMap & dof_map = system.get_dof_map();
300 
301  // Get a constant reference to the Finite Element type
302  // for the first (and only) variable in the system.
303  FEType fe_type = dof_map.variable_type(0);
304 
305  // Build a Finite Element object of the specified type. Since the
306  // FEBase::build() member dynamically creates memory we will
307  // store the object as a std::unique_ptr<FEBase>. This can be thought
308  // of as a pointer that will clean up after itself.
309  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
310 
311  // A 5th order Gauss quadrature rule for numerical integration.
312  QGauss qrule (dim, FIFTH);
313 
314  // Tell the finite element object to use our quadrature rule.
315  fe->attach_quadrature_rule (&qrule);
316 
317  // Declare a special finite element object for
318  // boundary integration.
319  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
320 
321  // Boundary integration requires one quadrature rule,
322  // with dimensionality one less than the dimensionality
323  // of the element.
324  QGauss qface(dim-1, FIFTH);
325 
326  // Tell the finite element object to use our
327  // quadrature rule.
328  fe_face->attach_quadrature_rule (&qface);
329 
330  // Here we define some references to cell-specific data that
331  // will be used to assemble the linear system.
332  // We begin with the element Jacobian * quadrature weight at each
333  // integration point.
334  const std::vector<Real> & JxW = fe->get_JxW();
335 
336  // The physical XY locations of the quadrature points on the element.
337  // These might be useful for evaluating spatially varying material
338  // properties at the quadrature points.
339  const std::vector<Point> & q_point = fe->get_xyz();
340 
341  // The element shape functions evaluated at the quadrature points.
342  const std::vector<std::vector<Real>> & phi = fe->get_phi();
343 
344  // The element shape function gradients evaluated at the quadrature
345  // points.
346  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
347 
348  // Define data structures to contain the element matrix
349  // and right-hand-side vector contribution. Following
350  // basic finite element terminology we will denote these
351  // "Ke" and "Fe". More detail is in example 3.
354 
355  // This vector will hold the degree of freedom indices for
356  // the element. These define where in the global system
357  // the element degrees of freedom get mapped.
358  std::vector<dof_id_type> dof_indices, dof_indices2;
359 
360  // The global system matrix
361  SparseMatrix<Number> & matrix = system.get_system_matrix();
362 
363  // Now we will loop over all the "local" elements in the mesh. We
364  // will compute the element matrix and right-hand-side contribution.
365  // See example 3 for a discussion of the element iterators. Here we
366  // only want to loop over elements that are owned by the local
367  // processor. This allows each processor to compute its components
368  // of the global matrix.
369  //
370  // "PARALLEL CHANGE"
371  for (const auto & elem : as_range(mesh.local_elements_begin(),
372  mesh.local_elements_end()))
373  {
374  // Start logging the shape function initialization.
375  // This is done through a simple function call with
376  // the name of the event to log.
377  perf_log.push("elem init");
378 
379  // Get the degree of freedom indices for the
380  // current element. These define where in the global
381  // matrix and right-hand-side this element will
382  // contribute to.
383  dof_map.dof_indices (elem, dof_indices, 0);
384  dof_map.dof_indices (elem, dof_indices2, 1);
385 
386  // libMesh::out << "dof_indices.size()="
387  // << dof_indices.size()
388  // << ", dof_indices2.size()="
389  // << dof_indices2.size()
390  // << std::endl;
391 
392  // Compute the element-specific data for the current
393  // element. This involves computing the location of the
394  // quadrature points (q_point) and the shape functions
395  // (phi, dphi) for the current element.
396  fe->reinit (elem);
397 
398  // Zero the element matrix and right-hand side before
399  // summing them. We use the resize member here because
400  // the number of degrees of freedom might have changed from
401  // the last element. Note that this will be the case if the
402  // element type is different (i.e. the last element was a
403  // triangle, now we are on a quadrilateral).
404  Ke.resize (std::max(dof_indices.size(), dof_indices2.size()),
405  std::max(dof_indices.size(), dof_indices2.size()));
406 
407  Fe.resize (std::max(dof_indices.size(), dof_indices2.size()));
408 
409  // Stop logging the shape function initialization.
410  // If you forget to stop logging an event the PerfLog
411  // object will probably catch the error and abort.
412  perf_log.pop("elem init");
413 
414  // Now we will build the element matrix. This involves
415  // a double loop to integrate the test functions (i) against
416  // the trial functions (j).
417  //
418  // We have split the numeric integration into two loops
419  // so that we can log the matrix and right-hand-side
420  // computation separately.
421  //
422  // Now start logging the element matrix computation
423  perf_log.push ("Ke");
424 
425  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
426  for (std::size_t i=0; i<phi.size(); i++)
427  for (std::size_t j=0; j<phi.size(); j++)
428  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
429 
430 
431  // Stop logging the matrix computation
432  perf_log.pop ("Ke");
433 
434  // Now we build the element right-hand-side contribution.
435  // This involves a single loop in which we integrate the
436  // "forcing function" in the PDE against the test functions.
437  //
438  // Start logging the right-hand-side computation
439  perf_log.push ("Fe");
440 
441  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
442  {
443  // fxy is the forcing function for the Poisson equation.
444  // In this case we set fxy to be a finite difference
445  // Laplacian approximation to the (known) exact solution.
446  //
447  // We will use the second-order accurate FD Laplacian
448  // approximation, which in 2D on a structured grid is
449  //
450  // u_xx + u_yy = (u(i-1,j) + u(i+1,j) +
451  // u(i,j-1) + u(i,j+1) +
452  // -4*u(i,j))/h^2
453  //
454  // Since the value of the forcing function depends only
455  // on the location of the quadrature point (q_point[qp])
456  // we will compute it here, outside of the i-loop
457  const Real x = q_point[qp](0);
458 #if LIBMESH_DIM > 1
459  const Real y = q_point[qp](1);
460 #else
461  const Real y = 0;
462 #endif
463 #if LIBMESH_DIM > 2
464  const Real z = q_point[qp](2);
465 #else
466  const Real z = 0;
467 #endif
468  const Real eps = 1.e-3;
469 
470  const Real uxx = (exact_solution(x-eps, y, z) +
471  exact_solution(x+eps, y, z) +
472  -2.*exact_solution(x, y, z))/eps/eps;
473 
474  const Real uyy = (exact_solution(x, y-eps, z) +
475  exact_solution(x, y+eps, z) +
476  -2.*exact_solution(x, y, z))/eps/eps;
477 
478  const Real uzz = (exact_solution(x, y, z-eps) +
479  exact_solution(x, y, z+eps) +
480  -2.*exact_solution(x, y, z))/eps/eps;
481 
482  Real fxy;
483  if (dim==1)
484  {
485  // In 1D, compute the rhs by differentiating the
486  // exact solution twice.
487  const Real pi = libMesh::pi;
488  fxy = (0.25*pi*pi)*sin(.5*pi*x);
489  }
490  else
491  {
492  fxy = - (uxx + uyy + ((dim==2) ? 0. : uzz));
493  }
494 
495  // Add the RHS contribution
496  for (std::size_t i=0; i<phi.size(); i++)
497  Fe(i) += JxW[qp]*fxy*phi[i][qp];
498  }
499 
500  // Stop logging the right-hand-side computation
501  perf_log.pop ("Fe");
502 
503  // At this point the interior element integration has
504  // been completed. However, we have not yet addressed
505  // boundary conditions. For this example we will only
506  // consider simple Dirichlet boundary conditions imposed
507  // via the penalty method. This is discussed at length in
508  // example 3.
509  {
510  // Start logging the boundary condition computation. We use a
511  // macro to log everything in this scope.
512  LOG_SCOPE_WITH("BCs", "", perf_log);
513 
514  // The following loops over the sides of the element.
515  // If the element has no neighbor on a side then that
516  // side MUST live on a boundary of the domain.
517  for (auto side : elem->side_index_range())
518  if ((elem->neighbor_ptr(side) == nullptr) ||
519  (elem->neighbor_ptr(side)->subdomain_id() != elem->subdomain_id()))
520  {
521 
522  // The penalty value. \frac{1}{\epsilon}
523  // in the discussion above.
524  const Real penalty = 1.e10;
525 
526  // The value of the shape functions at the quadrature
527  // points.
528  const std::vector<std::vector<Real>> & phi_face = fe_face->get_phi();
529 
530  // The Jacobian * Quadrature Weight at the quadrature
531  // points on the face.
532  const std::vector<Real> & JxW_face = fe_face->get_JxW();
533 
534  // The XYZ locations (in physical space) of the
535  // quadrature points on the face. This is where
536  // we will interpolate the boundary value function.
537  const std::vector<Point> & qface_point = fe_face->get_xyz();
538 
539  // Compute the shape function values on the element
540  // face.
541  fe_face->reinit(elem, side);
542 
543  // Loop over the face quadrature points for integration.
544  for (unsigned int qp=0; qp<qface.n_points(); qp++)
545  {
546  // The location on the boundary of the current
547  // face quadrature point.
548  const Real xf = qface_point[qp](0);
549 #if LIBMESH_DIM > 1
550  const Real yf = qface_point[qp](1);
551 #else
552  const Real yf = 0.;
553 #endif
554 #if LIBMESH_DIM > 2
555  const Real zf = qface_point[qp](2);
556 #else
557  const Real zf = 0.;
558 #endif
559 
560 
561  // The boundary value.
562  const Real value = exact_solution(xf, yf, zf);
563 
564  // Matrix contribution of the L2 projection.
565  for (std::size_t i=0; i<phi_face.size(); i++)
566  for (std::size_t j=0; j<phi_face.size(); j++)
567  Ke(i,j) += JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
568 
569  // Right-hand-side contribution of the L2
570  // projection.
571  for (std::size_t i=0; i<phi_face.size(); i++)
572  Fe(i) += JxW_face[qp]*penalty*value*phi_face[i][qp];
573  }
574  }
575  }
576 
577 
578  // The element matrix and right-hand-side are now built
579  // for this element. Add them to the global matrix and
580  // right-hand-side vector. The PetscMatrixBase::add_matrix()
581  // and PetscVector::add_vector() members do this for us.
582  // Start logging the insertion of the local (element)
583  // matrix and vector into the global matrix and vector
584  LOG_SCOPE_WITH("matrix insertion", "", perf_log);
585 
586  if (dof_indices.size())
587  {
588  matrix.add_matrix (Ke, dof_indices);
589  system.rhs->add_vector (Fe, dof_indices);
590  }
591 
592  if (dof_indices2.size())
593  {
594  matrix.add_matrix (Ke, dof_indices2);
595  system.rhs->add_vector (Fe, dof_indices2);
596  }
597  }
598 
599  // That's it. We don't need to do anything else to the
600  // PerfLog. When it goes out of scope (at this function return)
601  // it will print its log to the screen. Pretty easy, huh?
602 }
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.
const T_sys & get_system(std::string_view name) const
This is the MeshBase class.
Definition: mesh_base.h:75
The PerfLog class allows monitoring of specific events.
Definition: perf_log.h:145
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.
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
Definition: simple_range.h:57
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
Real exact_solution(const Real x, const Real y=0., 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
const Real pi
.
Definition: libmesh.h:299

◆ 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:299

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 92 of file subdomains_ex2.C.

References libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_poisson(), libMesh::System::attach_assemble_function(), libMesh::MeshTools::Generation::build_cube(), libMesh::LinearImplicitSystem::clear(), libMesh::command_line_next(), libMesh::default_solver_package(), dim, libMesh::EDGE2, libMesh::EDGE3, libMesh::EquationSystems::get_system(), libMesh::GnuPlotIO::GRID_ON, libMesh::HEX27, libMesh::HEX8, libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::INVALID_SOLVER_PACKAGE, mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::ParallelObject::processor_id(), libMesh::QUAD4, libMesh::QUAD9, libMesh::Real, libMesh::MeshOutput< MT >::write_equation_systems(), and libMesh::ExodusII_IO::write_equation_systems().

93 {
94  // Initialize libMesh and any dependent libraries, like in example 2.
95  LibMeshInit init (argc, argv);
96 
97  // This example requires a linear solver package.
98  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
99  "--enable-petsc, --enable-trilinos, or --enable-eigen");
100 
101  // Declare a performance log for the main program
102  // PerfLog perf_main("Main Program");
103 
104  // Create a GetPot object to parse the command line
105  GetPot command_line (argc, argv);
106 
107  // Check for proper calling arguments.
108  libmesh_error_msg_if(argc < 3, "Usage:\n" << "\t " << argv[0] << " -d 2(3)" << " -n 15");
109 
110  // Brief message to the user regarding the program name
111  // and command line arguments.
112  libMesh::out << "Running " << argv[0];
113 
114  for (int i=1; i<argc; i++)
115  libMesh::out << " " << argv[i];
116 
117  libMesh::out << std::endl << std::endl;
118 
119  // Read problem dimension from command line. Use int
120  // instead of unsigned since the GetPot overload is ambiguous
121  // otherwise.
122  const int dim = libMesh::command_line_next("-d", 2);
123 
124  // Skip higher-dimensional examples on a lower-dimensional libMesh build
125  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
126 
127  // Create a mesh with user-defined dimension on the default MPI
128  // communicator.
129  Mesh mesh (init.comm(), dim);
130 
131  // Read number of elements from command line
132  const int ps = libMesh::command_line_next("-n", 15);
133 
134  // Read FE order from command line
135  std::string order = "SECOND";
136  order = libMesh::command_line_next("-o", order);
137  order = libMesh::command_line_next("-Order", order);
138 
139  // Read FE Family from command line
140  std::string family = "LAGRANGE";
141  family = libMesh::command_line_next("-f", family);
142  family = libMesh::command_line_next("-FEFamily", family);
143 
144  // Cannot use discontinuous basis.
145  libmesh_error_msg_if(((family == "MONOMIAL") || (family == "XYZ")) &&
146  mesh.processor_id() == 0,
147  "This example requires a C^0 (or higher) FE basis.");
148 
149  // Use the MeshTools::Generation mesh generator to create a uniform
150  // grid on the square [-1,1]^D. We instruct the mesh generator
151  // to build a mesh of 8x8 Quad9 elements in 2D, or Hex27
152  // elements in 3D. Building these higher-order elements allows
153  // us to use higher-order approximation, as in example 3.
154 
155  Real halfwidth = dim > 1 ? 1. : 0.;
156  Real halfheight = dim > 2 ? 1. : 0.;
157 
158  if ((family == "LAGRANGE") && (order == "FIRST"))
159  {
160  // No reason to use high-order geometric elements if we are
161  // solving with low-order finite elements.
163  ps,
164  (dim>1) ? ps : 0,
165  (dim>2) ? ps : 0,
166  -1., 1.,
167  -halfwidth, halfwidth,
168  -halfheight, halfheight,
169  (dim==1) ? EDGE2 :
170  ((dim == 2) ? QUAD4 : HEX8));
171  }
172 
173  else
174  {
176  ps,
177  (dim>1) ? ps : 0,
178  (dim>2) ? ps : 0,
179  -1., 1.,
180  -halfwidth, halfwidth,
181  -halfheight, halfheight,
182  (dim==1) ? EDGE3 :
183  ((dim == 2) ? QUAD9 : HEX27));
184  }
185 
186  for (auto & elem : mesh.element_ptr_range())
187  {
188  const Point cent = elem->vertex_average();
189  if (dim > 1)
190  {
191  if ((cent(0) > 0) == (cent(1) > 0))
192  elem->subdomain_id() = 1;
193  }
194  else if (cent(0) > 0)
195  elem->subdomain_id() = 1;
196  }
197 
198  // Print information about the mesh to the screen.
199  mesh.print_info();
200 
201  // Create an equation systems object.
202  EquationSystems equation_systems (mesh);
203 
204  // Declare the system and its variables.
205  // Create a system named "Poisson"
206  LinearImplicitSystem & system =
207  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
208 
209 
210  std::set<subdomain_id_type> active_subdomains;
211 
212 
213  // Add the variable "u" to "Poisson". "u"
214  // will be approximated using second-order approximation.
215  active_subdomains.clear(); active_subdomains.insert(0);
216  system.add_variable("u",
217  Utility::string_to_enum<Order> (order),
218  Utility::string_to_enum<FEFamily>(family),
219  &active_subdomains);
220 
221  // Add the variable "v" to "Poisson". "v"
222  // will be approximated using second-order approximation.
223  active_subdomains.clear(); active_subdomains.insert(1);
224  system.add_variable("v",
225  Utility::string_to_enum<Order> (order),
226  Utility::string_to_enum<FEFamily>(family),
227  &active_subdomains);
228 
229  // Give the system a pointer to the matrix assembly
230  // function.
232 
233  // Initialize the data structures for the equation system.
234  equation_systems.init();
235 
236  // Print information about the system to the screen.
237  equation_systems.print_info();
238  mesh.print_info();
239 
240  // Solve the system "Poisson", just like example 2.
241  equation_systems.get_system("Poisson").solve();
242 
243  // After solving the system write the solution
244  // to a GMV-formatted plot file.
245  if (dim == 1)
246  {
247  GnuPlotIO plot(mesh, "Subdomains Example 2, 1D", GnuPlotIO::GRID_ON);
248  plot.write_equation_systems("gnuplot_script", equation_systems);
249  }
250  else
251  {
252 #ifdef LIBMESH_HAVE_EXODUS_API
254  "out_3.e" : "out_2.e", equation_systems);
255 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
256  }
257 
258  // All done.
259  return 0;
260 }
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 ...
MeshBase & mesh
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
SolverPackage default_solver_package()
Definition: libmesh.C:1117
This class implements writing meshes using GNUplot, designed for use only with 1D meshes...
Definition: gnuplot_io.h:43
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.
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 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
virtual void clear() override
Clear all the data structures associated with the system.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
processor_id_type processor_id() const
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
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)