34 #include "libmesh/getpot.h"
35 #include "libmesh/mesh_base.h"
36 #include "libmesh/error_vector.h"
41 : _error_estimator_options(input),
42 _mesh_adaptivity_options(input),
43 _refinement_type(INVALID)
70 if( refinement_stategy == std::string(
"error_tolerance") )
75 else if( refinement_stategy == std::string(
"nelem_target") )
79 if( !input.have_variable(
"MeshAdaptivity/nelem_target") )
81 std::cerr <<
"Error: Must specify nelem_target for" << std::endl
82 <<
" for error_tolerance refinement strategy."
89 else if( refinement_stategy == std::string(
"error_fraction") )
94 else if( refinement_stategy == std::string(
"elem_fraction") )
99 else if( refinement_stategy == std::string(
"mean_std_dev") )
106 std::cerr <<
"Error: Invalid refinement strategy " << refinement_stategy << std::endl
107 <<
" Valid refinement strategy options are: absolute_global_tolerance" << std::endl
108 <<
" error_tolerance" << std::endl
109 <<
" nelem_target" << std::endl
110 <<
" error_fraction" << std::endl
111 <<
" elem_fraction" << std::endl
112 <<
" mean_std_dev" << std::endl;
121 const libMesh::ErrorVector& error )
const
123 std::cout <<
"==========================================================" << std::endl
124 <<
"Checking convergence" << std::endl
125 <<
"==========================================================" << std::endl;
127 bool converged =
false;
129 libMesh::Real error_estimate = 0.0;
133 error_estimate = std::accumulate( error.begin(), error.end(), 0.0 );
137 error_estimate = error.l2_norm();
140 std::cout <<
"==========================================================" << std::endl
141 <<
"Error estimate = " << error_estimate << std::endl
142 <<
"==========================================================" << std::endl;
189 std::cerr <<
"Error: Invalid refinement option!" << std::endl;
207 std::cout <<
"==========================================================" << std::endl
208 <<
"Estimating error" << std::endl
209 <<
"==========================================================" << std::endl;
225 std::cout <<
"==========================================================" << std::endl
226 <<
"Performing Mesh Refinement" << std::endl
227 <<
"==========================================================" << std::endl;
236 std::cout <<
"==========================================================" << std::endl
237 <<
"Refined mesh to " << std::setw(12) << mesh.n_active_elem()
238 <<
" active elements" << std::endl
239 <<
" " << std::setw(16) << context.
system->n_active_dofs()
240 <<
" active dofs" << std::endl
241 <<
"==========================================================" << std::endl;
libMesh::Real coarsen_fraction() const
const std::string & refinement_strategy() const
unsigned int nelem_target() const
void set_refinement_type(const GetPot &input, const MeshAdaptivityOptions &mesh_adaptivity_options, RefinementFlaggingType &refinement_type)
bool coarsen_by_parents() const
SharedPtr< libMesh::EquationSystems > equation_system
unsigned int edge_level_mismatch_limit() const
bool plot_cell_errors() const
void build_mesh_refinement(libMesh::MeshBase &mesh)
libMesh::Real absolute_global_tolerance() const
unsigned int face_level_mismatch_limit() const
void perform_amr(SolverContext &context, const libMesh::ErrorVector &error)
void estimate_error_for_amr(SolverContext &context, libMesh::ErrorVector &error)
bool check_for_convergence(SolverContext &context, const libMesh::ErrorVector &error) const
SharedPtr< libMesh::ErrorEstimator > error_estimator
libMesh::Real coarsen_threshold() const
void flag_elements_for_refinement(const libMesh::ErrorVector &error)
libMesh::Real refine_fraction() const
Container for mesh adaptivity options.
unsigned int max_h_level() const
GRINS::MultiphysicsSystem * system
bool enforce_mismatch_limit_prior_to_refinement() const
unsigned int node_level_mismatch_limit() const
RefinementFlaggingType _refinement_type
Simple class to hold objects passed to Solver::solve.
libMesh::UniquePtr< libMesh::MeshRefinement > _mesh_refinement
const std::string & error_plot_prefix() const
MeshAdaptivityOptions _mesh_adaptivity_options