12 #include "libmesh/libmesh_config.h" 14 #ifdef LIBMESH_ENABLE_AMR 22 #include "libmesh/parallel_object.h" 25 #include "libmesh/mesh_refinement.h" 63 void init(
unsigned int steps,
unsigned int initial_steps);
72 void setParam(
const std::string & param_name,
const T & param_value);
145 bool adaptMesh(std::string marker_name = std::string());
244 ErrorVector &
getErrorVector(
const std::string & indicator_field);
318 template <
typename T>
322 if (param_name ==
"refine fraction")
328 else if (param_name ==
"coarsen fraction")
334 else if (param_name ==
"max h-level")
340 else if (param_name ==
"cycles_per_step")
342 else if (param_name ==
"recompute_markers_during_cycles")
345 mooseError(
"Invalid Param in adaptivity object");
347 #endif // LIBMESH_ENABLE_AMR bool adaptMesh(std::string marker_name=std::string())
Adapts the mesh based on the error estimator used.
bool initialAdaptMesh()
Used during initial adaptivity.
void doingPRefinement(bool doing_p_refinement, const MultiMooseEnum &disable_p_refinement_for_families)
Indicate whether the kind of adaptivity we're doing is p-refinement.
void setPrintMeshChanged(bool state=true)
std::shared_ptr< DisplacedProblem > _displaced_problem
Real _stop_time
When adaptivity stops.
std::string _marker_variable_name
Name of the marker variable if using the new adaptivity system.
std::unique_ptr< MeshRefinement > _displaced_mesh_refinement
A mesh refinement object for displaced mesh.
void setInterval(unsigned int interval)
Set the interval (number of timesteps) between refinement steps.
const unsigned int invalid_uint
Class for stuff related to variables.
std::unique_ptr< ErrorEstimator > _error_estimator
Error estimator to be used by the apps.
bool _initialized
on/off flag reporting if the adaptivity system has been initialized
static void uniformRefine(MooseMesh *mesh, unsigned int level=libMesh::invalid_uint)
Performs uniform refinement of the passed Mesh object.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
void setRecomputeMarkersFlag(const bool flag)
Set the flag to recompute markers during adaptivity cycles.
void setCyclesPerStep(const unsigned int &num)
Set the number of cycles_per_step.
unsigned int getCyclesPerStep() const
Pull out the number of cycles_per_step previously set through the AdaptivityAction.
bool _recompute_markers_during_cycles
Whether or not to recompute markers during adaptivity cycles.
void updateErrorVectors()
Update the ErrorVectors that have been requested through calls to getErrorVector().
void init(unsigned int steps, unsigned int initial_steps)
Initialize and turn on adaptivity for the simulation.
void setTimeActive(Real start_time, Real stop_time)
Sets the time when the adaptivity is active.
std::unique_ptr< ErrorVector > _error
Error vector for use with the error estimator.
void setMaxHLevel(unsigned int level)
Set the maximum refinement level (for the new Adaptivity system).
void setUseNewSystem()
Tells this object we're using the "new" adaptivity system.
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
MooseVariableFE< Real > MooseVariable
bool isOn()
Is adaptivity on?
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
std::string _initial_marker_variable_name
Name of the initial marker variable if using the new adaptivity system.
std::map< std::string, std::unique_ptr< ErrorVector > > _indicator_field_to_error_vector
Stores pointers to ErrorVectors associated with indicator field names.
bool _print_mesh_changed
True if we want to print out info when mesh has changed.
std::unique_ptr< MeshRefinement > _mesh_refinement
A mesh refinement object to be used either with initial refinement or with Adaptivity.
bool getRecomputeMarkersFlag() const
Pull out the _recompute_markers_during_cycles flag previously set through the AdaptivityAction.
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Adaptivity(FEProblemBase &fe_problem)
An inteface for the _console for outputting to the Console object.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
void setAdaptivityOn(bool state)
Allow adaptivity to be toggled programatically.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
void setErrorEstimator(const MooseEnum &error_estimator_name)
Set the error estimator.
unsigned int _steps
steps of adaptivity to perform
Interface for objects interacting with the PerfGraph.
void uniformRefineWithProjection()
Performs uniform refinement on the meshes in the current object.
unsigned int getMaxHLevel()
Return the maximum h-level.
unsigned int _initial_steps
the number of adaptivity steps to do at the beginning of simulation
bool _use_new_system
Whether or not to use the "new" adaptivity system.
unsigned int _interval
intreval between adaptivity runs
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Takes care of everything related to mesh adaptivity.
unsigned int getInitialSteps() const
Pull out the number of initial steps previously set by calling init()
ErrorVector & getErrorVector(const std::string &indicator_field)
Get an ErrorVector that will be filled up with values corresponding to the indicator field name passe...
unsigned int getSteps() const
Pull out the number of steps previously set by calling init()
unsigned int _cycles_per_step
The number of adaptivity cycles per step.
Real _start_time
When adaptivity start.
void setErrorNorm(SystemNorm &sys_norm)
Set the error norm (FIXME: improve description)
MooseVariableFE< VectorValue< Real > > VectorMooseVariable
bool isAdaptivityDue()
Query if an adaptivity step should be performed at the current time / time step.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
void setMarkerVariableName(std::string marker_field)
Sets the name of the field variable to actually use to flag elements for refinement / coarsening...
unsigned int _max_h_level
The maximum number of refinement levels.
void setParam(const std::string ¶m_name, const T ¶m_value)
Set adaptivity parameter.
void setInitialMarkerVariableName(std::string marker_field)
Sets the name of the field variable to actually use to flag elements for initial refinement / coarsen...
FEProblemBase & _fe_problem
bool isInitialized()
Returns whether or not Adaptivity::init() has ran.