Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #pragma once 11 : 12 : #include "libmesh/libmesh_config.h" 13 : 14 : #ifdef LIBMESH_ENABLE_AMR 15 : 16 : #include "Moose.h" 17 : #include "MooseError.h" 18 : #include "ConsoleStreamInterface.h" 19 : #include "MooseTypes.h" 20 : #include "PerfGraphInterface.h" 21 : 22 : #include "libmesh/parallel_object.h" 23 : 24 : // libMesh 25 : #include "libmesh/mesh_refinement.h" 26 : #include "libmesh/hp_coarsentest.h" 27 : #include "libmesh/sibling_coupling.h" 28 : 29 : class FEProblemBase; 30 : class MooseMesh; 31 : class DisplacedProblem; 32 : template <typename> 33 : class MooseVariableFE; 34 : typedef MooseVariableFE<Real> MooseVariable; 35 : typedef MooseVariableFE<libMesh::VectorValue<Real>> VectorMooseVariable; 36 : class MooseEnum; 37 : class MultiMooseEnum; 38 : 39 : // Forward declare classes in libMesh 40 : namespace libMesh 41 : { 42 : class SystemNorm; 43 : class ErrorVector; 44 : class ErrorEstimator; 45 : } 46 : 47 : /** 48 : * Defines types of mesh adaptivity options available 49 : */ 50 : enum class AdaptivityType 51 : { 52 : H = 0, 53 : P = 1, 54 : HP = 2 55 : }; 56 : 57 : /** 58 : * Takes care of everything related to mesh adaptivity 59 : * 60 : */ 61 : class Adaptivity : public ConsoleStreamInterface, 62 : public PerfGraphInterface, 63 : public libMesh::ParallelObject 64 : { 65 : public: 66 : Adaptivity(FEProblemBase & fe_problem); 67 : virtual ~Adaptivity(); 68 : 69 : /** 70 : * Initialize and turn on adaptivity for the simulation. initial_steps specifies the number of 71 : * adaptivity cycles to perform before the simulation starts and steps indicates the 72 : * number of adaptivity cycles to run during a steady (not transient) solve. steps is not used 73 : * for transient or eigen solves. 74 : * p_refinement indicates whether the refinement will be p-refinement or h-refinement. 75 : */ 76 : void init(const unsigned int steps, 77 : const unsigned int initial_steps, 78 : const AdaptivityType adaptivity_type); 79 : 80 : /** 81 : * Set adaptivity parameter 82 : * 83 : * @param param_name the name of the parameter 84 : * @param param_value the value of parameter 85 : */ 86 : template <typename T> 87 : void setParam(const std::string & param_name, const T & param_value); 88 : 89 : /** 90 : * Set the error estimator 91 : * 92 : * @param error_estimator_name the name of the error estimator (currently: Laplacian, Kelly, and 93 : * PatchRecovery) 94 : */ 95 : void setErrorEstimator(const MooseEnum & error_estimator_name); 96 : 97 : /** 98 : * Set the error norm (FIXME: improve description) 99 : */ 100 : void setErrorNorm(libMesh::SystemNorm & sys_norm); 101 : 102 : /** 103 : * 104 : */ 105 422 : void setPrintMeshChanged(bool state = true) { _print_mesh_changed = state; } 106 : 107 : /** 108 : * Pull out the number of initial steps previously set by calling init() 109 : * 110 : * @return the number of initial steps 111 : */ 112 109844 : unsigned int getInitialSteps() const { return _initial_steps; } 113 : 114 : /** 115 : * Pull out the number of steps previously set by calling init() 116 : * 117 : * @return the number of steps 118 : */ 119 27354 : unsigned int getSteps() const { return _steps; } 120 : 121 : /** 122 : * Pull out the number of cycles_per_step previously set through the AdaptivityAction 123 : * 124 : * @return the number of cycles per step 125 : */ 126 4555 : unsigned int getCyclesPerStep() const { return _cycles_per_step; } 127 : 128 : /** 129 : * Set the number of cycles_per_step 130 : * @param num The number of cycles per step to execute 131 : */ 132 1825 : void setCyclesPerStep(const unsigned int & num) { _cycles_per_step = num; } 133 : 134 : /** 135 : * Pull out the _recompute_markers_during_cycles flag previously set through the AdaptivityAction 136 : * 137 : * @return the flag to recompute markers during adaptivity cycles 138 : */ 139 4707 : bool getRecomputeMarkersFlag() const { return _recompute_markers_during_cycles; } 140 : 141 : /** 142 : * Set the flag to recompute markers during adaptivity cycles 143 : * @param flag The flag to recompute markers 144 : */ 145 1825 : void setRecomputeMarkersFlag(const bool flag) { _recompute_markers_during_cycles = flag; } 146 : 147 : /** 148 : * Adapts the mesh based on the error estimator used 149 : * 150 : * This method calls DisplacedProblem::undisplaceMesh, so it is necessary to 151 : * update the displaced problem's mesh (undoing the undisplacement) before any 152 : * objects relying on the displaced mesh are executed. However, projection of 153 : * the displacement variables is required before doing so; therefore, the 154 : * DisplacedProblem::updateMesh call needs to happen outside of this method. 155 : * @see FEProblemBase::adaptMesh and FEProblemMase::meshChangedHelper. 156 : * 157 : * @return a boolean that indicates whether the mesh was changed 158 : */ 159 : bool adaptMesh(std::string marker_name = std::string()); 160 : 161 : /** 162 : * Used during initial adaptivity. 163 : * 164 : * @return a boolean that indicates whether the mesh was changed 165 : */ 166 : bool initialAdaptMesh(); 167 : 168 : /** 169 : * Performs uniform refinement of the passed Mesh object. The 170 : * number of levels of refinement performed is stored in the 171 : * MooseMesh object. No solution projection is performed in this 172 : * version. 173 : */ 174 : static void uniformRefine(MooseMesh * mesh, unsigned int level = libMesh::invalid_uint); 175 : 176 : /** 177 : * Performs uniform refinement on the meshes in the current 178 : * object. Projections are performed of the solution vectors. 179 : */ 180 : void uniformRefineWithProjection(); 181 : 182 : /** 183 : * Allow adaptivity to be toggled programatically. 184 : * @param state The adaptivity state (on/off). 185 : */ 186 : void setAdaptivityOn(bool state); 187 : 188 : /** 189 : * Is adaptivity on? 190 : * 191 : * @return true if mesh adaptivity is on, otherwise false 192 : */ 193 180908 : bool isOn() { return _mesh_refinement_on; } 194 : 195 : /** 196 : * Returns whether or not Adaptivity::init() has ran. Can 197 : * be used to indicate if mesh adaptivity is available. 198 : * 199 : * @return true if the Adaptivity system is ready to be used, otherwise false 200 : */ 201 : bool isInitialized() { return _initialized; } 202 : 203 : /** 204 : * Sets the time when the adaptivity is active 205 : * @param start_time The time when adaptivity starts 206 : * @param stop_time The time when adaptivity stops 207 : */ 208 : void setTimeActive(Real start_time, Real stop_time); 209 : 210 : /** 211 : * Tells this object we're using the "new" adaptivity system. 212 : */ 213 : void setUseNewSystem(); 214 : 215 : /** 216 : * Sets the name of the field variable to actually use to flag elements for refinement / 217 : * coarsening. 218 : * This must be a CONSTANT, MONOMIAL Auxiliary Variable Name that contains values 219 : * corresponding to libMesh::Elem::RefinementState. 220 : * 221 : * @param marker_field The name of the field to use for refinement / coarsening. 222 : */ 223 : void setMarkerVariableName(std::string marker_field); 224 : 225 : /** 226 : * Sets the name of the field variable to actually use to flag elements for initial refinement / 227 : * coarsening. 228 : * This must be a CONSTANT, MONOMIAL Auxiliary Variable Name that contains values 229 : * corresponding to libMesh::Elem::RefinementState. 230 : * 231 : * @param marker_field The name of the field to use for refinement / coarsening. 232 : */ 233 : void setInitialMarkerVariableName(std::string marker_field); 234 : 235 : /** 236 : * Set the maximum refinement level (for the new Adaptivity system). 237 : */ 238 1825 : void setMaxHLevel(unsigned int level) { _max_h_level = level; } 239 : 240 : /** 241 : * Return the maximum h-level. 242 : */ 243 : unsigned int getMaxHLevel() { return _max_h_level; } 244 : 245 : /** 246 : * Set the interval (number of timesteps) between refinement steps. 247 : */ 248 2247 : void setInterval(unsigned int interval) { _interval = interval; } 249 : 250 : /** 251 : * Get an ErrorVector that will be filled up with values corresponding to the 252 : * indicator field name passed in. 253 : * 254 : * Note that this returns a reference... and the return value should be stored as a reference! 255 : * 256 : * @param indicator_field The name of the field to get an ErrorVector for. 257 : */ 258 : libMesh::ErrorVector & getErrorVector(const std::string & indicator_field); 259 : 260 : /** 261 : * Update the ErrorVectors that have been requested through calls to getErrorVector(). 262 : */ 263 : void updateErrorVectors(); 264 : 265 : /** 266 : * Query if an adaptivity step should be performed at the current time / time step 267 : */ 268 : bool isAdaptivityDue(); 269 : 270 : protected: 271 : FEProblemBase & _fe_problem; 272 : MooseMesh & _mesh; 273 : 274 : /// on/off flag reporting if the adaptivity is being used 275 : bool _mesh_refinement_on; 276 : /// on/off flag reporting if the adaptivity system has been initialized 277 : bool _initialized; 278 : /// A mesh refinement object to be used either with initial refinement or with Adaptivity. 279 : std::unique_ptr<libMesh::MeshRefinement> _mesh_refinement; 280 : /// Error estimator to be used by the apps. 281 : std::unique_ptr<libMesh::ErrorEstimator> _error_estimator; 282 : /// Error vector for use with the error estimator. 283 : std::unique_ptr<libMesh::ErrorVector> _error; 284 : 285 : std::shared_ptr<DisplacedProblem> _displaced_problem; 286 : 287 : /// A mesh refinement object for displaced mesh 288 : std::unique_ptr<libMesh::MeshRefinement> _displaced_mesh_refinement; 289 : 290 : /// the number of adaptivity steps to do at the beginning of simulation 291 : unsigned int _initial_steps; 292 : /// steps of adaptivity to perform 293 : unsigned int _steps; 294 : 295 : /// True if we want to print out info when mesh has changed 296 : bool _print_mesh_changed; 297 : 298 : /// Time 299 : Real & _t; 300 : /// Time Step 301 : int & _step; 302 : /// intreval between adaptivity runs 303 : unsigned int _interval; 304 : /// When adaptivity start 305 : Real _start_time; 306 : /// When adaptivity stops 307 : Real _stop_time; 308 : /// The number of adaptivity cycles per step 309 : unsigned int _cycles_per_step; 310 : 311 : /// Whether or not to use the "new" adaptivity system 312 : bool _use_new_system; 313 : 314 : /// Name of the marker variable if using the new adaptivity system 315 : std::string _marker_variable_name; 316 : 317 : /// Type of mesh adaptivity 318 : AdaptivityType _adaptivity_type; 319 : 320 : /// Name of the initial marker variable if using the new adaptivity system 321 : std::string _initial_marker_variable_name; 322 : 323 : /// The maximum number of refinement levels 324 : unsigned int _max_h_level; 325 : 326 : /// Whether or not to recompute markers during adaptivity cycles 327 : bool _recompute_markers_during_cycles; 328 : 329 : /// Sibling coupling object for HP adaptivity for evaluating data on elements' siblings in 330 : /// HPCoarsenTest 331 : std::unique_ptr<libMesh::SiblingCoupling> _sibling_coupling; 332 : 333 : /// Object for HP adaptivity 334 : std::unique_ptr<libMesh::HPCoarsenTest> _hp_coarsen_test; 335 : 336 : /// Stores pointers to ErrorVectors associated with indicator field names 337 : std::map<std::string, std::unique_ptr<libMesh::ErrorVector>> _indicator_field_to_error_vector; 338 : }; 339 : 340 : template <typename T> 341 : void 342 2110 : Adaptivity::setParam(const std::string & param_name, const T & param_value) 343 : { 344 2110 : if (param_name == "refine fraction") 345 : { 346 422 : _mesh_refinement->refine_fraction() = param_value; 347 422 : if (_displaced_mesh_refinement) 348 30 : _displaced_mesh_refinement->refine_fraction() = param_value; 349 : } 350 1688 : else if (param_name == "coarsen fraction") 351 : { 352 422 : _mesh_refinement->coarsen_fraction() = param_value; 353 422 : if (_displaced_mesh_refinement) 354 30 : _displaced_mesh_refinement->coarsen_fraction() = param_value; 355 : } 356 1266 : else if (param_name == "max h-level") 357 : { 358 422 : _mesh_refinement->max_h_level() = param_value; 359 422 : if (_displaced_mesh_refinement) 360 30 : _displaced_mesh_refinement->max_h_level() = param_value; 361 : } 362 844 : else if (param_name == "cycles_per_step") 363 422 : _cycles_per_step = param_value; 364 422 : else if (param_name == "recompute_markers_during_cycles") 365 422 : _recompute_markers_during_cycles = param_value; 366 : else 367 0 : mooseError("Invalid Param in adaptivity object"); 368 2110 : } 369 : #endif // LIBMESH_ENABLE_AMR