Go to the documentation of this file.
   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"  
   47 #include "libmesh/enum_solver_package.h" 
   48 #include "libmesh/enum_norm_type.h" 
   51 #include "heatsystem.h" 
   52 #include "libmesh/diff_solver.h" 
   53 #include "libmesh/euler_solver.h" 
   54 #include "libmesh/steady_solver.h" 
   60 int main (
int argc, 
char ** argv)
 
   67                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   69 #ifndef LIBMESH_ENABLE_AMR 
   70   libmesh_example_requires(
false, 
"--enable-amr");
 
   74 #ifndef LIBMESH_ENABLE_DIRICHLET 
   75   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
   82   libmesh_example_requires(
sizeof(
Real) > 4, 
"--disable-singleprecision");
 
   85   GetPot infile(
"fem_system_ex4.in");
 
   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);
 
   97   libmesh_example_requires(
dim <= LIBMESH_DIM, 
"2D/3D support");
 
  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;
 
  146     std::set<boundary_id_type> bcids;
 
  157     if (elem->dim() < 
dim)
 
  158       elem->subdomain_id() = 1;
 
  160   mesh_refinement.uniformly_refine(coarserefinements);
 
  173   system.
time_solver = libmesh_make_unique<SteadySolver>(system);
 
  176   equation_systems.
init ();
 
  183   solver.
quiet = infile(
"solver_quiet", 
true);
 
  198   unsigned int a_step = 0;
 
  199   for (; a_step != max_adaptivesteps; ++a_step)
 
  207       std::unique_ptr<ErrorEstimator> error_estimator;
 
  211       if (global_tolerance != 0.)
 
  215           libmesh_assert_equal_to (nelem_target, 0);
 
  224           error_estimator.reset(u);
 
  230           libmesh_assert_greater (nelem_target, 0);
 
  236           error_estimator = libmesh_make_unique<KellyErrorEstimator>();
 
  239       error_estimator->estimate_error(system, error);
 
  248       if (global_tolerance != 0.)
 
  253       if (global_tolerance != 0.)
 
  260       if (global_tolerance != 0.)
 
  264           if (global_error < global_tolerance)
 
  266           mesh_refinement.flag_elements_by_error_tolerance(error);
 
  272           if (mesh_refinement.flag_elements_by_nelem_target(error))
 
  274               mesh_refinement.refine_and_coarsen_elements();
 
  275               equation_systems.
reinit();
 
  276               a_step = max_adaptivesteps;
 
  282       mesh_refinement.refine_and_coarsen_elements();
 
  283       equation_systems.
reinit();
 
  287                    << 
" active elements and " 
  293   if (a_step == max_adaptivesteps)
 
  301 #ifdef LIBMESH_HAVE_EXODUS_API 
  303     (
"out.e", equation_systems);
 
  304 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  306 #ifdef LIBMESH_HAVE_GMV 
  308     (
"out.gmv", equation_systems);
 
  309 #endif // #ifdef LIBMESH_HAVE_GMV 
  311 #ifdef LIBMESH_HAVE_FPARSER 
  314   const std::string exact_str = (
dim == 2) ?
 
  315     "sin(pi*x)*sin(pi*y)" : 
"sin(pi*x)*sin(pi*y)*sin(pi*z)";
 
  325   libmesh_assert_less(
err, 2e-3);
 
  327 #endif // #ifdef LIBMESH_HAVE_FPARSER 
  329 #endif // #ifndef LIBMESH_ENABLE_AMR 
  
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
 
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
 
A Function generated (via FParser) by parsing a mathematical expression.
 
virtual void solve() override
Invokes the solver associated with the system.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
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.
 
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
 
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
 
virtual Real l2_norm() const
 
The libMesh namespace provides an interface to certain functionality in the library.
 
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
 
SolverPackage default_solver_package()
 
This class implements writing meshes in the GMV format.
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
virtual SimpleRange< element_iterator > element_ptr_range()=0
 
Implements (adaptive) mesh refinement algorithms for a MeshBase.
 
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...
 
std::size_t n_active_dofs() const
 
Real relative_step_tolerance
 
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false.
 
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
 
virtual void init()
Initialize all the systems.
 
The ErrorVector is a specialization of the StatisticsVector for error data computed on a finite eleme...
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
Real l2_error(const std::string &sys_name, const std::string &unknown_name)
 
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(...
 
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...
 
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
 
This is the EquationSystems class.
 
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 ...
 
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
 
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
 
Real relative_residual_tolerance
 
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
 
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
 
virtual Real mean() const override
 
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
 
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.
 
int main(int argc, char **argv)
 
virtual dof_id_type n_active_elem() const =0
 
virtual T maximum() const
 
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...