1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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/replicated_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 ReplicatedMesh & mesh1,
66                       const ReplicatedMesh & mesh2);
67 void assemble_poisson(EquationSystems & es,
68                       const std::string & system_name);
69 void assemble_and_solve(MeshBase &,
70                         EquationSystems &);
71 
main(int argc,char ** argv)72 int main (int argc, char ** argv)
73 {
74   START_LOG("Initialize and create cubes", "main");
75   LibMeshInit init (argc, argv);
76 
77   // This example requires a linear solver package.
78   libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
79                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
80 
81   // Create a GetPot object to parse the command line
82   GetPot command_line (argc, argv);
83 
84   // Check for proper calling arguments.
85   libmesh_error_msg_if(argc < 3, "Usage:\n" << "\t " << argv[0] << " -n 15");
86 
87   // Brief message to the user regarding the program name
88   // and command line arguments.
89   libMesh::out << "Running " << argv[0];
90 
91   for (int i=1; i<argc; i++)
92     libMesh::out << " " << argv[i];
93 
94   libMesh::out << std::endl << std::endl;
95 
96   // This is 3D-only problem
97   const unsigned int dim = 3;
98 
99   // Skip higher-dimensional examples on a lower-dimensional libMesh build
100   libmesh_example_requires(dim <= LIBMESH_DIM, "3D support");
101 
102   // We use Dirichlet boundary conditions here
103 #ifndef LIBMESH_ENABLE_DIRICHLET
104   libmesh_example_requires(false, "--enable-dirichlet");
105 #endif
106 
107   // Read number of elements used in each cube from command line
108   int ps = 10;
109   if (command_line.search(1, "-n"))
110     ps = command_line.next(ps);
111 
112   // Generate eight meshes that will be stitched
113   ReplicatedMesh mesh (init.comm());
114   ReplicatedMesh mesh1(init.comm());
115   ReplicatedMesh mesh2(init.comm());
116   ReplicatedMesh mesh3(init.comm());
117   ReplicatedMesh mesh4(init.comm());
118   ReplicatedMesh mesh5(init.comm());
119   ReplicatedMesh mesh6(init.comm());
120   ReplicatedMesh mesh7(init.comm());
121   MeshTools::Generation::build_cube (mesh, ps, ps, ps, -1,    0,    0,  1,  0, 1, HEX8);
122   MeshTools::Generation::build_cube (mesh1, ps, ps, ps,    0,  1,    0,  1,  0, 1, HEX8);
123   MeshTools::Generation::build_cube (mesh2, ps, ps, ps, -1,    0, -1,    0,  0, 1, HEX8);
124   MeshTools::Generation::build_cube (mesh3, ps, ps, ps,    0,  1, -1,    0,  0, 1, HEX8);
125   MeshTools::Generation::build_cube (mesh4, ps, ps, ps, -1,    0,    0,  1, -1, 0, HEX8);
126   MeshTools::Generation::build_cube (mesh5, ps, ps, ps,    0,  1,    0,  1, -1, 0, HEX8);
127   MeshTools::Generation::build_cube (mesh6, ps, ps, ps, -1,    0, -1,    0, -1, 0, HEX8);
128   MeshTools::Generation::build_cube (mesh7, ps, ps, ps,    0,  1, -1,    0, -1, 0, HEX8);
129 
130   // Generate a single unstitched reference mesh
131   ReplicatedMesh nostitch_mesh(init.comm());
132   MeshTools::Generation::build_cube (nostitch_mesh, ps*2, ps*2, ps*2, -1, 1, -1, 1, -1, 1, HEX8);
133   STOP_LOG("Initialize and create cubes", "main");
134 
135   START_LOG("Stitching", "main");
136   // We stitch the meshes in a hierarchical way.
137   mesh.stitch_meshes(mesh1, 2, 4, TOLERANCE, true, true, false, false);
138   mesh2.stitch_meshes(mesh3, 2, 4, TOLERANCE, true, true, false, false);
139   mesh.stitch_meshes(mesh2, 1, 3, TOLERANCE, true, true, false, false);
140   mesh4.stitch_meshes(mesh5, 2, 4, TOLERANCE, true, true, false, false);
141   mesh6.stitch_meshes(mesh7, 2, 4, TOLERANCE, true, true, false, false);
142   mesh4.stitch_meshes(mesh6, 1, 3, TOLERANCE, true, true, false, false);
143   mesh.stitch_meshes(mesh4, 0, 5, TOLERANCE, true, true, false, false);
144   STOP_LOG("Stitching", "main");
145 
146   START_LOG("Initialize and solve systems", "main");
147   EquationSystems equation_systems_stitch (mesh);
148   assemble_and_solve(mesh, equation_systems_stitch);
149 
150   EquationSystems equation_systems_nostitch (nostitch_mesh);
151   assemble_and_solve(nostitch_mesh, equation_systems_nostitch);
152   STOP_LOG("Initialize and solve systems", "main");
153 
154   START_LOG("Result comparison", "main");
155   ExactSolution comparison(equation_systems_stitch);
156   comparison.attach_reference_solution(&equation_systems_nostitch);
157   comparison.compute_error("Poisson", "u");
158   Real error = comparison.l2_error("Poisson", "u");
159   libmesh_assert_less(error, TOLERANCE*sqrt(TOLERANCE));
160   libMesh::out << "L2 error between stitched and non-stitched cases: " << error << std::endl;
161   STOP_LOG("Result comparison", "main");
162 
163   START_LOG("Output", "main");
164 #ifdef LIBMESH_HAVE_EXODUS_API
165   ExodusII_IO(mesh).write_equation_systems("solution_stitch.exo",
166                                            equation_systems_stitch);
167   ExodusII_IO(nostitch_mesh).write_equation_systems("solution_nostitch.exo",
168                                                     equation_systems_nostitch);
169 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
170   STOP_LOG("Output", "main");
171 
172   return 0;
173 }
174 
assemble_and_solve(MeshBase & mesh,EquationSystems & equation_systems)175 void assemble_and_solve(MeshBase & mesh,
176                         EquationSystems & equation_systems)
177 {
178   mesh.print_info();
179 
180   LinearImplicitSystem & system =
181     equation_systems.add_system<LinearImplicitSystem> ("Poisson");
182 
183 #ifdef LIBMESH_ENABLE_DIRICHLET
184   unsigned int u_var = system.add_variable("u", FIRST, LAGRANGE);
185 
186   system.attach_assemble_function (assemble_poisson);
187 
188   // the cube has boundaries IDs 0, 1, 2, 3, 4 and 5
189   std::set<boundary_id_type> boundary_ids;
190   for (int j = 0; j<6; ++j)
191     boundary_ids.insert(j);
192 
193   // Create a vector storing the variable numbers which the BC applies to
194   std::vector<unsigned int> variables(1);
195   variables[0] = u_var;
196 
197   ZeroFunction<> zf;
198 
199   // Most DirichletBoundary users will want to supply a "locally
200   // indexed" functor
201   DirichletBoundary dirichlet_bc(boundary_ids, variables, zf,
202                                  LOCAL_VARIABLE_ORDER);
203   system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
204 #endif // LIBMESH_ENABLE_DIRICHLET
205 
206   equation_systems.init();
207   equation_systems.print_info();
208 
209 #ifdef LIBMESH_ENABLE_AMR
210   MeshRefinement mesh_refinement(mesh);
211 
212   mesh_refinement.refine_fraction()  = 0.7;
213   mesh_refinement.coarsen_fraction() = 0.3;
214   mesh_refinement.max_h_level()      = 5;
215 
216   const unsigned int max_r_steps = 2;
217 
218   for (unsigned int r_step=0; r_step<=max_r_steps; r_step++)
219     {
220       system.solve();
221       if (r_step != max_r_steps)
222         {
223           ErrorVector error;
224           KellyErrorEstimator error_estimator;
225 
226           error_estimator.estimate_error(system, error);
227 
228           libMesh::out << "Error estimate\nl2 norm = "
229                        << error.l2_norm()
230                        << "\nmaximum = "
231                        << error.maximum()
232                        << std::endl;
233 
234           mesh_refinement.flag_elements_by_error_fraction (error);
235 
236           mesh_refinement.refine_and_coarsen_elements();
237 
238           equation_systems.reinit();
239         }
240     }
241 #else
242   system.solve();
243 #endif
244 }
245 
assemble_poisson(EquationSystems & es,const std::string & libmesh_dbg_var (system_name))246 void assemble_poisson(EquationSystems & es,
247                       const std::string & libmesh_dbg_var(system_name))
248 {
249   libmesh_assert_equal_to (system_name, "Poisson");
250 
251   const MeshBase & mesh = es.get_mesh();
252   const unsigned int dim = mesh.mesh_dimension();
253   LinearImplicitSystem & system = es.get_system<LinearImplicitSystem>("Poisson");
254 
255   const DofMap & dof_map = system.get_dof_map();
256 
257   FEType fe_type = dof_map.variable_type(0);
258   std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
259   QGauss qrule (dim, FIFTH);
260   fe->attach_quadrature_rule (&qrule);
261   std::unique_ptr<FEBase> fe_face (FEBase::build(dim, fe_type));
262   QGauss qface(dim-1, FIFTH);
263   fe_face->attach_quadrature_rule (&qface);
264 
265   const std::vector<Real> & JxW = fe->get_JxW();
266   const std::vector<std::vector<Real>> & phi = fe->get_phi();
267   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
268 
269   DenseMatrix<Number> Ke;
270   DenseVector<Number> Fe;
271 
272   std::vector<dof_id_type> dof_indices;
273 
274   for (const auto & elem : mesh.active_local_element_ptr_range())
275     {
276       dof_map.dof_indices (elem, dof_indices);
277 
278       fe->reinit (elem);
279 
280       Ke.resize (dof_indices.size(),
281                  dof_indices.size());
282 
283       Fe.resize (dof_indices.size());
284 
285       for (unsigned int qp=0; qp<qrule.n_points(); qp++)
286         {
287           for (std::size_t i=0; i<phi.size(); i++)
288             {
289               Fe(i) += JxW[qp]*phi[i][qp];
290               for (std::size_t j=0; j<phi.size(); j++)
291                 Ke(i,j) += JxW[qp]*(dphi[i][qp]*dphi[j][qp]);
292             }
293         }
294 
295       dof_map.constrain_element_matrix_and_vector (Ke, Fe, dof_indices);
296 
297       system.matrix->add_matrix (Ke, dof_indices);
298       system.rhs->add_vector    (Fe, dof_indices);
299     }
300 }
301