libMesh
fem_system_ex4.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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>FEMSystem Example 4 - Mixed-dimension heat transfer equation
21 // with FEMSystem</h1>
22 // \author Roy Stogner
23 // \date 2015
24 //
25 // This example shows how elements of multiple dimensions can be
26 // linked and computed upon simultaneously within the
27 // DifferentiableSystem class framework
28 
29 // C++ includes
30 #include <iomanip>
31 
32 // Basic include files
33 #include "libmesh/elem.h"
34 #include "libmesh/equation_systems.h"
35 #include "libmesh/error_vector.h"
36 #include "libmesh/exact_solution.h"
37 #include "libmesh/getpot.h"
38 #include "libmesh/gmv_io.h"
39 #include "libmesh/exodusII_io.h"
40 #include "libmesh/kelly_error_estimator.h"
41 #include "libmesh/mesh.h"
42 #include "libmesh/mesh_generation.h"
43 #include "libmesh/mesh_refinement.h"
44 #include "libmesh/parsed_function.h"
45 #include "libmesh/uniform_refinement_estimator.h"
46 #include "libmesh/auto_ptr.h" // libmesh_make_unique
47 #include "libmesh/enum_solver_package.h"
48 #include "libmesh/enum_norm_type.h"
49 
50 // The systems and solvers we may use
51 #include "heatsystem.h"
52 #include "libmesh/diff_solver.h"
53 #include "libmesh/euler_solver.h"
54 #include "libmesh/steady_solver.h"
55 
56 // Bring in everything from the libMesh namespace
57 using namespace libMesh;
58 
59 // The main program.
60 int main (int argc, char ** argv)
61 {
62  // Initialize libMesh.
63  LibMeshInit init (argc, argv);
64 
65  // This example requires a linear solver package.
66  libmesh_example_requires(libMesh::default_solver_package() != INVALID_SOLVER_PACKAGE,
67  "--enable-petsc, --enable-trilinos, or --enable-eigen");
68 
69 #ifndef LIBMESH_ENABLE_AMR
70  libmesh_example_requires(false, "--enable-amr");
71 #else
72 
73  // We use Dirichlet boundary conditions here
74 #ifndef LIBMESH_ENABLE_DIRICHLET
75  libmesh_example_requires(false, "--enable-dirichlet");
76 #endif
77 
78  // This doesn't converge with Eigen BICGSTAB for some reason...
79  libmesh_example_requires(libMesh::default_solver_package() != EIGEN_SOLVERS, "--enable-petsc");
80 
81  // This doesn't converge without at least double precision
82  libmesh_example_requires(sizeof(Real) > 4, "--disable-singleprecision");
83 
84  // Parse the input file
85  GetPot infile("fem_system_ex4.in");
86 
87  // Read in parameters from the input file
88  const Real global_tolerance = infile("global_tolerance", 0.);
89  const unsigned int nelem_target = infile("n_elements", 400);
90  const Real deltat = infile("deltat", 0.005);
91  const unsigned int coarsegridsize = infile("coarsegridsize", 20);
92  const unsigned int coarserefinements = infile("coarserefinements", 0);
93  const unsigned int max_adaptivesteps = infile("max_adaptivesteps", 10);
94  const unsigned int dim = infile("dimension", 2);
95 
96  // Skip higher-dimensional examples on a lower-dimensional libMesh build
97  libmesh_example_requires(dim <= LIBMESH_DIM, "2D/3D support");
98 
99  // We have only defined 2 and 3 dimensional problems
100  libmesh_assert (dim == 2 || dim == 3);
101 
102  // Create a mesh, with dimension to be overridden later, distributed
103  // across the default MPI communicator.
104  Mesh mesh(init.comm());
105 
106  // And an object to refine it
107  MeshRefinement mesh_refinement(mesh);
108  mesh_refinement.coarsen_by_parents() = true;
109  mesh_refinement.absolute_global_tolerance() = global_tolerance;
110  mesh_refinement.nelem_target() = nelem_target;
111  mesh_refinement.refine_fraction() = 0.3;
112  mesh_refinement.coarsen_fraction() = 0.3;
113  mesh_refinement.coarsen_threshold() = 0.1;
114 
115  // Use the MeshTools::Generation mesh generator to create a uniform
116  // grid on the square or cube. We crop the domain at y=2/3 to allow
117  // for a homogeneous Neumann BC in our benchmark there.
118  boundary_id_type bcid = 3; // +y in 3D
119  if (dim == 2)
120  {
122  (mesh,
123  coarsegridsize,
124  coarsegridsize*2/3, // int arithmetic best we can do here
125  0., 1.,
126  0., 2./3.,
127  QUAD9);
128  bcid = 2; // +y in 2D
129  }
130  else if (dim == 3)
131  {
133  (mesh,
134  coarsegridsize,
135  coarsegridsize*2/3,
136  coarsegridsize,
137  0., 1.,
138  0., 2./3.,
139  0., 1.,
140  HEX27);
141  }
142 
143  {
144  // Add boundary elements corresponding to the +y boundary of our
145  // volume mesh
146  std::set<boundary_id_type> bcids;
147  bcids.insert(bcid);
150  }
151 
152  // To work around ExodusII file format limitations, we need elements
153  // of different dimensionality to belong to different subdomains.
154  // Our interior elements defaulted to subdomain id 0, so we'll set
155  // boundary elements to subdomain 1.
156  for (auto & elem : mesh.element_ptr_range())
157  if (elem->dim() < dim)
158  elem->subdomain_id() = 1;
159 
160  mesh_refinement.uniformly_refine(coarserefinements);
161 
162  // Print information about the mesh to the screen.
163  mesh.print_info();
164 
165  // Create an equation systems object.
166  EquationSystems equation_systems (mesh);
167 
168  // Declare the system "Heat" and its variables.
169  HeatSystem & system =
170  equation_systems.add_system<HeatSystem> ("Heat");
171 
172  // Solve this as a steady system
173  system.time_solver = libmesh_make_unique<SteadySolver>(system);
174 
175  // Initialize the system
176  equation_systems.init ();
177 
178  // Set the time stepping options
179  system.deltat = deltat;
180 
181  // And the nonlinear solver options
182  DiffSolver & solver = *(system.time_solver->diff_solver().get());
183  solver.quiet = infile("solver_quiet", true);
184  solver.verbose = !solver.quiet;
185  solver.max_nonlinear_iterations = infile("max_nonlinear_iterations", 15);
186  solver.relative_step_tolerance = infile("relative_step_tolerance", 1.e-3);
187  solver.relative_residual_tolerance = infile("relative_residual_tolerance", 0.0);
188  solver.absolute_residual_tolerance = infile("absolute_residual_tolerance", 0.0);
189 
190  // And the linear solver options
191  solver.max_linear_iterations = infile("max_linear_iterations", 50000);
192  solver.initial_linear_tolerance = infile("initial_linear_tolerance", 1.e-3);
193 
194  // Print information about the system to the screen.
195  equation_systems.print_info();
196 
197  // Adaptively solve the steady solution
198  unsigned int a_step = 0;
199  for (; a_step != max_adaptivesteps; ++a_step)
200  {
201  system.solve();
202 
203  system.postprocess();
204 
205  ErrorVector error;
206 
207  std::unique_ptr<ErrorEstimator> error_estimator;
208 
209  // To solve to a tolerance in this problem we
210  // need a better estimator than Kelly
211  if (global_tolerance != 0.)
212  {
213  // We can't adapt to both a tolerance and a mesh
214  // size at once
215  libmesh_assert_equal_to (nelem_target, 0);
216 
219 
220  // The lid-driven cavity problem isn't in H1, so
221  // lets estimate L2 error
222  u->error_norm = L2;
223 
224  error_estimator.reset(u);
225  }
226  else
227  {
228  // If we aren't adapting to a tolerance we need a
229  // target mesh size
230  libmesh_assert_greater (nelem_target, 0);
231 
232  // Kelly is a lousy estimator to use for a problem
233  // not in H1 - if we were doing more than a few
234  // timesteps we'd need to turn off or limit the
235  // maximum level of our adaptivity eventually
236  error_estimator = libmesh_make_unique<KellyErrorEstimator>();
237  }
238 
239  error_estimator->estimate_error(system, error);
240 
241  // Print out status at each adaptive step.
242  Real global_error = error.l2_norm();
243  libMesh::out << "Adaptive step "
244  << a_step
245  << ": "
246  << std::endl;
247 
248  if (global_tolerance != 0.)
249  libMesh::out << "Global_error = "
250  << global_error
251  << std::endl;
252 
253  if (global_tolerance != 0.)
254  libMesh::out << "Worst element error = "
255  << error.maximum()
256  << ", mean = "
257  << error.mean()
258  << std::endl;
259 
260  if (global_tolerance != 0.)
261  {
262  // If we've reached our desired tolerance, we
263  // don't need any more adaptive steps
264  if (global_error < global_tolerance)
265  break;
266  mesh_refinement.flag_elements_by_error_tolerance(error);
267  }
268  else
269  {
270  // If flag_elements_by_nelem_target returns true, this
271  // should be our last adaptive step.
272  if (mesh_refinement.flag_elements_by_nelem_target(error))
273  {
274  mesh_refinement.refine_and_coarsen_elements();
275  equation_systems.reinit();
276  a_step = max_adaptivesteps;
277  break;
278  }
279  }
280 
281  // Carry out the adaptive mesh refinement/coarsening
282  mesh_refinement.refine_and_coarsen_elements();
283  equation_systems.reinit();
284 
285  libMesh::out << "Refined mesh to "
286  << mesh.n_active_elem()
287  << " active elements and "
288  << equation_systems.n_active_dofs()
289  << " active dofs."
290  << std::endl;
291  }
292  // Do one last solve if necessary
293  if (a_step == max_adaptivesteps)
294  {
295  system.solve();
296 
297  system.postprocess();
298  }
299 
300 
301 #ifdef LIBMESH_HAVE_EXODUS_API
303  ("out.e", equation_systems);
304 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
305 
306 #ifdef LIBMESH_HAVE_GMV
308  ("out.gmv", equation_systems);
309 #endif // #ifdef LIBMESH_HAVE_GMV
310 
311 #ifdef LIBMESH_HAVE_FPARSER
312  // Check that we got close to the analytic solution
313  ExactSolution exact_sol(equation_systems);
314  const std::string exact_str = (dim == 2) ?
315  "sin(pi*x)*sin(pi*y)" : "sin(pi*x)*sin(pi*y)*sin(pi*z)";
316  ParsedFunction<Number> exact_func(exact_str);
317  exact_sol.attach_exact_value(0, &exact_func);
318  exact_sol.compute_error("Heat", "T");
319 
320  Real err = exact_sol.l2_error("Heat", "T");
321 
322  // Print out the error value
323  libMesh::out << "L2-Error is: " << err << std::endl;
324 
325  libmesh_assert_less(err, 2e-3);
326 
327 #endif // #ifdef LIBMESH_HAVE_FPARSER
328 
329 #endif // #ifndef LIBMESH_ENABLE_AMR
330 
331  // All done.
332  return 0;
333 }
libMesh::DifferentiableSystem::deltat
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
Definition: diff_system.h:260
libMesh::ExactSolution
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
Definition: exact_solution.h:73
libMesh::ParsedFunction
A Function generated (via FParser) by parsing a mathematical expression.
Definition: parsed_function.h:60
libMesh::FEMSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1046
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::MeshBase::get_boundary_info
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
Definition: mesh_base.h:132
libMesh::ErrorEstimator::error_norm
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
Definition: error_estimator.h:164
libMesh::StatisticsVector::l2_norm
virtual Real l2_norm() const
Definition: statistics.C:37
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::MeshTools::Generation::build_cube
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.
Definition: mesh_generation.C:298
libMesh::DiffSolver
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
Definition: diff_solver.h:70
libMesh::default_solver_package
SolverPackage default_solver_package()
Definition: libmesh.C:993
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::GMVIO
This class implements writing meshes in the GMV format.
Definition: gmv_io.h:54
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::boundary_id_type
int8_t boundary_id_type
Definition: id_types.h:51
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DifferentiableSystem::time_solver
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
Definition: diff_system.h:233
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
libMesh::MeshBase::element_ptr_range
virtual SimpleRange< element_iterator > element_ptr_range()=0
libMesh::MeshRefinement
Implements (adaptive) mesh refinement algorithms for a MeshBase.
Definition: mesh_refinement.h:61
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::MeshTools::Generation::build_square
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
Definition: mesh_generation.C:1501
libMesh::ExactSolution::attach_exact_value
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
Clone and attach an arbitrary functor which computes the exact value of the system sys_num solution a...
Definition: exact_solution.C:123
libMesh::HEX27
Definition: enum_elem_type.h:49
HeatSystem
Definition: heatsystem.h:30
libMesh::EquationSystems::n_active_dofs
std::size_t n_active_dofs() const
Definition: equation_systems.C:1362
libMesh::DiffSolver::relative_step_tolerance
Real relative_step_tolerance
Definition: diff_solver.h:204
libMesh::DiffSolver::verbose
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false.
Definition: diff_solver.h:168
libMesh::FEMSystem::postprocess
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
Definition: fem_system.C:1097
libMesh::INVALID_SOLVER_PACKAGE
Definition: enum_solver_package.h:43
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::ErrorVector
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
Definition: error_vector.h:50
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::UniformRefinementEstimator
This class implements a `‘brute force’' error estimator which integrates differences between the curr...
Definition: uniform_refinement_estimator.h:45
libMesh::ExactSolution::l2_error
Real l2_error(const std::string &sys_name, const std::string &unknown_name)
Definition: exact_solution.C:333
libMesh::ExactSolution::compute_error
void compute_error(const std::string &sys_name, const std::string &unknown_name)
Computes and stores the error in the solution value e = u-u_h, the gradient grad(e) = grad(u) - grad(...
Definition: exact_solution.C:227
libMesh::BoundaryInfo::add_elements
void add_elements(const std::set< boundary_id_type > &requested_boundary_ids, UnstructuredMesh &boundary_mesh)
Generates elements along the boundary of our _mesh, which use pre-existing nodes on the boundary_mesh...
Definition: boundary_info.C:370
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
libMesh::DiffSolver::quiet
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
Definition: diff_solver.h:162
libMesh::EquationSystems::reinit
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
Definition: equation_systems.C:121
libMesh::DiffSolver::relative_residual_tolerance
Real relative_residual_tolerance
Definition: diff_solver.h:192
libMesh::DiffSolver::absolute_residual_tolerance
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
Definition: diff_solver.h:191
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
libMesh::ErrorVector::mean
virtual Real mean() const override
Definition: error_vector.C:69
libMesh::DiffSolver::max_linear_iterations
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
Definition: diff_solver.h:148
libMesh::L2
Definition: enum_norm_type.h:36
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::err
OStreamProxy err
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::DiffSolver::max_nonlinear_iterations
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
Definition: diff_solver.h:156
libMesh::MeshBase::prepare_for_use
void prepare_for_use(const bool skip_renumber_nodes_and_elements=false, const bool skip_find_neighbors=false)
Prepare a newly ecreated (or read) mesh for use.
Definition: mesh_base.C:318
main
int main(int argc, char **argv)
Definition: fem_system_ex4.C:60
libMesh::out
OStreamProxy out
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::MeshBase::n_active_elem
virtual dof_id_type n_active_elem() const =0
libMesh::StatisticsVector::maximum
virtual T maximum() const
Definition: statistics.C:62
libMesh::DiffSolver::initial_linear_tolerance
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...
Definition: diff_solver.h:210