libMesh
Functions
miscellaneous_ex10.C File Reference

Go to the source code of this file.

Functions

bool compare_elements (const UnstructuredMesh &mesh1, const UnstructuredMesh &mesh2)
 
void assemble_poisson (EquationSystems &es, const std::string &system_name)
 
void assemble_and_solve (MeshBase &, EquationSystems &)
 
int main (int argc, char **argv)
 
void assemble_poisson (EquationSystems &es, const std::string &libmesh_dbg_var(system_name))
 

Function Documentation

◆ assemble_and_solve()

void assemble_and_solve ( MeshBase mesh,
EquationSystems equation_systems 
)

Definition at line 174 of file miscellaneous_ex10.C.

References libMesh::DofMap::add_dirichlet_boundary(), libMesh::EquationSystems::add_system(), libMesh::System::add_variable(), assemble_poisson(), libMesh::System::attach_assemble_function(), libMesh::MeshRefinement::coarsen_fraction(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::FIRST, libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::System::get_dof_map(), libMesh::EquationSystems::init(), libMesh::StatisticsVector< T >::l2_norm(), libMesh::LAGRANGE, libMesh::LOCAL_VARIABLE_ORDER, libMesh::MeshRefinement::max_h_level(), libMesh::StatisticsVector< T >::maximum(), mesh, libMesh::out, libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::MeshRefinement::refine_fraction(), libMesh::EquationSystems::reinit(), and libMesh::LinearImplicitSystem::solve().

Referenced by main().

176 {
177  mesh.print_info();
178 
179  LinearImplicitSystem & system =
180  equation_systems.add_system<LinearImplicitSystem> ("Poisson");
181 
182 #ifdef LIBMESH_ENABLE_DIRICHLET
183  unsigned int u_var = system.add_variable("u", FIRST, LAGRANGE);
184 
186 
187  ZeroFunction<> zf;
188 
189  // the cube has boundaries IDs 0, 1, 2, 3, 4 and 5
190 
191  // Most DirichletBoundary users will want to supply a "locally
192  // indexed" functor
193  DirichletBoundary dirichlet_bc({0,1,2,3,4,5}, {u_var}, zf,
195  system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
196 #endif // LIBMESH_ENABLE_DIRICHLET
197 
198  equation_systems.init();
199  equation_systems.print_info();
200 
201 #ifdef LIBMESH_ENABLE_AMR
202  MeshRefinement mesh_refinement(mesh);
203 
204  mesh_refinement.refine_fraction() = 0.7;
205  mesh_refinement.coarsen_fraction() = 0.3;
206  mesh_refinement.max_h_level() = 5;
207 
208  const unsigned int max_r_steps = 2;
209 
210  for (unsigned int r_step=0; r_step<=max_r_steps; r_step++)
211  {
212  system.solve();
213  if (r_step != max_r_steps)
214  {
215  ErrorVector error;
216  KellyErrorEstimator error_estimator;
217 
218  error_estimator.estimate_error(system, error);
219 
220  libMesh::out << "Error estimate\nl2 norm = "
221  << error.l2_norm()
222  << "\nmaximum = "
223  << error.maximum()
224  << std::endl;
225 
226  mesh_refinement.flag_elements_by_error_fraction (error);
227 
228  mesh_refinement.refine_and_coarsen_elements();
229 
230  equation_systems.reinit();
231  }
232  }
233 #else
234  system.solve();
235 #endif
236 }
virtual T maximum() const
Definition: statistics.C:62
ConstFunction that simply returns 0.
Definition: zero_function.h:38
virtual Real l2_norm() const
Definition: statistics.C:37
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
virtual void solve() override
Assembles & solves the linear system A*x=b.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void assemble_poisson(EquationSystems &es, const std::string &system_name)
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
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
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
This class implements the Kelly error indicator which is based on the flux jumps between elements...
OStreamProxy out
void add_dirichlet_boundary(const DirichletBoundary &dirichlet_boundary)
Adds a copy of the specified Dirichlet boundary to the system.
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
This function uses the derived class&#39;s jump error estimate formula to estimate the error on each cell...
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
const DofMap & get_dof_map() const
Definition: system.h:2374

◆ 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 assemble_and_solve().

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 238 of file miscellaneous_ex10.C.

References libMesh::SparseMatrix< T >::add_matrix(), libMesh::NumericVector< T >::add_vector(), libMesh::FEGenericBase< OutputType >::build(), dim, 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::DenseVector< T >::resize(), libMesh::DenseMatrix< T >::resize(), and libMesh::ExplicitSystem::rhs.

240 {
241  libmesh_assert_equal_to (system_name, "Poisson");
242 
243  const MeshBase & mesh = es.get_mesh();
244  const unsigned int dim = mesh.mesh_dimension();
245  LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
246 
247  const DofMap & dof_map = system.get_dof_map();
248 
249  FEType fe_type = dof_map.variable_type(0);
250  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
251  QGauss qrule (dim, FIFTH);
252  fe->attach_quadrature_rule (&qrule);
253  std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
254  QGauss qface(dim-1, FIFTH);
255  fe_face->attach_quadrature_rule (&qface);
256 
257  const std::vector<Real> & JxW = fe->get_JxW();
258  const std::vector<std::vector<Real>> & phi = fe->get_phi();
259  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
260 
263 
264  std::vector<dof_id_type> dof_indices;
265  SparseMatrix<Number> & matrix = system.get_system_matrix();
266 
267  for (const auto & elem : mesh.active_local_element_ptr_range())
268  {
269  dof_map.dof_indices (elem, dof_indices);
270 
271  fe->reinit (elem);
272 
273  Ke.resize (dof_indices.size(),
274  dof_indices.size());
275 
276  Fe.resize (dof_indices.size());
277 
278  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
279  {
280  for (std::size_t i=0; i<phi.size(); i++)
281  {
282  Fe(i) += JxW[qp]*phi[i][qp];
283  for (std::size_t j=0; j<phi.size(); j++)
284  Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
285  }
286  }
287 
288  dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
289 
290  matrix.add_matrix (Ke, dof_indices);
291  system.rhs->add_vector (Fe, dof_indices);
292  }
293 }
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
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.
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
const DofMap & get_dof_map() const
Definition: system.h:2374
const SparseMatrix< Number > & get_system_matrix() const

◆ compare_elements()

bool compare_elements ( const UnstructuredMesh mesh1,
const UnstructuredMesh mesh2 
)

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 72 of file miscellaneous_ex10.C.

References assemble_and_solve(), libMesh::ExactSolution::attach_reference_solution(), libMesh::MeshTools::Generation::build_cube(), libMesh::command_line_next(), libMesh::ExactSolution::compute_error(), libMesh::default_solver_package(), dim, libMesh::HEX8, libMesh::TriangleWrapper::init(), libMesh::INVALID_SOLVER_PACKAGE, libMesh::ExactSolution::l2_error(), mesh, libMesh::out, libMesh::Real, libMesh::TOLERANCE, and libMesh::ExodusII_IO::write_equation_systems().

73 {
74  LibMeshInit init (argc, argv);
75 
76  // This example requires a linear solver package.
77  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
78  "--enable-petsc, --enable-trilinos, or --enable-eigen");
79 
80  // Check for proper calling arguments.
81  libmesh_error_msg_if(argc < 3, "Usage:\n" << "\t " << argv[0] << " -n 15");
82 
83  // Brief message to the user regarding the program name
84  // and command line arguments.
85  libMesh::out << "Running " << argv[0];
86 
87  for (int i=1; i<argc; i++)
88  libMesh::out << " " << argv[i];
89 
90  libMesh::out << std::endl << std::endl;
91 
92  // This is 3D-only problem
93  const unsigned int dim = 3;
94 
95  // Skip higher-dimensional examples on a lower-dimensional libMesh build
96  libmesh_example_requires(dim <= LIBMESH_DIM, "3D support");
97 
98  // We use Dirichlet boundary conditions here
99 #ifndef LIBMESH_ENABLE_DIRICHLET
100  libmesh_example_requires(false, "--enable-dirichlet");
101 #endif
102 
103  // Read number of elements used in each cube from command line
104  const int ps = libMesh::command_line_next("-n", 10);
105 
106  // Generate eight meshes that will be stitched
107  Mesh mesh (init.comm());
108  Mesh mesh1(init.comm());
109  Mesh mesh2(init.comm());
110  Mesh mesh3(init.comm());
111  Mesh mesh4(init.comm());
112  Mesh mesh5(init.comm());
113  Mesh mesh6(init.comm());
114  Mesh mesh7(init.comm());
115  Mesh nostitch_mesh(init.comm());
116  {
117  LOG_SCOPE("Initialize and create cubes", "main");
118  MeshTools::Generation::build_cube (mesh, ps, ps, ps, -1, 0, 0, 1, 0, 1, HEX8);
119  MeshTools::Generation::build_cube (mesh1, ps, ps, ps, 0, 1, 0, 1, 0, 1, HEX8);
120  MeshTools::Generation::build_cube (mesh2, ps, ps, ps, -1, 0, -1, 0, 0, 1, HEX8);
121  MeshTools::Generation::build_cube (mesh3, ps, ps, ps, 0, 1, -1, 0, 0, 1, HEX8);
122  MeshTools::Generation::build_cube (mesh4, ps, ps, ps, -1, 0, 0, 1, -1, 0, HEX8);
123  MeshTools::Generation::build_cube (mesh5, ps, ps, ps, 0, 1, 0, 1, -1, 0, HEX8);
124  MeshTools::Generation::build_cube (mesh6, ps, ps, ps, -1, 0, -1, 0, -1, 0, HEX8);
125  MeshTools::Generation::build_cube (mesh7, ps, ps, ps, 0, 1, -1, 0, -1, 0, HEX8);
126 
127  // Generate a single unstitched reference mesh
128  MeshTools::Generation::build_cube (nostitch_mesh, ps*2, ps*2, ps*2, -1, 1, -1, 1, -1, 1, HEX8);
129  }
130 
131  // We stitch the meshes in a hierarchical way.
132  {
133  LOG_SCOPE("Stitching", "main");
134  mesh.stitch_meshes(mesh1, 2, 4, TOLERANCE, true, true, false, false);
135  mesh2.stitch_meshes(mesh3, 2, 4, TOLERANCE, true, true, false, false);
136  mesh.stitch_meshes(mesh2, 1, 3, TOLERANCE, true, true, false, false);
137  mesh4.stitch_meshes(mesh5, 2, 4, TOLERANCE, true, true, false, false);
138  mesh6.stitch_meshes(mesh7, 2, 4, TOLERANCE, true, true, false, false);
139  mesh4.stitch_meshes(mesh6, 1, 3, TOLERANCE, true, true, false, false);
140  mesh.stitch_meshes(mesh4, 0, 5, TOLERANCE, true, true, false, false);
141  }
142 
143  EquationSystems equation_systems_stitch (mesh);
144  EquationSystems equation_systems_nostitch (nostitch_mesh);
145  {
146  LOG_SCOPE("Initialize and solve systems", "main");
147  assemble_and_solve(mesh, equation_systems_stitch);
148  assemble_and_solve(nostitch_mesh, equation_systems_nostitch);
149  }
150 
151  {
152  LOG_SCOPE("Result comparison", "main");
153  ExactSolution comparison(equation_systems_stitch);
154  comparison.attach_reference_solution(&equation_systems_nostitch);
155  comparison.compute_error("Poisson", "u");
156  Real error = comparison.l2_error("Poisson", "u");
157  libmesh_assert_less(error, TOLERANCE*sqrt(TOLERANCE));
158  libMesh::out << "L2 error between stitched and non-stitched cases: " << error << std::endl;
159  }
160 
161 #ifdef LIBMESH_HAVE_EXODUS_API
162  {
163  LOG_SCOPE("Output", "main");
164  ExodusII_IO(mesh).write_equation_systems("solution_stitch.exo",
165  equation_systems_stitch);
166  ExodusII_IO(nostitch_mesh).write_equation_systems("solution_nostitch.exo",
167  equation_systems_nostitch);
168  }
169 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
170 
171  return 0;
172 }
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
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.
static constexpr Real TOLERANCE
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
MeshBase & mesh
void assemble_and_solve(MeshBase &, EquationSystems &)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
SolverPackage default_solver_package()
Definition: libmesh.C:1117
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 init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
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.