libMesh
miscellaneous_ex10.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1>Miscellaneous Example 10 - Stitching meshes </h1>
21 // \author Dana Christen
22 // \date 2014
23 //
24 // This example shows how to generate a domain by stitching 8 cubic meshes
25 // together. Then a Poisson problem is solved on the stitched domain,
26 // and compared to the Poisson problem on a reference (unstitched) mesh
27 // to verify that the results match.
28 
29 
30 // C++ include files that we need
31 #include <iostream>
32 #include <algorithm>
33 #include <math.h>
34 #include <set>
35 #include <sstream>
36 #include <fstream>
37 
38 // libMesh includes
39 #include "libmesh/libmesh.h"
40 #include "libmesh/mesh.h"
41 #include "libmesh/mesh_generation.h"
42 #include "libmesh/linear_implicit_system.h"
43 #include "libmesh/equation_systems.h"
44 #include "libmesh/exact_solution.h"
45 #include "libmesh/exodusII_io.h"
46 #include "libmesh/fe.h"
47 #include "libmesh/quadrature_gauss.h"
48 #include "libmesh/dof_map.h"
49 #include "libmesh/sparse_matrix.h"
50 #include "libmesh/numeric_vector.h"
51 #include "libmesh/dense_matrix.h"
52 #include "libmesh/dense_vector.h"
53 #include "libmesh/elem.h"
54 #include "libmesh/dirichlet_boundaries.h"
55 #include "libmesh/zero_function.h"
56 #include "libmesh/libmesh_logging.h"
57 #include "libmesh/getpot.h"
58 #include "libmesh/error_vector.h"
59 #include "libmesh/kelly_error_estimator.h"
60 #include "libmesh/mesh_refinement.h"
61 #include "libmesh/enum_solver_package.h"
62 
63 using namespace libMesh;
64 
65 bool compare_elements(const UnstructuredMesh & mesh1,
66  const UnstructuredMesh & mesh2);
68  const std::string & system_name);
70  EquationSystems &);
71 
72 int main (int argc, char ** argv)
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 }
173 
175  EquationSystems & equation_systems)
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 }
237 
239  const std::string & libmesh_dbg_var(system_name))
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 }
virtual T maximum() const
Definition: statistics.C:62
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
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.
ConstFunction that simply returns 0.
Definition: zero_function.h:38
virtual Real l2_norm() const
Definition: statistics.C:37
static constexpr Real TOLERANCE
unsigned int dim
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
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 ...
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]...
Real & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
NumericVector< Number > * rhs
The system matrix.
This class allows one to associate Dirichlet boundary values with a given set of mesh boundary ids an...
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
void assemble_and_solve(MeshBase &, EquationSystems &)
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
The libMesh namespace provides an interface to certain functionality in the library.
virtual void solve() override
Assembles & solves the linear system A*x=b.
const T_sys & get_system(std::string_view name) const
unsigned int & max_h_level()
The max_h_level is the greatest refinement level an element should reach.
This is the MeshBase class.
Definition: mesh_base.h:75
SolverPackage default_solver_package()
Definition: libmesh.C:1117
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.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
void assemble_poisson(EquationSystems &es, const std::string &system_name)
The UnstructuredMesh class is derived from the MeshBase class.
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 compute_error(std::string_view sys_name, std::string_view unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
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
int main(int argc, char **argv)
unsigned int n_points() const
Definition: quadrature.h:131
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
Real l2_error(std::string_view sys_name, std::string_view unknown_name)
This class implements the Kelly error indicator which is based on the flux jumps between elements...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
OStreamProxy out
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
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
void flag_elements_by_error_fraction(const ErrorVector &error_per_cell, const Real refine_fraction=0.3, const Real coarsen_fraction=0.0, const unsigned int max_level=libMesh::invalid_uint)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell...
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.
void attach_reference_solution(const EquationSystems *es_fine)
Attach function similar to system.h which allows the user to attach a second EquationSystems object w...
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.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
const DofMap & get_dof_map() const
Definition: system.h:2374
bool compare_elements(const UnstructuredMesh &mesh1, const UnstructuredMesh &mesh2)
const SparseMatrix< Number > & get_system_matrix() const
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.