Go to the documentation of this file.
33 #include "libmesh/equation_systems.h"
34 #include "libmesh/error_vector.h"
35 #include "libmesh/getpot.h"
36 #include "libmesh/exodusII_io.h"
37 #include "libmesh/kelly_error_estimator.h"
38 #include "libmesh/mesh.h"
39 #include "libmesh/replicated_mesh.h"
40 #include "libmesh/distributed_mesh.h"
41 #include "libmesh/mesh_generation.h"
42 #include "libmesh/mesh_refinement.h"
43 #include "libmesh/uniform_refinement_estimator.h"
44 #include "libmesh/auto_ptr.h"
45 #include "libmesh/enum_solver_package.h"
46 #include "libmesh/enum_norm_type.h"
50 #include "libmesh/diff_solver.h"
51 #include "libmesh/petsc_diff_solver.h"
52 #include "libmesh/newton_solver.h"
53 #include "libmesh/euler_solver.h"
54 #include "libmesh/steady_solver.h"
55 #include "libmesh/partitioner.h"
61 int main (
int argc,
char ** argv)
68 "--enable-petsc, --enable-trilinos, or --enable-eigen");
71 #ifdef LIBMESH_DEFAULT_SINGLE_PRECISION
72 libmesh_example_requires(
false,
"--disable-singleprecision");
75 #ifndef LIBMESH_ENABLE_AMR
76 libmesh_example_requires(
false,
"--enable-amr");
80 #ifndef LIBMESH_ENABLE_DIRICHLET
81 libmesh_example_requires(
false,
"--enable-dirichlet");
90 GetPot infile(
"fem_system_ex1.in");
93 infile.parse_command_line(argc, argv);
96 const Real global_tolerance = infile(
"global_tolerance", 0.);
97 const unsigned int nelem_target = infile(
"n_elements", 400);
98 const bool transient = infile(
"transient",
true);
99 const Real deltat = infile(
"deltat", 0.005);
100 unsigned int n_timesteps = infile(
"n_timesteps", 20);
101 const unsigned int coarsegridsize = infile(
"coarsegridsize", 1);
102 const unsigned int coarserefinements = infile(
"coarserefinements", 0);
103 const unsigned int max_adaptivesteps = infile(
"max_adaptivesteps", 10);
104 const unsigned int dim = infile(
"dimension", 2);
105 const std::string slvr_type = infile(
"solver_type",
"newton");
106 const std::string mesh_type = infile(
"mesh_type" ,
"replicated");
108 #ifdef LIBMESH_HAVE_EXODUS_API
109 const unsigned int write_interval = infile(
"write_interval", 5);
113 libmesh_example_requires(
dim <= LIBMESH_DIM,
"2D/3D support");
120 std::shared_ptr<UnstructuredMesh>
mesh;
122 if (mesh_type ==
"distributed")
123 mesh = std::make_shared<DistributedMesh>(
init.comm());
124 else if (mesh_type ==
"replicated")
125 mesh = std::make_shared<ReplicatedMesh>(
init.comm());
127 libmesh_error_msg(
"Error: specified mesh_type not understood");
160 if (slvr_type ==
"petscdiff")
181 system.
time_solver = libmesh_make_unique<EulerSolver>(system);
184 system.
time_solver = libmesh_make_unique<SteadySolver>(system);
185 libmesh_assert_equal_to (n_timesteps, 1);
189 equation_systems.
init ();
195 if (slvr_type ==
"newton")
196 system.
time_solver->diff_solver() = libmesh_make_unique<NewtonSolver>(system);
197 else if (slvr_type ==
"petscdiff")
198 #if defined(LIBMESH_HAVE_PETSC) && defined(LIBMESH_HAVE_METAPHYSICL)
199 system.
time_solver->diff_solver() = libmesh_make_unique<PetscDiffSolver>(system);
201 libmesh_example_requires(
false,
"--enable-petsc --enable-metaphysicl-required");
204 libmesh_error_msg(
"Error: specified solver_type not understood");
209 solver.
quiet = infile(
"solver_quiet",
true);
212 infile(
"max_nonlinear_iterations", 15);
214 infile(
"relative_step_tolerance", 1.e-3);
216 infile(
"relative_residual_tolerance", 0.0);
218 infile(
"absolute_residual_tolerance", 0.0);
222 infile(
"max_linear_iterations", 50000);
224 infile(
"initial_linear_tolerance", 1.e-3);
231 for (
unsigned int t_step=0; t_step != n_timesteps; ++t_step)
241 unsigned int a_step = 0;
242 for (; a_step != max_adaptivesteps; ++a_step)
250 std::unique_ptr<ErrorEstimator> error_estimator;
254 if (global_tolerance != 0.)
258 libmesh_assert_equal_to (nelem_target, 0);
260 error_estimator = libmesh_make_unique<UniformRefinementEstimator>();
264 error_estimator->error_norm =
L2;
270 libmesh_assert_greater (nelem_target, 0);
276 error_estimator = libmesh_make_unique<KellyErrorEstimator>();
280 std::vector<Real> weights(2,1.0);
282 weights.push_back(1.0);
283 weights.push_back(0.0);
285 std::vector<FEMNormType>
286 norms(1, error_estimator->error_norm.type(0));
287 error_estimator->error_norm =
SystemNorm(norms, weights);
289 error_estimator->estimate_error(system, error);
298 if (global_tolerance != 0.)
303 if (global_tolerance != 0.)
310 if (global_tolerance != 0.)
314 if (global_error < global_tolerance)
325 equation_systems.
reinit();
326 a_step = max_adaptivesteps;
333 equation_systems.
reinit();
337 <<
" active elements and "
343 if (a_step == max_adaptivesteps)
353 #ifdef LIBMESH_HAVE_EXODUS_API
355 if ((t_step+1)%write_interval == 0)
357 std::ostringstream file_name;
373 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
375 #endif // #ifndef LIBMESH_ENABLE_AMR
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
bool refine_and_coarsen_elements()
Refines and coarsens user-requested elements.
virtual void solve() override
Invokes the solver associated with the system.
Real & coarsen_threshold()
The coarsen_threshold provides hysteresis in AMR/C strategies.
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.
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()
void write_timestep(const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr)
Writes out the solution at a specific timestep.
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
bool flag_elements_by_nelem_target(const ErrorVector &error_per_cell)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell.
virtual void init()
The initialization function.
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.
Implements (adaptive) mesh refinement algorithms for a MeshBase.
std::size_t n_active_dofs() const
int main(int argc, char **argv)
virtual void postprocess()
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
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.
void flag_elements_by_error_tolerance(const ErrorVector &error_per_cell)
Flags elements for coarsening and refinement based on the computed error passed in error_per_cell.
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.
dof_id_type & nelem_target()
If nelem_target is set to a nonzero value, methods like flag_elements_by_nelem_target() will attempt ...
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.
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
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 & coarsen_fraction()
The coarsen_fraction sets either a desired target or a desired maximum number of elements to flag for...
Real relative_residual_tolerance
Real & refine_fraction()
The refine_fraction sets either a desired target or a desired maximum number of elements to flag for ...
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.
Real & absolute_global_tolerance()
If absolute_global_tolerance is set to a nonzero value, methods like flag_elements_by_global_toleranc...
bool & coarsen_by_parents()
If coarsen_by_parents is true, complete groups of sibling elements (elements with the same parent) wi...
This class defines a norm/seminorm to be applied to a NumericVector which contains coefficients in a ...
void allow_remote_element_removal(bool allow)
If false is passed in then this mesh will no longer have remote elements deleted when being prepared ...
void allow_renumbering(bool allow)
If false is passed in then this mesh will no longer be renumbered when being prepared for use.
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 uniformly_refine(unsigned int n=1)
Uniformly refines the mesh n times.
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 ...
virtual std::unique_ptr< Partitioner > & partitioner()
A partitioner to use at each prepare_for_use()