www.mooseframework.org
Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Adaptivity Class Reference

Takes care of everything related to mesh adaptivity. More...

#include <Adaptivity.h>

Inheritance diagram for Adaptivity:
[legend]

Public Member Functions

 Adaptivity (FEProblemBase &fe_problem)
 
virtual ~Adaptivity ()
 
void init (unsigned int steps, unsigned int initial_steps)
 Initialize and turn on adaptivity for the simulation. More...
 
template<typename T >
void setParam (const std::string &param_name, const T &param_value)
 Set adaptivity parameter. More...
 
void setErrorEstimator (const MooseEnum &error_estimator_name)
 Set the error estimator. More...
 
void setErrorNorm (SystemNorm &sys_norm)
 Set the error norm (FIXME: improve description) More...
 
void setPrintMeshChanged (bool state=true)
 
unsigned int getInitialSteps () const
 Pull out the number of initial steps previously set by calling init() More...
 
unsigned int getSteps () const
 Pull out the number of steps previously set by calling init() More...
 
unsigned int getCyclesPerStep () const
 Pull out the number of cycles_per_step previously set through the AdaptivityAction. More...
 
void setCyclesPerStep (const unsigned int &num)
 Set the number of cycles_per_step. More...
 
bool getRecomputeMarkersFlag () const
 Pull out the _recompute_markers_during_cycles flag previously set through the AdaptivityAction. More...
 
void setRecomputeMarkersFlag (const bool flag)
 Set the flag to recompute markers during adaptivity cycles. More...
 
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. More...
 
bool adaptMesh (std::string marker_name=std::string())
 Adapts the mesh based on the error estimator used. More...
 
bool initialAdaptMesh ()
 Used during initial adaptivity. More...
 
void uniformRefineWithProjection ()
 Performs uniform refinement on the meshes in the current object. More...
 
void setAdaptivityOn (bool state)
 Allow adaptivity to be toggled programatically. More...
 
bool isOn ()
 Is adaptivity on? More...
 
bool isInitialized ()
 Returns whether or not Adaptivity::init() has ran. More...
 
void setTimeActive (Real start_time, Real stop_time)
 Sets the time when the adaptivity is active. More...
 
void setUseNewSystem ()
 Tells this object we're using the "new" adaptivity system. More...
 
void setMarkerVariableName (std::string marker_field)
 Sets the name of the field variable to actually use to flag elements for refinement / coarsening. More...
 
void setInitialMarkerVariableName (std::string marker_field)
 Sets the name of the field variable to actually use to flag elements for initial refinement / coarsening. More...
 
void setMaxHLevel (unsigned int level)
 Set the maximum refinement level (for the new Adaptivity system). More...
 
unsigned int getMaxHLevel ()
 Return the maximum h-level. More...
 
void setInterval (unsigned int interval)
 Set the interval (number of timesteps) between refinement steps. More...
 
ErrorVector & getErrorVector (const std::string &indicator_field)
 Get an ErrorVector that will be filled up with values corresponding to the indicator field name passed in. More...
 
void updateErrorVectors ()
 Update the ErrorVectors that have been requested through calls to getErrorVector(). More...
 
bool isAdaptivityDue ()
 Query if an adaptivity step should be performed at the current time / time step. More...
 
PerfGraphperfGraph ()
 Get the PerfGraph. More...
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static void uniformRefine (MooseMesh *mesh, unsigned int level=libMesh::invalid_uint)
 Performs uniform refinement of the passed Mesh object. More...
 
static InputParameters validParams ()
 

Public Attributes

const ConsoleStream _console
 An instance of helper class to write streams to the Console objects. More...
 

Protected Member Functions

PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 Call to register a named section for timing. More...
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 Call to register a named section for timing. More...
 
std::string timedSectionName (const std::string &section_name) const
 

Protected Attributes

FEProblemBase_fe_problem
 
MooseMesh_mesh
 
bool _mesh_refinement_on
 on/off flag reporting if the adaptivity is being used More...
 
bool _initialized
 on/off flag reporting if the adaptivity system has been initialized More...
 
std::unique_ptr< MeshRefinement > _mesh_refinement
 A mesh refinement object to be used either with initial refinement or with Adaptivity. More...
 
std::unique_ptr< ErrorEstimator > _error_estimator
 Error estimator to be used by the apps. More...
 
std::unique_ptr< ErrorVector > _error
 Error vector for use with the error estimator. More...
 
std::shared_ptr< DisplacedProblem_displaced_problem
 
std::unique_ptr< MeshRefinement > _displaced_mesh_refinement
 A mesh refinement object for displaced mesh. More...
 
unsigned int _initial_steps
 the number of adaptivity steps to do at the beginning of simulation More...
 
unsigned int _steps
 steps of adaptivity to perform More...
 
bool _print_mesh_changed
 True if we want to print out info when mesh has changed. More...
 
Real_t
 Time. More...
 
int_step
 Time Step. More...
 
unsigned int _interval
 intreval between adaptivity runs More...
 
Real _start_time
 When adaptivity start. More...
 
Real _stop_time
 When adaptivity stops. More...
 
unsigned int _cycles_per_step
 The number of adaptivity cycles per step. More...
 
bool _use_new_system
 Whether or not to use the "new" adaptivity system. More...
 
std::string _marker_variable_name
 Name of the marker variable if using the new adaptivity system. More...
 
std::string _initial_marker_variable_name
 Name of the initial marker variable if using the new adaptivity system. More...
 
unsigned int _max_h_level
 The maximum number of refinement levels. More...
 
bool _recompute_markers_during_cycles
 Whether or not to recompute markers during adaptivity cycles. More...
 
std::map< std::string, std::unique_ptr< ErrorVector > > _indicator_field_to_error_vector
 Stores pointers to ErrorVectors associated with indicator field names. More...
 
bool _p_refinement_flag = false
 
MooseApp_pg_moose_app
 The MooseApp that owns the PerfGraph. More...
 
const std::string _prefix
 A prefix to use for all sections. More...
 
const Parallel::Communicator & _communicator
 

Detailed Description

Takes care of everything related to mesh adaptivity.

Definition at line 49 of file Adaptivity.h.

Constructor & Destructor Documentation

◆ Adaptivity()

Adaptivity::Adaptivity ( FEProblemBase fe_problem)

Definition at line 31 of file Adaptivity.C.

32  : ConsoleStreamInterface(fe_problem.getMooseApp()),
33  PerfGraphInterface(fe_problem.getMooseApp().perfGraph(), "Adaptivity"),
34  ParallelObject(fe_problem.getMooseApp()),
35  _fe_problem(fe_problem),
37  _mesh_refinement_on(false),
38  _initialized(false),
39  _initial_steps(0),
40  _steps(0),
41  _print_mesh_changed(false),
42  _t(_fe_problem.time()),
44  _interval(1),
48  _use_new_system(false),
49  _max_h_level(0),
51 {
52 }
ParallelObject(const Parallel::Communicator &comm_in)
Real _stop_time
When adaptivity stops.
Definition: Adaptivity.h:293
virtual Real & time() const
bool _initialized
on/off flag reporting if the adaptivity system has been initialized
Definition: Adaptivity.h:263
bool _recompute_markers_during_cycles
Whether or not to recompute markers during adaptivity cycles.
Definition: Adaptivity.h:310
MooseApp & getMooseApp() const
Get the MooseApp this class is associated with.
Definition: MooseBase.h:45
MooseMesh & _mesh
Definition: Adaptivity.h:258
auto max(const L &left, const R &right)
Real & _t
Time.
Definition: Adaptivity.h:285
bool _print_mesh_changed
True if we want to print out info when mesh has changed.
Definition: Adaptivity.h:282
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Definition: Adaptivity.h:261
ConsoleStreamInterface(MooseApp &app)
A class for providing a helper stream object for writting message to all the Output objects...
unsigned int _steps
steps of adaptivity to perform
Definition: Adaptivity.h:279
virtual int & timeStep() const
unsigned int _initial_steps
the number of adaptivity steps to do at the beginning of simulation
Definition: Adaptivity.h:277
bool _use_new_system
Whether or not to use the "new" adaptivity system.
Definition: Adaptivity.h:298
unsigned int _interval
intreval between adaptivity runs
Definition: Adaptivity.h:289
virtual MooseMesh & mesh() override
unsigned int _cycles_per_step
The number of adaptivity cycles per step.
Definition: Adaptivity.h:295
Real _start_time
When adaptivity start.
Definition: Adaptivity.h:291
int & _step
Time Step.
Definition: Adaptivity.h:287
PerfGraphInterface(const MooseObject *moose_object)
For objects that are MooseObjects with a default prefix of type()
unsigned int _max_h_level
The maximum number of refinement levels.
Definition: Adaptivity.h:307
PerfGraph & perfGraph()
Get the PerfGraph for this app.
Definition: MooseApp.h:133
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257

◆ ~Adaptivity()

Adaptivity::~Adaptivity ( )
virtual

Definition at line 54 of file Adaptivity.C.

54 {}

Member Function Documentation

◆ adaptMesh()

bool Adaptivity::adaptMesh ( std::string  marker_name = std::string())

Adapts the mesh based on the error estimator used.

Returns
a boolean that indicates whether the mesh was changed

Definition at line 123 of file Adaptivity.C.

Referenced by FEProblemBase::adaptMesh(), and initialAdaptMesh().

124 {
125  TIME_SECTION("adaptMesh", 3, "Adapting Mesh");
126 
127  // If the marker name is supplied, use it. Otherwise, use the one in _marker_variable_name
128  if (marker_name.empty())
129  marker_name = _marker_variable_name;
130 
131  bool mesh_changed = false;
132 
133  // If mesh adaptivity is carried out in a distributed (scalable) way
134  bool distributed_adaptivity = false;
135 
136  if (_use_new_system)
137  {
138  if (!marker_name.empty()) // Only flag if a marker variable name has been set
139  {
140  _mesh_refinement->clean_refinement_flags();
141 
142  std::vector<Number> serialized_solution;
143 
144  auto distributed_mesh = dynamic_cast<DistributedMesh *>(&_fe_problem.mesh().getMesh());
145 
146  // Element range
147  std::unique_ptr<ConstElemRange> all_elems;
148  // If the mesh is distributed and we do not do "gather to zero" or "allgather".
149  // Then it is safe to not serialize solution.
150  // Some output idiom (Exodus) will do "gather to zero". That being said,
151  // if you have exodus output on, mesh adaptivty is not scalable.
152  if (distributed_mesh && !distributed_mesh->is_serial_on_zero())
153  {
154  // We update here to make sure local solution is up-to-date
156  distributed_adaptivity = true;
157 
158  // We can not assume that geometric and algebraic ghosting functors cover
159  // the same set of elements/nodes. That being said, in general,
160  // we would expect G(e) > A(e). Here G(e) is the set of elements reserved
161  // by the geometric ghosting functors, and A(e) corresponds to
162  // the one covered by the algebraic ghosting functors.
163  // Therefore, we have to work only on local elements instead of
164  // ghosted + local elements. The ghosted solution might not be enough
165  // for ghosted+local elements. But it is always sufficient for local elements.
166  // After we set markers for all local elements, we will do a global
167  // communication to sync markers for ghosted elements from their owners.
168  all_elems = std::make_unique<ConstElemRange>(
169  _fe_problem.mesh().getMesh().active_local_elements_begin(),
170  _fe_problem.mesh().getMesh().active_local_elements_end());
171  }
172  else // This is not scalable but it might be useful for small-size problems
173  {
175  _fe_problem.getAuxiliarySystem().solution().localize(serialized_solution);
176  distributed_adaptivity = false;
177 
178  // For a replicated mesh or a serialized distributed mesh, the solution
179  // is serialized to everyone. Then we update markers for all active elements.
180  // In this case, we can avoid a global communication to update mesh.
181  // I do not know if it is a good idea, but it the old code behavior.
182  // We might not care about much since a replicated mesh
183  // or a serialized distributed mesh is not scalable anyway.
184  all_elems =
185  std::make_unique<ConstElemRange>(_fe_problem.mesh().getMesh().active_elements_begin(),
186  _fe_problem.mesh().getMesh().active_elements_end());
187  }
188 
189  FlagElementsThread fet(
190  _fe_problem, serialized_solution, _max_h_level, marker_name, !distributed_adaptivity);
191  Threads::parallel_reduce(*all_elems, fet);
193  }
194  }
195  else
196  {
197  // Compute the error for each active element
198  _error_estimator->estimate_error(_fe_problem.getNonlinearSystemBase(/*nl_sys=*/0).system(),
199  *_error);
200 
201  // Flag elements to be refined and coarsened
202  _mesh_refinement->flag_elements_by_error_fraction(*_error);
203 
204  if (_displaced_problem)
205  // Reuse the error vector and refine the displaced mesh
206  _displaced_mesh_refinement->flag_elements_by_error_fraction(*_error);
207  }
208 
209  // If the DisplacedProblem is active, undisplace the DisplacedMesh
210  // in preparation for refinement. We can't safely refine the
211  // DisplacedMesh directly, since the Hilbert keys computed on the
212  // inconsistenly-displaced Mesh are different on different
213  // processors, leading to inconsistent Hilbert keys. We must do
214  // this before the undisplaced Mesh is refined, so that the
215  // element and node numbering is still consistent.
216  if (_displaced_problem)
217  _displaced_problem->undisplaceMesh();
218 
219  // If markers are added to only local elements,
220  // we sync them here.
221  if (distributed_adaptivity)
222  _mesh_refinement->make_flags_parallel_consistent();
223 
224  if (_p_refinement_flag)
225  _mesh_refinement->switch_h_to_p_refinement();
226 
227  // Perform refinement and coarsening
228  mesh_changed = _mesh_refinement->refine_and_coarsen_elements();
229 
230  if (_displaced_problem && mesh_changed)
231  {
232  // If markers are added to only local elements,
233  // we sync them here.
234  if (distributed_adaptivity)
235  _displaced_mesh_refinement->make_flags_parallel_consistent();
236 
237  if (_p_refinement_flag)
238  _displaced_mesh_refinement->switch_h_to_p_refinement();
239 
240 #ifndef NDEBUG
241  bool displaced_mesh_changed =
242 #endif
243  _displaced_mesh_refinement->refine_and_coarsen_elements();
244 
245  // Since the undisplaced mesh changed, the displaced mesh better have changed!
246  mooseAssert(displaced_mesh_changed, "Undisplaced mesh changed, but displaced mesh did not!");
247  }
248 
249  if (mesh_changed && _print_mesh_changed)
250  {
251  _console << "\nMesh Changed:\n";
252  _mesh.printInfo();
253  _console << std::flush;
254  }
255 
256  return mesh_changed;
257 }
virtual void update(bool update_libmesh_system=true)
Update the system (doing libMesh magic)
Definition: SystemBase.C:1211
std::shared_ptr< DisplacedProblem > _displaced_problem
Definition: Adaptivity.h:271
std::string _marker_variable_name
Name of the marker variable if using the new adaptivity system.
Definition: Adaptivity.h:301
std::unique_ptr< MeshRefinement > _displaced_mesh_refinement
A mesh refinement object for displaced mesh.
Definition: Adaptivity.h:274
std::unique_ptr< ErrorEstimator > _error_estimator
Error estimator to be used by the apps.
Definition: Adaptivity.h:267
bool _p_refinement_flag
Definition: Adaptivity.h:315
NumericVector< Number > & solution()
Definition: SystemBase.h:176
std::unique_ptr< ErrorVector > _error
Error vector for use with the error estimator.
Definition: Adaptivity.h:269
MooseMesh & _mesh
Definition: Adaptivity.h:258
bool _print_mesh_changed
True if we want to print out info when mesh has changed.
Definition: Adaptivity.h:282
std::unique_ptr< MeshRefinement > _mesh_refinement
A mesh refinement object to be used either with initial refinement or with Adaptivity.
Definition: Adaptivity.h:265
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3198
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
AuxiliarySystem & getAuxiliarySystem()
virtual void close()=0
void printInfo(std::ostream &os=libMesh::out, const unsigned int verbosity=0) const
Calls print_info() on the underlying Mesh.
Definition: MooseMesh.C:3212
bool _use_new_system
Whether or not to use the "new" adaptivity system.
Definition: Adaptivity.h:298
virtual System & system() override
Get the reference to the libMesh system.
virtual MooseMesh & mesh() override
const ConsoleStream _console
An instance of helper class to write streams to the Console objects.
unsigned int _max_h_level
The maximum number of refinement levels.
Definition: Adaptivity.h:307
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257
virtual void localize(std::vector< Number > &v_local) const =0

◆ doingPRefinement()

void Adaptivity::doingPRefinement ( bool  doing_p_refinement,
const MultiMooseEnum disable_p_refinement_for_families 
)

Indicate whether the kind of adaptivity we're doing is p-refinement.

Parameters
doing_p_refinementWhether we're doing p-refinement
disable_p_refinement_for_familiesFamilies to disable p-refinement for

Definition at line 390 of file Adaptivity.C.

Referenced by AdaptivityAction::act(), and SetAdaptivityOptionsAction::act().

392 {
393  _p_refinement_flag = doing_p_refinement;
394  _fe_problem.doingPRefinement(doing_p_refinement, disable_p_refinement_for_families);
395 }
bool _p_refinement_flag
Definition: Adaptivity.h:315
virtual void doingPRefinement(bool doing_p_refinement, const MultiMooseEnum &disable_p_refinement_for_families) override
Indicate whether the kind of adaptivity we&#39;re doing is p-refinement.
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257

◆ getCyclesPerStep()

unsigned int Adaptivity::getCyclesPerStep ( ) const
inline

Pull out the number of cycles_per_step previously set through the AdaptivityAction.

Returns
the number of cycles per step

Definition at line 111 of file Adaptivity.h.

Referenced by FEProblemBase::adaptMesh().

111 { return _cycles_per_step; }
unsigned int _cycles_per_step
The number of adaptivity cycles per step.
Definition: Adaptivity.h:295

◆ getErrorVector()

ErrorVector & Adaptivity::getErrorVector ( const std::string &  indicator_field)

Get an ErrorVector that will be filled up with values corresponding to the indicator field name passed in.

Note that this returns a reference... and the return value should be stored as a reference!

Parameters
indicator_fieldThe name of the field to get an ErrorVector for.

Definition at line 354 of file Adaptivity.C.

Referenced by Marker::getErrorVector().

355 {
356  // Insert or retrieve error vector
357  auto insert_pair = moose_try_emplace(
358  _indicator_field_to_error_vector, indicator_field, std::make_unique<ErrorVector>());
359  return *insert_pair.first->second;
360 }
std::pair< typename M::iterator, bool > moose_try_emplace(M &m, const typename M::key_type &k, Args &&... args)
Function to mirror the behavior of the C++17 std::map::try_emplace() method (no hint).
Definition: Moose.h:94
std::map< std::string, std::unique_ptr< ErrorVector > > _indicator_field_to_error_vector
Stores pointers to ErrorVectors associated with indicator field names.
Definition: Adaptivity.h:313

◆ getInitialSteps()

unsigned int Adaptivity::getInitialSteps ( ) const
inline

Pull out the number of initial steps previously set by calling init()

Returns
the number of initial steps

Definition at line 97 of file Adaptivity.h.

Referenced by FEProblemBase::hasInitialAdaptivity(), FEProblemBase::initialAdaptMesh(), and FEProblemBase::initialSetup().

97 { return _initial_steps; }
unsigned int _initial_steps
the number of adaptivity steps to do at the beginning of simulation
Definition: Adaptivity.h:277

◆ getMaxHLevel()

unsigned int Adaptivity::getMaxHLevel ( )
inline

Return the maximum h-level.

Definition at line 229 of file Adaptivity.h.

229 { return _max_h_level; }
unsigned int _max_h_level
The maximum number of refinement levels.
Definition: Adaptivity.h:307

◆ getRecomputeMarkersFlag()

bool Adaptivity::getRecomputeMarkersFlag ( ) const
inline

Pull out the _recompute_markers_during_cycles flag previously set through the AdaptivityAction.

Returns
the flag to recompute markers during adaptivity cycles

Definition at line 124 of file Adaptivity.h.

Referenced by FEProblemBase::adaptMesh().

bool _recompute_markers_during_cycles
Whether or not to recompute markers during adaptivity cycles.
Definition: Adaptivity.h:310

◆ getSteps()

unsigned int Adaptivity::getSteps ( ) const
inline

Pull out the number of steps previously set by calling init()

Returns
the number of steps

Definition at line 104 of file Adaptivity.h.

Referenced by Steady::execute(), and Eigenvalue::execute().

104 { return _steps; }
unsigned int _steps
steps of adaptivity to perform
Definition: Adaptivity.h:279

◆ init()

void Adaptivity::init ( unsigned int  steps,
unsigned int  initial_steps 
)

Initialize and turn on adaptivity for the simulation.

initial_steps specifies the number of adaptivity cycles to perform before the simulation starts and steps indicates the number of adaptivity cycles to run during a steady (not transient) solve. steps is not used for transient or eigen solves.

Definition at line 57 of file Adaptivity.C.

Referenced by AdaptivityAction::act(), and SetAdaptivityOptionsAction::act().

58 {
59  // Get the pointer to the DisplacedProblem, this cannot be done at construction because
60  // DisplacedProblem
61  // does not exist at that point.
63 
64  _mesh_refinement = std::make_unique<MeshRefinement>(_mesh);
65  _error = std::make_unique<ErrorVector>();
66 
67  EquationSystems & es = _fe_problem.es();
68  es.parameters.set<bool>("adaptivity") = true;
69 
70  _initial_steps = initial_steps;
71  _steps = steps;
72  _mesh_refinement_on = true;
73 
74  _mesh_refinement->set_periodic_boundaries_ptr(
75  _fe_problem.getNonlinearSystemBase(/*nl_sys=*/0).dofMap().get_periodic_boundaries());
76 
77  // displaced problem
78  if (_displaced_problem != nullptr)
79  {
80  EquationSystems & displaced_es = _displaced_problem->es();
81  displaced_es.parameters.set<bool>("adaptivity") = true;
82 
84  _displaced_mesh_refinement = std::make_unique<MeshRefinement>(_displaced_problem->mesh());
85 
86  // The periodic boundaries pointer allows the MeshRefinement
87  // object to determine elements which are "topological" neighbors,
88  // i.e. neighbors across periodic boundaries, for the purposes of
89  // refinement.
90  _displaced_mesh_refinement->set_periodic_boundaries_ptr(
91  _fe_problem.getNonlinearSystemBase(/*nl_sys=*/0).dofMap().get_periodic_boundaries());
92 
93  // TODO: This is currently an empty function on the DisplacedProblem... could it be removed?
94  _displaced_problem->initAdaptivity();
95  }
96 
97  // indicate the Adaptivity system has been initialized
98  _initialized = true;
99 }
std::shared_ptr< DisplacedProblem > _displaced_problem
Definition: Adaptivity.h:271
std::unique_ptr< MeshRefinement > _displaced_mesh_refinement
A mesh refinement object for displaced mesh.
Definition: Adaptivity.h:274
bool _initialized
on/off flag reporting if the adaptivity system has been initialized
Definition: Adaptivity.h:263
std::unique_ptr< ErrorVector > _error
Error vector for use with the error estimator.
Definition: Adaptivity.h:269
MooseMesh & _mesh
Definition: Adaptivity.h:258
virtual EquationSystems & es() override
virtual DofMap & dofMap()
Gets writeable reference to the dof map.
Definition: SystemBase.C:1131
std::unique_ptr< MeshRefinement > _mesh_refinement
A mesh refinement object to be used either with initial refinement or with Adaptivity.
Definition: Adaptivity.h:265
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Definition: Adaptivity.h:261
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
unsigned int _steps
steps of adaptivity to perform
Definition: Adaptivity.h:279
unsigned int _initial_steps
the number of adaptivity steps to do at the beginning of simulation
Definition: Adaptivity.h:277
virtual std::shared_ptr< const DisplacedProblem > getDisplacedProblem() const
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257

◆ initialAdaptMesh()

bool Adaptivity::initialAdaptMesh ( )

Used during initial adaptivity.

Returns
a boolean that indicates whether the mesh was changed

Definition at line 260 of file Adaptivity.C.

Referenced by FEProblemBase::initialAdaptMesh().

261 {
263 }
bool adaptMesh(std::string marker_name=std::string())
Adapts the mesh based on the error estimator used.
Definition: Adaptivity.C:123
std::string _initial_marker_variable_name
Name of the initial marker variable if using the new adaptivity system.
Definition: Adaptivity.h:304

◆ isAdaptivityDue()

bool Adaptivity::isAdaptivityDue ( )

Query if an adaptivity step should be performed at the current time / time step.

Definition at line 384 of file Adaptivity.C.

Referenced by FEProblemBase::adaptMesh().

385 {
386  return _mesh_refinement_on && (_start_time <= _t && _t < _stop_time) && _step % _interval == 0;
387 }
Real _stop_time
When adaptivity stops.
Definition: Adaptivity.h:293
Real & _t
Time.
Definition: Adaptivity.h:285
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Definition: Adaptivity.h:261
unsigned int _interval
intreval between adaptivity runs
Definition: Adaptivity.h:289
Real _start_time
When adaptivity start.
Definition: Adaptivity.h:291
int & _step
Time Step.
Definition: Adaptivity.h:287

◆ isInitialized()

bool Adaptivity::isInitialized ( )
inline

Returns whether or not Adaptivity::init() has ran.

Can be used to indicate if mesh adaptivity is available.

Returns
true if the Adaptivity system is ready to be used, otherwise false

Definition at line 187 of file Adaptivity.h.

187 { return _initialized; }
bool _initialized
on/off flag reporting if the adaptivity system has been initialized
Definition: Adaptivity.h:263

◆ isOn()

bool Adaptivity::isOn ( )
inline

Is adaptivity on?

Returns
true if mesh adaptivity is on, otherwise false

Definition at line 179 of file Adaptivity.h.

Referenced by FEProblemBase::addAnyRedistributers(), FEProblemBase::checkProblemIntegrity(), and FEProblemBase::initialSetup().

179 { return _mesh_refinement_on; }
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Definition: Adaptivity.h:261

◆ perfGraph()

PerfGraph & PerfGraphInterface::perfGraph ( )
inherited

Get the PerfGraph.

Definition at line 78 of file PerfGraphInterface.C.

Referenced by CommonOutputAction::act(), PerfGraphData::finalize(), and PerfGraphOutput::output().

79 {
80  return _pg_moose_app.perfGraph();
81 }
MooseApp & _pg_moose_app
The MooseApp that owns the PerfGraph.
PerfGraph & perfGraph()
Get the PerfGraph for this app.
Definition: MooseApp.h:133

◆ registerTimedSection() [1/2]

PerfID PerfGraphInterface::registerTimedSection ( const std::string &  section_name,
const unsigned int  level 
) const
protectedinherited

Call to register a named section for timing.

Parameters
section_nameThe name of the code section to be timed
levelThe importance of the timer - lower is more important (0 will always come out)
Returns
The ID of the section - use when starting timing

Definition at line 53 of file PerfGraphInterface.C.

55 {
56  const auto timed_section_name = timedSectionName(section_name);
57  if (!moose::internal::getPerfGraphRegistry().sectionExists(timed_section_name))
58  return moose::internal::getPerfGraphRegistry().registerSection(timed_section_name, level);
59  else
60  return moose::internal::getPerfGraphRegistry().sectionID(timed_section_name);
61 }
PerfID registerSection(const std::string &section_name, const unsigned int level)
Call to register a named section for timing.
std::string timedSectionName(const std::string &section_name) const
PerfID sectionID(const std::string &section_name) const
Given a name return the PerfID The name of the section.
PerfGraphRegistry & getPerfGraphRegistry()
Get the global PerfGraphRegistry singleton.

◆ registerTimedSection() [2/2]

PerfID PerfGraphInterface::registerTimedSection ( const std::string &  section_name,
const unsigned int  level,
const std::string &  live_message,
const bool  print_dots = true 
) const
protectedinherited

Call to register a named section for timing.

Parameters
section_nameThe name of the code section to be timed
levelThe importance of the timer - lower is more important (0 will always come out)
live_messageThe message to be printed to the screen during execution
print_dotsWhether or not progress dots should be printed for this section
Returns
The ID of the section - use when starting timing

Definition at line 64 of file PerfGraphInterface.C.

68 {
69  const auto timed_section_name = timedSectionName(section_name);
70  if (!moose::internal::getPerfGraphRegistry().sectionExists(timed_section_name))
72  timedSectionName(section_name), level, live_message, print_dots);
73  else
74  return moose::internal::getPerfGraphRegistry().sectionID(timed_section_name);
75 }
PerfID registerSection(const std::string &section_name, const unsigned int level)
Call to register a named section for timing.
std::string timedSectionName(const std::string &section_name) const
PerfID sectionID(const std::string &section_name) const
Given a name return the PerfID The name of the section.
PerfGraphRegistry & getPerfGraphRegistry()
Get the global PerfGraphRegistry singleton.

◆ setAdaptivityOn()

void Adaptivity::setAdaptivityOn ( bool  state)

Allow adaptivity to be toggled programatically.

Parameters
stateThe adaptivity state (on/off).

Definition at line 319 of file Adaptivity.C.

320 {
321  // check if Adaptivity has been initialized before turning on
322  if (state == true && !_initialized)
323  mooseError("Mesh adaptivity system not available");
324 
325  _mesh_refinement_on = state;
326 }
bool _initialized
on/off flag reporting if the adaptivity system has been initialized
Definition: Adaptivity.h:263
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:284
bool _mesh_refinement_on
on/off flag reporting if the adaptivity is being used
Definition: Adaptivity.h:261

◆ setCyclesPerStep()

void Adaptivity::setCyclesPerStep ( const unsigned int num)
inline

Set the number of cycles_per_step.

Parameters
numThe number of cycles per step to execute

Definition at line 117 of file Adaptivity.h.

Referenced by SetAdaptivityOptionsAction::act().

117 { _cycles_per_step = num; }
unsigned int _cycles_per_step
The number of adaptivity cycles per step.
Definition: Adaptivity.h:295

◆ setErrorEstimator()

void Adaptivity::setErrorEstimator ( const MooseEnum error_estimator_name)

Set the error estimator.

Parameters
error_estimator_namethe name of the error estimator (currently: Laplacian, Kelly, and PatchRecovery)

Definition at line 102 of file Adaptivity.C.

Referenced by AdaptivityAction::act().

103 {
104  if (error_estimator_name == "KellyErrorEstimator")
105  _error_estimator = std::make_unique<KellyErrorEstimator>();
106  else if (error_estimator_name == "LaplacianErrorEstimator")
107  _error_estimator = std::make_unique<LaplacianErrorEstimator>();
108  else if (error_estimator_name == "PatchRecoveryErrorEstimator")
109  _error_estimator = std::make_unique<PatchRecoveryErrorEstimator>();
110  else
111  mooseError(std::string("Unknown error_estimator selection: ") +
112  std::string(error_estimator_name));
113 }
std::unique_ptr< ErrorEstimator > _error_estimator
Error estimator to be used by the apps.
Definition: Adaptivity.h:267
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:284

◆ setErrorNorm()

void Adaptivity::setErrorNorm ( SystemNorm &  sys_norm)

Set the error norm (FIXME: improve description)

Definition at line 116 of file Adaptivity.C.

Referenced by AdaptivityAction::act().

117 {
118  mooseAssert(_error_estimator, "error_estimator not initialized. Did you call init_adaptivity()?");
119  _error_estimator->error_norm = sys_norm;
120 }
std::unique_ptr< ErrorEstimator > _error_estimator
Error estimator to be used by the apps.
Definition: Adaptivity.h:267

◆ setInitialMarkerVariableName()

void Adaptivity::setInitialMarkerVariableName ( std::string  marker_field)

Sets the name of the field variable to actually use to flag elements for initial refinement / coarsening.

This must be a CONSTANT, MONOMIAL Auxiliary Variable Name that contains values corresponding to libMesh::Elem::RefinementState.

Parameters
marker_fieldThe name of the field to use for refinement / coarsening.

Definition at line 348 of file Adaptivity.C.

Referenced by SetAdaptivityOptionsAction::act().

349 {
350  _initial_marker_variable_name = marker_field;
351 }
std::string _initial_marker_variable_name
Name of the initial marker variable if using the new adaptivity system.
Definition: Adaptivity.h:304

◆ setInterval()

void Adaptivity::setInterval ( unsigned int  interval)
inline

Set the interval (number of timesteps) between refinement steps.

Definition at line 234 of file Adaptivity.h.

Referenced by AdaptivityAction::act(), and SetAdaptivityOptionsAction::act().

234 { _interval = interval; }
unsigned int _interval
intreval between adaptivity runs
Definition: Adaptivity.h:289

◆ setMarkerVariableName()

void Adaptivity::setMarkerVariableName ( std::string  marker_field)

Sets the name of the field variable to actually use to flag elements for refinement / coarsening.

This must be a CONSTANT, MONOMIAL Auxiliary Variable Name that contains values corresponding to libMesh::Elem::RefinementState.

Parameters
marker_fieldThe name of the field to use for refinement / coarsening.

Definition at line 342 of file Adaptivity.C.

Referenced by SetAdaptivityOptionsAction::act().

343 {
344  _marker_variable_name = marker_field;
345 }
std::string _marker_variable_name
Name of the marker variable if using the new adaptivity system.
Definition: Adaptivity.h:301

◆ setMaxHLevel()

void Adaptivity::setMaxHLevel ( unsigned int  level)
inline

Set the maximum refinement level (for the new Adaptivity system).

Definition at line 224 of file Adaptivity.h.

Referenced by SetAdaptivityOptionsAction::act().

224 { _max_h_level = level; }
unsigned int _max_h_level
The maximum number of refinement levels.
Definition: Adaptivity.h:307

◆ setParam()

template<typename T >
void Adaptivity::setParam ( const std::string &  param_name,
const T &  param_value 
)

Set adaptivity parameter.

Parameters
param_namethe name of the parameter
param_valuethe value of parameter

Definition at line 320 of file Adaptivity.h.

Referenced by AdaptivityAction::act().

321 {
322  if (param_name == "refine fraction")
323  {
324  _mesh_refinement->refine_fraction() = param_value;
326  _displaced_mesh_refinement->refine_fraction() = param_value;
327  }
328  else if (param_name == "coarsen fraction")
329  {
330  _mesh_refinement->coarsen_fraction() = param_value;
332  _displaced_mesh_refinement->coarsen_fraction() = param_value;
333  }
334  else if (param_name == "max h-level")
335  {
336  _mesh_refinement->max_h_level() = param_value;
338  _displaced_mesh_refinement->max_h_level() = param_value;
339  }
340  else if (param_name == "cycles_per_step")
341  _cycles_per_step = param_value;
342  else if (param_name == "recompute_markers_during_cycles")
343  _recompute_markers_during_cycles = param_value;
344  else
345  mooseError("Invalid Param in adaptivity object");
346 }
std::unique_ptr< MeshRefinement > _displaced_mesh_refinement
A mesh refinement object for displaced mesh.
Definition: Adaptivity.h:274
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:284
bool _recompute_markers_during_cycles
Whether or not to recompute markers during adaptivity cycles.
Definition: Adaptivity.h:310
std::unique_ptr< MeshRefinement > _mesh_refinement
A mesh refinement object to be used either with initial refinement or with Adaptivity.
Definition: Adaptivity.h:265
unsigned int _cycles_per_step
The number of adaptivity cycles per step.
Definition: Adaptivity.h:295

◆ setPrintMeshChanged()

void Adaptivity::setPrintMeshChanged ( bool  state = true)
inline

Definition at line 90 of file Adaptivity.h.

Referenced by AdaptivityAction::act().

90 { _print_mesh_changed = state; }
bool _print_mesh_changed
True if we want to print out info when mesh has changed.
Definition: Adaptivity.h:282

◆ setRecomputeMarkersFlag()

void Adaptivity::setRecomputeMarkersFlag ( const bool  flag)
inline

Set the flag to recompute markers during adaptivity cycles.

Parameters
flagThe flag to recompute markers

Definition at line 130 of file Adaptivity.h.

Referenced by SetAdaptivityOptionsAction::act().

bool _recompute_markers_during_cycles
Whether or not to recompute markers during adaptivity cycles.
Definition: Adaptivity.h:310

◆ setTimeActive()

void Adaptivity::setTimeActive ( Real  start_time,
Real  stop_time 
)

Sets the time when the adaptivity is active.

Parameters
start_timeThe time when adaptivity starts
stop_timeThe time when adaptivity stops

Definition at line 329 of file Adaptivity.C.

Referenced by AdaptivityAction::act(), and SetAdaptivityOptionsAction::act().

330 {
331  _start_time = start_time;
332  _stop_time = stop_time;
333 }
Real _stop_time
When adaptivity stops.
Definition: Adaptivity.h:293
Real _start_time
When adaptivity start.
Definition: Adaptivity.h:291

◆ setUseNewSystem()

void Adaptivity::setUseNewSystem ( )

Tells this object we're using the "new" adaptivity system.

Definition at line 336 of file Adaptivity.C.

Referenced by SetAdaptivityOptionsAction::act().

337 {
338  _use_new_system = true;
339 }
bool _use_new_system
Whether or not to use the "new" adaptivity system.
Definition: Adaptivity.h:298

◆ timedSectionName()

std::string PerfGraphInterface::timedSectionName ( const std::string &  section_name) const
protectedinherited
Returns
The name of the timed section with the name section_name.

Optionally adds a prefix if one is defined.

Definition at line 47 of file PerfGraphInterface.C.

Referenced by PerfGraphInterface::registerTimedSection().

48 {
49  return _prefix.empty() ? "" : (_prefix + "::") + section_name;
50 }
const std::string _prefix
A prefix to use for all sections.

◆ uniformRefine()

void Adaptivity::uniformRefine ( MooseMesh mesh,
unsigned int  level = libMesh::invalid_uint 
)
static

Performs uniform refinement of the passed Mesh object.

The number of levels of refinement performed is stored in the MooseMesh object. No solution projection is performed in this version.

Definition at line 266 of file Adaptivity.C.

Referenced by SetupMeshCompleteAction::act(), and FEProblemBase::uniformRefine().

267 {
268  mooseAssert(mesh, "Mesh pointer must not be NULL");
269 
270  // NOTE: we are using a separate object here, since adaptivity may not be on, but we need to be
271  // able to do refinements
272  MeshRefinement mesh_refinement(*mesh);
273  if (level == libMesh::invalid_uint)
274  level = mesh->uniformRefineLevel();
275 
276  // Skip deletion and repartition will make uniform refinements run more
277  // efficiently, but at the same time, there might be extra ghosting elements.
278  // The number of layers of additional ghosting elements depends on the number
279  // of uniform refinement levels. This should happen only when you have a "fine enough"
280  // coarse mesh and want to refine the mesh by a few levels. Otherwise, it might
281  // introduce an unbalanced workload and too large ghosting domain.
282  if (mesh->skipDeletionRepartitionAfterRefine())
283  {
284  mesh->getMesh().skip_partitioning(true);
285  mesh->getMesh().allow_remote_element_removal(false);
286  mesh->needsRemoteElemDeletion(false);
287  }
288 
289  mesh_refinement.uniformly_refine(level);
290 }
const unsigned int invalid_uint
MeshBase & mesh

◆ uniformRefineWithProjection()

void Adaptivity::uniformRefineWithProjection ( )

Performs uniform refinement on the meshes in the current object.

Projections are performed of the solution vectors.

Definition at line 293 of file Adaptivity.C.

Referenced by FEProblemBase::initialSetup().

294 {
295  TIME_SECTION("uniformRefineWithProjection", 2, "Uniformly Refining and Reprojecting");
296 
297  // NOTE: we are using a separate object here, since adaptivity may not be on, but we need to be
298  // able to do refinements
299  MeshRefinement mesh_refinement(_mesh);
300  unsigned int level = _mesh.uniformRefineLevel();
301  MeshRefinement displaced_mesh_refinement(_displaced_problem ? _displaced_problem->mesh() : _mesh);
302 
303  // we have to go step by step so EquationSystems::reinit() won't freak out
304  for (unsigned int i = 0; i < level; i++)
305  {
306  // See comment above about why refining the displaced mesh is potentially unsafe.
307  if (_displaced_problem)
308  _displaced_problem->undisplaceMesh();
309 
310  mesh_refinement.uniformly_refine(1);
311 
312  if (_displaced_problem)
313  displaced_mesh_refinement.uniformly_refine(1);
315  }
316 }
std::shared_ptr< DisplacedProblem > _displaced_problem
Definition: Adaptivity.h:271
MooseMesh & _mesh
Definition: Adaptivity.h:258
virtual void meshChanged() override
Update data after a mesh change.
unsigned int uniformRefineLevel() const
Returns the level of uniform refinement requested (zero if AMR is disabled).
Definition: MooseMesh.C:2966
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257

◆ updateErrorVectors()

void Adaptivity::updateErrorVectors ( )

Update the ErrorVectors that have been requested through calls to getErrorVector().

Definition at line 363 of file Adaptivity.C.

Referenced by FEProblemBase::computeMarkers().

364 {
365  TIME_SECTION("updateErrorVectors", 5, "Updating Error Vectors");
366 
367  // Resize all of the ErrorVectors in case the mesh has changed
368  for (const auto & it : _indicator_field_to_error_vector)
369  {
370  ErrorVector & vec = *(it.second);
371  vec.assign(_mesh.getMesh().max_elem_id(), 0);
372  }
373 
374  // Fill the vectors with the local contributions
376  Threads::parallel_reduce(*_mesh.getActiveLocalElementRange(), uevt);
377 
378  // Now sum across all processors
379  for (const auto & it : _indicator_field_to_error_vector)
380  _fe_problem.comm().sum((std::vector<float> &)*(it.second));
381 }
ConstElemRange * getActiveLocalElementRange()
Return pointers to range objects for various types of ranges (local nodes, boundary elems...
Definition: MooseMesh.C:1040
const Parallel::Communicator & comm() const
MooseMesh & _mesh
Definition: Adaptivity.h:258
std::map< std::string, std::unique_ptr< ErrorVector > > _indicator_field_to_error_vector
Stores pointers to ErrorVectors associated with indicator field names.
Definition: Adaptivity.h:313
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3198
FEProblemBase & _fe_problem
Definition: Adaptivity.h:257

◆ validParams()

InputParameters PerfGraphInterface::validParams ( )
staticinherited

Definition at line 16 of file PerfGraphInterface.C.

17 {
19  return params;
20 }
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
InputParameters emptyInputParameters()

Member Data Documentation

◆ _console

const ConsoleStream ConsoleStreamInterface::_console
inherited

An instance of helper class to write streams to the Console objects.

Definition at line 31 of file ConsoleStreamInterface.h.

Referenced by IterationAdaptiveDT::acceptStep(), MeshOnlyAction::act(), SetupDebugAction::act(), MaterialOutputAction::act(), adaptMesh(), FEProblemBase::adaptMesh(), PerfGraph::addToExecutionList(), SimplePredictor::apply(), SystemBase::applyScalingFactors(), MultiApp::backup(), FEProblemBase::backupMultiApps(), CoarsenedPiecewiseLinear::buildCoarsenedGrid(), MeshDiagnosticsGenerator::checkElementOverlap(), MeshDiagnosticsGenerator::checkElementTypes(), MeshDiagnosticsGenerator::checkElementVolumes(), FEProblemBase::checkExceptionAndStopSolve(), MeshDiagnosticsGenerator::checkLocalJacobians(), MeshDiagnosticsGenerator::checkNonConformalMesh(), MeshDiagnosticsGenerator::checkNonConformalMeshFromAdaptivity(), MeshDiagnosticsGenerator::checkNonPlanarSides(), FEProblemBase::checkProblemIntegrity(), ReferenceResidualProblem::checkRelativeConvergence(), MeshDiagnosticsGenerator::checkSidesetsOrientation(), IterationAdaptiveDT::computeAdaptiveDT(), Transient::computeConstrainedDT(), FixedPointSolve::computeCustomConvergencePostprocessor(), NonlinearSystemBase::computeDamping(), IterationAdaptiveDT::computeDT(), IterationAdaptiveDT::computeFailedDT(), IterationAdaptiveDT::computeInitialDT(), IterationAdaptiveDT::computeInterpolationDT(), NonlinearSystemBase::computeScaling(), Problem::console(), IterationAdaptiveDT::constrainStep(), TimeStepper::constrainStep(), MultiApp::createApp(), FEProblemBase::execMultiApps(), FEProblemBase::execMultiAppTransfers(), MessageFromInput::execute(), Steady::execute(), Eigenvalue::execute(), ActionWarehouse::executeActionsWithAction(), ActionWarehouse::executeAllActions(), ElementQualityChecker::finalize(), FEProblemBase::finishMultiAppStep(), MeshRepairGenerator::fixOverlappingNodes(), CoarsenBlockGenerator::generate(), MeshGenerator::generateInternal(), VariableCondensationPreconditioner::getDofToCondense(), InversePowerMethod::init(), NonlinearEigen::init(), FEProblemBase::initialAdaptMesh(), EigenExecutionerBase::inversePowerIteration(), FEProblemBase::joinAndFinalize(), Transient::keepGoing(), IterationAdaptiveDT::limitDTByFunction(), IterationAdaptiveDT::limitDTToPostprocessorValue(), FEProblemBase::logAdd(), EigenExecutionerBase::makeBXConsistent(), Console::meshChanged(), MooseBaseErrorInterface::mooseDeprecated(), MooseBaseErrorInterface::mooseInfo(), MooseBaseErrorInterface::mooseWarning(), MooseBaseErrorInterface::mooseWarningNonPrefixed(), ReferenceResidualProblem::nonlinearConvergenceSetup(), ReporterDebugOutput::output(), PerfGraphOutput::output(), MaterialPropertyDebugOutput::output(), DOFMapOutput::output(), VariableResidualNormsDebugOutput::output(), Console::output(), ControlOutput::outputActiveObjects(), ControlOutput::outputChangedControls(), ControlOutput::outputControls(), Console::outputInput(), Console::outputPostprocessors(), PseudoTimestep::outputPseudoTimestep(), Console::outputReporters(), Console::outputScalarVariables(), Console::outputSystemInformation(), FEProblemBase::possiblyRebuildGeomSearchPatches(), EigenExecutionerBase::postExecute(), AB2PredictorCorrector::postSolve(), ActionWarehouse::printActionDependencySets(), SolutionInvalidity::printDebug(), EigenExecutionerBase::printEigenvalue(), SecantSolve::printFixedPointConvergenceHistory(), SteffensenSolve::printFixedPointConvergenceHistory(), PicardSolve::printFixedPointConvergenceHistory(), FixedPointSolve::printFixedPointConvergenceReason(), PerfGraphLivePrint::printLiveMessage(), MaterialPropertyDebugOutput::printMaterialMap(), PerfGraphLivePrint::printStats(), AutomaticMortarGeneration::projectPrimaryNodesSinglePair(), AutomaticMortarGeneration::projectSecondaryNodesSinglePair(), CoarsenBlockGenerator::recursiveCoarsen(), SolutionTimeAdaptiveDT::rejectStep(), MultiApp::restore(), FEProblemBase::restoreMultiApps(), SimplePredictor::shouldApply(), Checkpoint::shouldOutput(), SubProblem::showFunctorRequestors(), SubProblem::showFunctors(), FullSolveMultiApp::showStatusMessage(), FEProblemSolve::solve(), FixedPointSolve::solve(), NonlinearSystem::solve(), EigenProblem::solve(), LStableDirk2::solve(), LStableDirk3::solve(), ImplicitMidpoint::solve(), ExplicitTVDRK2::solve(), LStableDirk4::solve(), AStableDirk4::solve(), ExplicitRK2::solve(), TransientMultiApp::solveStep(), FixedPointSolve::solveStep(), PerfGraphLivePrint::start(), AB2PredictorCorrector::step(), NonlinearEigen::takeStep(), Transient::takeStep(), Console::writeTimestepInformation(), Console::writeVariableNorms(), and FEProblemBase::~FEProblemBase().

◆ _cycles_per_step

unsigned int Adaptivity::_cycles_per_step
protected

The number of adaptivity cycles per step.

Definition at line 295 of file Adaptivity.h.

Referenced by getCyclesPerStep(), setCyclesPerStep(), and setParam().

◆ _displaced_mesh_refinement

std::unique_ptr<MeshRefinement> Adaptivity::_displaced_mesh_refinement
protected

A mesh refinement object for displaced mesh.

Definition at line 274 of file Adaptivity.h.

Referenced by adaptMesh(), init(), and setParam().

◆ _displaced_problem

std::shared_ptr<DisplacedProblem> Adaptivity::_displaced_problem
protected

Definition at line 271 of file Adaptivity.h.

Referenced by adaptMesh(), init(), and uniformRefineWithProjection().

◆ _error

std::unique_ptr<ErrorVector> Adaptivity::_error
protected

Error vector for use with the error estimator.

Definition at line 269 of file Adaptivity.h.

Referenced by adaptMesh(), and init().

◆ _error_estimator

std::unique_ptr<ErrorEstimator> Adaptivity::_error_estimator
protected

Error estimator to be used by the apps.

Definition at line 267 of file Adaptivity.h.

Referenced by adaptMesh(), setErrorEstimator(), and setErrorNorm().

◆ _fe_problem

FEProblemBase& Adaptivity::_fe_problem
protected

◆ _indicator_field_to_error_vector

std::map<std::string, std::unique_ptr<ErrorVector> > Adaptivity::_indicator_field_to_error_vector
protected

Stores pointers to ErrorVectors associated with indicator field names.

Definition at line 313 of file Adaptivity.h.

Referenced by getErrorVector(), and updateErrorVectors().

◆ _initial_marker_variable_name

std::string Adaptivity::_initial_marker_variable_name
protected

Name of the initial marker variable if using the new adaptivity system.

Definition at line 304 of file Adaptivity.h.

Referenced by initialAdaptMesh(), and setInitialMarkerVariableName().

◆ _initial_steps

unsigned int Adaptivity::_initial_steps
protected

the number of adaptivity steps to do at the beginning of simulation

Definition at line 277 of file Adaptivity.h.

Referenced by getInitialSteps(), and init().

◆ _initialized

bool Adaptivity::_initialized
protected

on/off flag reporting if the adaptivity system has been initialized

Definition at line 263 of file Adaptivity.h.

Referenced by init(), isInitialized(), and setAdaptivityOn().

◆ _interval

unsigned int Adaptivity::_interval
protected

intreval between adaptivity runs

Definition at line 289 of file Adaptivity.h.

Referenced by isAdaptivityDue(), and setInterval().

◆ _marker_variable_name

std::string Adaptivity::_marker_variable_name
protected

Name of the marker variable if using the new adaptivity system.

Definition at line 301 of file Adaptivity.h.

Referenced by adaptMesh(), and setMarkerVariableName().

◆ _max_h_level

unsigned int Adaptivity::_max_h_level
protected

The maximum number of refinement levels.

Definition at line 307 of file Adaptivity.h.

Referenced by adaptMesh(), getMaxHLevel(), and setMaxHLevel().

◆ _mesh

MooseMesh& Adaptivity::_mesh
protected

Definition at line 258 of file Adaptivity.h.

Referenced by adaptMesh(), init(), uniformRefineWithProjection(), and updateErrorVectors().

◆ _mesh_refinement

std::unique_ptr<MeshRefinement> Adaptivity::_mesh_refinement
protected

A mesh refinement object to be used either with initial refinement or with Adaptivity.

Definition at line 265 of file Adaptivity.h.

Referenced by adaptMesh(), init(), and setParam().

◆ _mesh_refinement_on

bool Adaptivity::_mesh_refinement_on
protected

on/off flag reporting if the adaptivity is being used

Definition at line 261 of file Adaptivity.h.

Referenced by init(), isAdaptivityDue(), isOn(), and setAdaptivityOn().

◆ _p_refinement_flag

bool Adaptivity::_p_refinement_flag = false
protected

Definition at line 315 of file Adaptivity.h.

Referenced by adaptMesh(), and doingPRefinement().

◆ _pg_moose_app

MooseApp& PerfGraphInterface::_pg_moose_app
protectedinherited

The MooseApp that owns the PerfGraph.

Definition at line 124 of file PerfGraphInterface.h.

Referenced by PerfGraphInterface::perfGraph().

◆ _prefix

const std::string PerfGraphInterface::_prefix
protectedinherited

A prefix to use for all sections.

Definition at line 127 of file PerfGraphInterface.h.

Referenced by PerfGraphInterface::timedSectionName().

◆ _print_mesh_changed

bool Adaptivity::_print_mesh_changed
protected

True if we want to print out info when mesh has changed.

Definition at line 282 of file Adaptivity.h.

Referenced by adaptMesh(), and setPrintMeshChanged().

◆ _recompute_markers_during_cycles

bool Adaptivity::_recompute_markers_during_cycles
protected

Whether or not to recompute markers during adaptivity cycles.

Definition at line 310 of file Adaptivity.h.

Referenced by getRecomputeMarkersFlag(), setParam(), and setRecomputeMarkersFlag().

◆ _start_time

Real Adaptivity::_start_time
protected

When adaptivity start.

Definition at line 291 of file Adaptivity.h.

Referenced by isAdaptivityDue(), and setTimeActive().

◆ _step

int& Adaptivity::_step
protected

Time Step.

Definition at line 287 of file Adaptivity.h.

Referenced by isAdaptivityDue().

◆ _steps

unsigned int Adaptivity::_steps
protected

steps of adaptivity to perform

Definition at line 279 of file Adaptivity.h.

Referenced by getSteps(), and init().

◆ _stop_time

Real Adaptivity::_stop_time
protected

When adaptivity stops.

Definition at line 293 of file Adaptivity.h.

Referenced by isAdaptivityDue(), and setTimeActive().

◆ _t

Real& Adaptivity::_t
protected

Time.

Definition at line 285 of file Adaptivity.h.

Referenced by isAdaptivityDue().

◆ _use_new_system

bool Adaptivity::_use_new_system
protected

Whether or not to use the "new" adaptivity system.

Definition at line 298 of file Adaptivity.h.

Referenced by adaptMesh(), and setUseNewSystem().


The documentation for this class was generated from the following files: