https://mooseframework.inl.gov
ElementSubdomainModifierBase.C
Go to the documentation of this file.
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 
11 #include "DisplacedProblem.h"
12 #include "MaterialWarehouse.h"
13 
14 #include "libmesh/dof_map.h"
15 #include "libmesh/remote_elem.h"
16 #include "libmesh/parallel_ghost_sync.h"
17 #include "libmesh/numeric_vector.h"
18 #include "libmesh/parameters.h"
19 #include <iterator>
20 #include <unordered_set>
21 
24 {
26 
27  // ESMs only operate on the undisplaced mesh to gather subdomain and sideset information. This
28  // information is used to modify both the undisplaced and the displaced meshes. It is the
29  // developer's responsibility to make sure the element IDs and sideset info are consistent across
30  // both meshes.
31  params.set<bool>("use_displaced_mesh") = false;
32  params.suppressParameter<bool>("use_displaced_mesh");
33 
34  params.addParam<std::vector<BoundaryName>>(
35  "moving_boundaries",
36  {},
37  "Moving boundaries between subdomains. These boundaries (both sidesets and nodesets) will be "
38  "updated for elements that change subdomain. The subdomains that each moving "
39  "boundary lies between shall be specified using the parameter "
40  "'moving_boundary_subdomain_pairs'. If one boundary and multiple subdomain pairs are "
41  "specified, then it is assumed that the pairs all apply to the boundary. A boundary will be "
42  "created on the mesh if it does not already exist.");
43  params.addParam<std::vector<std::vector<SubdomainName>>>(
44  "moving_boundary_subdomain_pairs",
45  {},
46  "The subdomain pairs associated with each moving boundary. For each pair of subdomains, only "
47  "the element side from the first subdomain will be added to the moving boundary, i.e., the "
48  "side normal is pointing from the first subdomain to the second subdomain. The pairs shall "
49  "be delimited by ';'. If a pair only has one subdomain, the moving boundary is associated "
50  "with the subdomain's external boundary, i.e., when the elements have no neighboring "
51  "elements.");
52 
53  params.addParam<std::vector<SubdomainName>>(
54  "reinitialize_subdomains",
55  {"ANY_BLOCK_ID"},
56  "By default, any element which changes subdomain is reinitialized. If a list of subdomains "
57  "(IDs or names) is provided, then only elements whose new subdomain is in the list will be "
58  "reinitialized. If an empty list is set, then no elements will be reinitialized.");
59 
60  params.addParam<bool>(
61  "old_subdomain_reinitialized",
62  true,
63  "This parameter must be set with a non-empty list in 'reinitialize_subdomains'. When set to "
64  "the default true, the element's old subdomain is not considered when determining if an "
65  "element should be reinitialized. If set to false, only elements whose old subdomain was not "
66  "in 'reinitialize_subdomains' are reinitialized. ");
67 
68  params.addParam<std::vector<VariableName>>(
69  "reinitialize_variables", {}, "Which variables to reinitialize when subdomain changes.");
70  MooseEnum reinit_strategy("IC POLYNOMIAL_NEIGHBOR POLYNOMIAL_WHOLE POLYNOMIAL_NEARBY", "IC");
71  params.addParam<std::vector<MooseEnum>>(
72  "reinitialization_strategy",
73  {reinit_strategy},
74  "The strategy used to reinitialize the solution when elements change subdomain. If multiple "
75  "strategies are provided, each strategy will be applied to the corresponding variable. If "
76  "only one strategy is provided, it will be applied to all variables.");
77  params.addParam<std::vector<UserObjectName>>(
78  "polynomial_fitters",
79  {},
80  "List of NodalPatchRecovery UserObjects used for polynomial fitting during variable "
81  "reinitialization. Required only if 'reinitialization_strategy' includes "
82  "POLYNOMIAL_NEIGHBOR, POLYNOMIAL_WHOLE, or POLYNOMIAL_NEARBY.");
83  params.addParam<int>(
84  "nearby_kd_tree_leaf_max_size",
85  10,
86  "Maximum number of elements in a leaf node of the K-D tree used to search for nearby "
87  "elements. Only needed if 'reinitialization_strategy' is set to POLYNOMIAL_NEARBY.");
88  params.addParam<double>(
89  "nearby_distance_threshold",
90  -1.0,
91  "Threshold for considering elements as 'nearby' in the K-D tree search. Only elements within "
92  "this distance will be considered for polynomial fitting.");
93  params.addParam<std::vector<bool>>(
94  "restore_overridden_dofs",
95  {},
96  "A list of boolean flags, one for each variable in 'reinitialize_variables', specifying "
97  "whether overridden DOF values should be restored after reinitialization for each variable. "
98  "This is useful when the solved values on these DOFs should be preserved. If the list is "
99  "empty, overridden DOF values will NOT be restored for any variable by default.");
100 
101  params.addParam<bool>("skip_restore_subdomain_changes",
102  false,
103  "Skip restoring the subdomain changes if the timestep is not advanced.");
104 
105  params.registerBase("MeshModifier");
106 
107  return params;
108 }
109 
111  : ElementUserObject(parameters),
112  _displaced_problem(_fe_problem.getDisplacedProblem().get()),
113  _displaced_mesh(_displaced_problem ? &_displaced_problem->mesh() : nullptr),
114  _nl_sys(_fe_problem.getNonlinearSystemBase(systemNumber())),
115  _aux_sys(_fe_problem.getAuxiliarySystem()),
116  _t_step_old(declareRestartableData<int>("t_step_old", 0)),
117  _restep(false),
118  _skip_restore_subdomain_changes(getParam<bool>("skip_restore_subdomain_changes")),
119  _old_subdomain_reinitialized(getParam<bool>("old_subdomain_reinitialized")),
120  _pr_names(getParam<std::vector<UserObjectName>>("polynomial_fitters")),
121  _reinit_vars(getParam<std::vector<VariableName>>("reinitialize_variables")),
122  _leaf_max_size(getParam<int>("nearby_kd_tree_leaf_max_size")),
123  _nearby_distance_threshold(getParam<double>("nearby_distance_threshold"))
124 {
125 }
126 
127 void
129 {
130  const std::vector<SubdomainName> subdomain_names_to_reinitialize =
131  getParam<std::vector<SubdomainName>>("reinitialize_subdomains");
132 
133  if (std::find(subdomain_names_to_reinitialize.begin(),
134  subdomain_names_to_reinitialize.end(),
135  "ANY_BLOCK_ID") != subdomain_names_to_reinitialize.end())
137  else
138  _subdomain_ids_to_reinitialize = _mesh.getSubdomainIDs(subdomain_names_to_reinitialize);
139 
140  std::set<SubdomainID> set_subdomain_ids_to_reinitialize(_subdomain_ids_to_reinitialize.begin(),
142 
143  if (_old_subdomain_reinitialized == false &&
147  set_subdomain_ids_to_reinitialize == _mesh.meshSubdomains()))
148  paramError("old_subdomain_reinitialized",
149  "'old_subdomain_reinitialized' can only be set to false if "
150  "reinitialize_subdomains does "
151  "not cover the whole model, otherwise no elements will be reinitialized as it is "
152  "impossible for an element's old subdomain to not be in the list.");
154  paramError("old_subdomain_reinitialized",
155  "'old_subdomain_reinitialized' can only be set to false if "
156  "reinitialize_subdomains is set to a non-empty list of subdomains, otherwise no "
157  "elements will be reinitialized, as it is impossible for an element's new subdomain "
158  "to be in the list.");
159 
160  auto bnd_names = getParam<std::vector<BoundaryName>>("moving_boundaries");
161  auto bnd_ids = _mesh.getBoundaryIDs(bnd_names, true);
162  const auto bnd_subdomains =
163  getParam<std::vector<std::vector<SubdomainName>>>("moving_boundary_subdomain_pairs");
164 
165  if (bnd_names.size() == 1 && bnd_subdomains.size() > 1)
166  {
167  bnd_names.insert(bnd_names.end(), bnd_subdomains.size() - 1, bnd_names[0]);
168  bnd_ids.insert(bnd_ids.end(), bnd_subdomains.size() - 1, bnd_ids[0]);
169  }
170  else if (bnd_names.size() != bnd_subdomains.size())
171  paramError("moving_boundary_subdomain_pairs",
172  "Each moving boundary must correspond to a pair of subdomains. ",
173  bnd_names.size(),
174  " boundaries are specified by the parameter 'moving_boundaries', while ",
175  bnd_subdomains.size(),
176  " subdomain pairs are provided. Alternatively, if one boundary and multiple "
177  "subdomain pairs are provided, then the subdomain pairs all apply to one boundary.");
178 
179  for (auto i : index_range(bnd_names))
180  {
181  _moving_boundary_names[bnd_ids[i]] = bnd_names[i];
182 
183  if (bnd_subdomains[i].size() == 2)
184  _moving_boundaries[{_mesh.getSubdomainID(bnd_subdomains[i][0]),
185  _mesh.getSubdomainID(bnd_subdomains[i][1])}] = bnd_ids[i];
186  else if (bnd_subdomains[i].size() == 1)
188  bnd_ids[i];
189  else
190  paramError("moving_boundary_subdomain_pairs",
191  "Each subdomain pair must contain 1 or 2 subdomain names, but ",
192  bnd_subdomains[i].size(),
193  " are given.");
194  }
195 
196  // Check if the variables to reinitialize exist in the system
197  for (const auto & var_name : _reinit_vars)
198  if (!_nl_sys.hasVariable(var_name) && !_aux_sys.hasVariable(var_name))
199  paramError("reinitialize_variables", "Variable ", var_name, " does not exist.");
200 
201  // Determine the reinitialization strategy for each variable.
202  // (1) If they are of the same size, we perform a 1-to-1 mapping.
203  // (2) If only one strategy or restore flag is provided, it applies to all variables.
204  const auto reinit_strategy_in = getParam<std::vector<MooseEnum>>("reinitialization_strategy");
205  const auto restore_overridden_dofs_in = getParam<std::vector<bool>>("restore_overridden_dofs");
206 
207  if (std::any_of(reinit_strategy_in.begin(),
208  reinit_strategy_in.end(),
209  [](const MooseEnum & val) { return val == "POLYNOMIAL_NEARBY"; }) &&
210  !isParamSetByUser("nearby_distance_threshold"))
211  mooseError("The 'nearby_distance_threshold' parameter must be set when using the "
212  "POLYNOMIAL_NEARBY reinitialization strategy.");
213 
214  if (reinit_strategy_in.size() == 1)
215  _reinit_strategy.resize(_reinit_vars.size(), reinit_strategy_in[0].getEnum<ReinitStrategy>());
216  else if (reinit_strategy_in.size() == _reinit_vars.size())
217  for (const auto & e : reinit_strategy_in)
218  _reinit_strategy.push_back(e.getEnum<ReinitStrategy>());
219  else
220  paramError(
221  "reinitialization_strategy",
222  "The 'reinitialization_strategy' parameter must have either a single value or a number "
223  "of values equal to the number of 'reinitialize_variables'. "
224  "Got ",
225  reinit_strategy_in.size(),
226  " strategies for ",
227  _reinit_vars.size(),
228  " variables.");
229 
230  if (restore_overridden_dofs_in.size() == 1)
231  {
232  if (restore_overridden_dofs_in[0])
234  _reinit_vars; // Restore overridden DOFs for all reinitialized variables
235  }
236  else if (restore_overridden_dofs_in.size() == _reinit_vars.size())
237  {
238  for (auto i : index_range(_reinit_vars))
239  if (restore_overridden_dofs_in[i])
241  }
242  else
243  {
244  if (!restore_overridden_dofs_in.empty())
245  paramError(
246  "restore_overridden_dofs",
247  "The 'restore_overridden_dofs' parameter must have either a single value or a number "
248  "of values equal to the number of 'reinitialize_variables'. "
249  "Got ",
250  restore_overridden_dofs_in.size(),
251  " restore_overridden_dofs for ",
252  _reinit_vars.size(),
253  " variables.");
254  }
255 
256  // For all the other variables, we set the reinitialization strategy to IC
257  for (const auto & var_name : _fe_problem.getVariableNames())
258  if (std::find(_reinit_vars.begin(), _reinit_vars.end(), var_name) == _reinit_vars.end())
259  {
260  _reinit_vars.push_back(var_name);
262  }
263 
264  // Map variable names to the index of the nodal patch recovery user object
265  _pr.resize(_pr_names.size());
266  std::size_t pr_count = 0;
267  for (auto i : index_range(_reinit_vars))
271  {
272  _var_name_to_pr_idx[_reinit_vars[i]] = pr_count;
273  if (pr_count >= _pr_names.size())
274  paramError("polynomial_fitters",
275  "The number of polynomial fitters (",
276  _pr_names.size(),
277  ") is less than the number of variables to reinitialize with polynomial "
278  "extrapolation.");
279  _pr[pr_count] =
281  mooseAssert(_pr[pr_count]->variableName() == _reinit_vars[i],
282  "The patch recovery UserObject's variable name must match the variable being "
283  "reinitialized in ElementSubdomainModifierBase.");
284  _depend_uo.insert(_pr_names[pr_count]);
285  pr_count++;
286  }
287  if (_pr_names.size() != pr_count)
288  paramError("polynomial_fitters",
289  "Mismatch between number of reinitialization strategies using polynomial "
290  "extrapolation and polynomial fitters (expected: ",
291  pr_count,
292  ", given: ",
293  _pr_names.size(),
294  ").");
295 }
296 
297 void
299 {
301  {
302  mooseInfoRepeated(name(), ": Restoring element subdomain changes.");
303 
304  // Reverse the subdomain changes
305  auto moved_elem_reversed = _moved_elems;
306  for (auto & [elem_id, subdomain] : moved_elem_reversed)
307  std::swap(subdomain.first, subdomain.second);
308 
309  _restep = true;
310  modify(moved_elem_reversed);
311  _restep = false;
312  }
313 
315 }
316 
317 void
319  const std::unordered_map<dof_id_type, std::pair<SubdomainID, SubdomainID>> & moved_elems)
320 {
321  if (!_restep)
322  // Cache the moved elements for potential restore
323  _moved_elems = moved_elems;
324 
325  // If nothing need to change, just return.
326  // This will skip all mesh changes, and so no adaptivity mesh files will be written.
327  auto n_moved_elem = moved_elems.size();
328  gatherSum(n_moved_elem);
329  if (n_moved_elem == 0)
330  return;
331 
332  // Create moving boundaries on the undisplaced and displaced meshes
333  //
334  // Note: We do this _everytime_ because previous execution might have removed the sidesets and
335  // nodesets. Most of the moving boundary algorithms below assume that the moving sidesets and
336  // nodesets already exist on the mesh.
338  if (_displaced_mesh)
340 
341  // This has to be done *before* subdomain changes are applied
342  findReinitializedElemsAndNodes(moved_elems);
343 
344  // Apply cached subdomain changes
345  applySubdomainChanges(moved_elems, _mesh);
346  if (_displaced_mesh)
347  applySubdomainChanges(moved_elems, *_displaced_mesh);
348 
349  // Update moving boundaries
350  gatherMovingBoundaryChanges(moved_elems);
352  if (_displaced_mesh)
354 
355  // Some variable reinitialization strategies require patch elements to be gathered
356  // This has to be done *before* reinitializing the equation systems because we need to find
357  // currently evaluable elements
358  if (!_restep)
359  {
360  _evaluable_elems.clear();
361  _patch_elem_ids.clear();
362  for (auto i : index_range(_reinit_vars))
364  }
365 
366  // Reinit equation systems
368  /*intermediate_change=*/false, /*contract_mesh=*/false, /*clean_refinement_flags=*/false);
369 
370  // Initialize solution and stateful material properties
371  if (!_restep)
372  {
373  applyIC();
376  }
377 }
378 
379 void
381 {
382  auto & bnd_info = mesh.getMesh().get_boundary_info();
383  for (const auto & [bnd_id, bnd_name] : _moving_boundary_names)
384  {
385  bnd_info.sideset_name(bnd_id) = bnd_name;
386  bnd_info.nodeset_name(bnd_id) = bnd_name;
387  }
388 }
389 
390 void
392  const std::unordered_map<dof_id_type, std::pair<SubdomainID, SubdomainID>> & moved_elems,
393  MooseMesh & mesh)
394 {
395  for (const auto & [elem_id, subdomain] : moved_elems)
396  {
397  // Change the element's subdomain ID
398  auto elem = mesh.elemPtr(elem_id);
399  const auto & [from, to] = subdomain;
400  mooseAssert(elem->subdomain_id() == from, "Inconsistent element subdomain ID.");
401  elem->subdomain_id() = to;
402 
403  // Change the ancestors' (if any) subdomain ID
404  setAncestorsSubdomainIDs(elem, to);
405  }
406 
407  // Synchronize ghost element subdomain changes
408  libMesh::SyncSubdomainIds sync(mesh.getMesh());
409  Parallel::sync_dofobject_data_by_id(
410  mesh.getMesh().comm(), mesh.getMesh().elements_begin(), mesh.getMesh().elements_end(), sync);
411 }
412 
413 void
415 {
416  auto curr_elem = elem;
417 
418  for (unsigned int i = curr_elem->level(); i > 0; --i)
419  {
420  // Change the parent's subdomain
421  curr_elem = curr_elem->parent();
422  curr_elem->subdomain_id() = subdomain_id;
423  }
424 }
425 
426 void
428  const std::unordered_map<dof_id_type, std::pair<SubdomainID, SubdomainID>> & moved_elems)
429 {
430  // Clear moving boundary changes from last execution
431  _add_element_sides.clear();
432  _add_neighbor_sides.clear();
433  _remove_element_sides.clear();
434  _remove_neighbor_sides.clear();
435 
436  const auto & sidesets = _mesh.getMesh().get_boundary_info().get_sideset_map();
437 
438  for (const auto & [elem_id, subdomain_assignment] : moved_elems)
439  {
440  auto elem = _mesh.elemPtr(elem_id);
441 
442  // The existing moving boundaries on the element side should be removed
443  for (auto itr = sidesets.lower_bound(elem); itr != sidesets.upper_bound(elem); itr++)
444  if (_moving_boundary_names.count(itr->second.second))
445  _remove_element_sides[elem->id()].emplace(itr->second.first, itr->second.second);
446 
447  for (auto side : elem->side_index_range())
448  {
449  auto neigh = elem->neighbor_ptr(side);
450 
451  // Don't mess with remote element neighbor
452  if (neigh && neigh == libMesh::remote_elem)
453  continue;
454  // If neighbor doesn't exist
455  else if (!neigh)
456  gatherMovingBoundaryChangesHelper(elem, side, nullptr, 0);
457  // If neighbor exists
458  else
459  {
460  auto neigh_side = neigh->which_neighbor_am_i(elem);
461 
462  if (neigh->active())
463  gatherMovingBoundaryChangesHelper(elem, side, neigh, neigh_side);
464  else
465  {
466  // Find the active neighbors of the element
467  std::vector<const Elem *> active_neighs;
468  // Neighbor has active children, they are neighbors of the element along that side
469  mooseAssert(!neigh->subactive(),
470  "The case where the active neighbor is an ancestor of this neighbor is not "
471  "handled at this time.");
472  neigh->active_family_tree_by_neighbor(active_neighs, elem);
473 
474  for (auto active_neigh : active_neighs)
475  gatherMovingBoundaryChangesHelper(elem, side, active_neigh, neigh_side);
476  }
477  }
478  }
479  }
480 }
481 
482 void
484  unsigned short side,
485  const Elem * neigh,
486  unsigned short neigh_side)
487 {
488  const auto & sidesets = _mesh.getMesh().get_boundary_info().get_sideset_map();
489 
490  // Detect element side change
491  SubdomainPair subdomain_pair = {elem->subdomain_id(),
492  neigh ? neigh->subdomain_id() : Moose::INVALID_BLOCK_ID};
493  if (_moving_boundaries.count(subdomain_pair))
494  _add_element_sides[elem->id()].emplace(side, _moving_boundaries.at(subdomain_pair));
495 
496  if (neigh)
497  {
498  // The existing moving boundaries on the neighbor side should be removed
499  for (auto itr = sidesets.lower_bound(neigh); itr != sidesets.upper_bound(neigh); itr++)
500  if (itr->second.first == neigh_side && _moving_boundary_names.count(itr->second.second))
501  _remove_neighbor_sides[neigh->id()].emplace(itr->second.first, itr->second.second);
502 
503  // Detect neighbor side change (by reversing the subdomain pair)
504  subdomain_pair = {subdomain_pair.second, subdomain_pair.first};
505  if (_moving_boundaries.count(subdomain_pair))
506  _add_neighbor_sides[neigh->id()].emplace(neigh_side, _moving_boundaries.at(subdomain_pair));
507  }
508 }
509 
510 void
512 {
513  auto & bnd_info = mesh.getMesh().get_boundary_info();
514 
515  // Remove all boundary nodes from the previous moving boundaries
516  auto nodesets = bnd_info.get_nodeset_map();
517  for (const auto & [node_id, bnd] : nodesets)
518  if (_moving_boundary_names.count(bnd))
519  bnd_info.remove_node(node_id, bnd);
520 
521  // Keep track of ghost element changes
522  std::unordered_map<processor_id_type,
523  std::vector<std::tuple<dof_id_type, unsigned short, BoundaryID>>>
524  add_ghost_sides, remove_ghost_sides;
525 
526  // Remove element sides from moving boundaries
527  for (const auto & [elem_id, sides] : _remove_element_sides)
528  for (const auto & [side, bnd] : sides)
529  bnd_info.remove_side(mesh.elemPtr(elem_id), side, bnd);
530 
531  // Remove neighbor sides from moving boundaries
532  for (const auto & [elem_id, sides] : _remove_neighbor_sides)
533  {
534  auto elem = mesh.elemPtr(elem_id);
535  for (const auto & [side, bnd] : sides)
536  {
537  bnd_info.remove_side(elem, side, bnd);
538  // Keep track of changes to ghosted elements
539  if (elem->processor_id() != processor_id())
540  remove_ghost_sides[elem->processor_id()].push_back({elem_id, side, bnd});
541  }
542  }
543 
544  Parallel::push_parallel_vector_data(
545  bnd_info.comm(),
546  remove_ghost_sides,
547  [&mesh,
548  &bnd_info](processor_id_type,
549  const std::vector<std::tuple<dof_id_type, unsigned short, BoundaryID>> & received)
550  {
551  for (const auto & [elem_id, side, bnd] : received)
552  bnd_info.remove_side(mesh.elemPtr(elem_id), side, bnd);
553  });
554 
555  // Add element sides to moving boundaries
556  for (const auto & [elem_id, sides] : _add_element_sides)
557  for (const auto & [side, bnd] : sides)
558  bnd_info.add_side(mesh.elemPtr(elem_id), side, bnd);
559 
560  // Add neighbor sides to moving boundaries
561  for (const auto & [elem_id, sides] : _add_neighbor_sides)
562  {
563  auto elem = mesh.elemPtr(elem_id);
564  for (const auto & [side, bnd] : sides)
565  {
566  bnd_info.add_side(elem, side, bnd);
567  // Keep track of changes to ghosted elements
568  if (elem->processor_id() != processor_id())
569  add_ghost_sides[elem->processor_id()].push_back({elem_id, side, bnd});
570  }
571  }
572 
573  Parallel::push_parallel_vector_data(
574  bnd_info.comm(),
575  add_ghost_sides,
576  [&mesh,
577  &bnd_info](processor_id_type,
578  const std::vector<std::tuple<dof_id_type, unsigned short, BoundaryID>> & received)
579  {
580  for (const auto & [elem_id, side, bnd] : received)
581  bnd_info.add_side(mesh.elemPtr(elem_id), side, bnd);
582  });
583 
584  bnd_info.parallel_sync_side_ids();
585  bnd_info.parallel_sync_node_ids();
586  mesh.update();
587 }
588 
589 void
591  ReinitStrategy reinit_strategy)
592 {
593  switch (reinit_strategy)
594  {
595  case ReinitStrategy::IC:
596  // No additional preparation needed for IC
597  break;
601  {
602  if (_var_name_to_pr_idx.find(var_name) == _var_name_to_pr_idx.end())
603  return;
604  const int pr_idx = _var_name_to_pr_idx[var_name];
605  // The patch elements might be different for each variable
606  gatherPatchElements(var_name, reinit_strategy);
607 
608  // Notify the patch recovery user object about the patch elements
609  _pr[pr_idx]->sync(_patch_elem_ids[var_name]);
610 
611  break;
612  }
613  default:
614  mooseError("Unknown reinitialization strategy");
615  break;
616  }
617 }
618 
619 void
621 {
622  // Clear cached ranges
626 
628  if (_displaced_mesh)
630 }
631 
632 void
634 {
635  auto & bnd_info = mesh.getMesh().get_boundary_info();
636  auto sidesets = bnd_info.get_sideset_map();
637  for (const auto & i : sidesets)
638  {
639  auto elem = i.first;
640  auto side = i.second.first;
641  auto bnd = i.second.second;
642  if (_moving_boundary_names.count(bnd) && !elem->active())
643  {
644  bnd_info.remove_side(elem, side, bnd);
645 
646  std::vector<const Elem *> elem_family;
647  elem->active_family_tree_by_side(elem_family, side);
648  for (auto felem : elem_family)
649  bnd_info.add_side(felem, side, bnd);
650  }
651  }
652 
653  bnd_info.parallel_sync_side_ids();
654  bnd_info.parallel_sync_node_ids();
655 }
656 
657 void
659  const std::unordered_map<dof_id_type, std::pair<SubdomainID, SubdomainID>> & moved_elems)
660 {
661  // Clear cached element reinitialization data
662  _reinitialized_elems.clear();
663  _reinitialized_nodes.clear();
664 
665  // One more algorithm:
666  // (1) Loop over moved elements
667  // (2) If neighbor element processor ID is not the same as current processor ID (ghost element),
668  // push the moved element ID to the neighbor processor
669 
670  std::unordered_map<processor_id_type, std::unordered_set<dof_id_type>> push_data_set;
671  std::unordered_map<processor_id_type, std::vector<dof_id_type>> push_data;
672 
673  for (const auto & [elem_id, subdomain] : moved_elems)
674  {
675  mooseAssert(_mesh.elemPtr(elem_id)->active(), "Moved elements should be active");
676  // Default: any element that changes subdomain is reinitialized
680  _reinitialized_elems.insert(elem_id);
681  else // Reinitialize if new subdomain is in list of subdomains to be reinitialized
682  {
683  const auto & [from, to] = subdomain;
685  _reinitialized_elems.insert(elem_id);
686  // Only reinitialize if original subdomain is not in list of subdomains
689  _reinitialized_elems.insert(elem_id);
690  else // New subdomain is not in list of subdomains
691  continue;
692  }
693  const auto & elem = _mesh.elemPtr(elem_id);
694 
695  // (1) Loop over nodes of moved elements
696  // (2) node to element map is used to find neighbor elements
697  // (3) If neighbor element processor ID is not the same as current processor ID (means that the
698  // current element is ghosted element to the neighbor processor), push the moved element (or
699  // reinitialized, or newly-activated) ID to the neighbor processor
700  for (const auto & node : elem->node_ref_range())
701  for (const auto & neigh_id : _mesh.nodeToElemMap().at(node.id()))
702  if (neigh_id != elem_id) // Don't check the element itself
703  {
704  const auto neigh_elem = _mesh.elemPtr(neigh_id);
705  if (neigh_elem->processor_id() != processor_id())
706  push_data_set[neigh_elem->processor_id()].insert(elem_id);
707  }
708 
709  for (unsigned int i = 0; i < elem->n_nodes(); ++i)
710  if (nodeIsNewlyReinitialized(elem->node_id(i)))
711  _reinitialized_nodes.insert(elem->node_id(i));
712  }
713 
714  for (auto & [pid, s] : push_data_set)
715  push_data[pid] = {s.begin(), s.end()};
716 
718 
719  auto push_receiver =
720  [this](const processor_id_type, const std::vector<dof_id_type> & received_data)
721  {
722  for (const auto & id : received_data)
724  };
725 
726  Parallel::push_parallel_vector_data(_mesh.comm(), push_data, push_receiver);
727 }
728 
729 bool
731 {
732  // Default: any element that changes subdomain is reinitialized
736  return true;
737 
738  // Is subdomain in list of subdomains to be reinitialized
741  id) != _subdomain_ids_to_reinitialize.end();
742 }
743 
744 bool
746 {
747  // If any of the node neighbor elements has reinitialized, then the node is NOT newly
748  // reinitialized.
749  for (auto neighbor_elem_id : _mesh.nodeToElemMap().at(node_id))
750  if (subdomainIsReinitialized(_mesh.elemPtr(neighbor_elem_id)->subdomain_id()))
751  return false;
752  return true;
753 }
754 
755 void
757 {
758  // Before reinitializing variables, some DOFs may be overwritten.
759  // By default, these overwritten DOF values are NOT restored.
760  // If the user sets `restore_overridden_dofs` to true, we first save the current
761  // values of these DOFs, then restore them after reinitialization.
762  for (const auto & var_name : _vars_to_restore_overridden_dofs)
763  storeOverriddenDofValues(var_name);
764 
765  // ReinitStrategy::IC
766  std::set<VariableName> ic_vars;
767  for (auto i : index_range(_reinit_vars))
769  ic_vars.insert(_reinit_vars[i]);
770  if (!ic_vars.empty())
773 
774  // ReinitStrategy::POLYNOMIAL_NEIGHBOR, POLYNOMIAL_WHOLE, POLYNOMIAL_NEARBY
775  for (const auto & [var, patch] : _patch_elem_ids)
777 
778  // See the comment above, now we restore the values of the dofs that were overridden
779  for (const auto & var_name : _vars_to_restore_overridden_dofs)
780  restoreOverriddenDofValues(var_name);
781 
782  mooseAssert(_fe_problem.numSolverSystems() <= 1,
783  "This code was written for a single nonlinear system");
784  // Set old and older solutions on the reinitialized dofs to the reinitialized values
785  // note: from current -> old -> older
791 
792  // Note: Need method to handle solve failures at timesteps where subdomain changes. The old
793  // solutions are now set to the reinitialized values. Does this impact restoring solutions
794 }
795 
796 void
798 {
799  const auto & sys = _fe_problem.getSystem(var_name);
800  const auto & current_solution = *sys.current_local_solution;
801  const auto & dof_map = sys.get_dof_map();
802  const auto & var = _fe_problem.getStandardVariable(0, var_name);
803  const auto var_num = var.number();
804 
805  // Get the DOFs on the reinitialized elements
806  // Here we should loop over both ghosted and local reinitialized elements.
807  // The ghosted elements here can take care of DoFs that is belong to the reinitialized
808  // elements but are not on the current processor.
809  std::set<dof_id_type> reinitialized_dofs;
810  for (const auto & elem_id : _semi_local_reinitialized_elems)
811  {
812  const auto & elem = _mesh.elemPtr(elem_id);
813  std::vector<dof_id_type> elem_dofs;
814  dof_map.dof_indices(elem, elem_dofs, var_num);
815  reinitialized_dofs.insert(elem_dofs.begin(), elem_dofs.end());
816  }
817 
818  // Get existing DOFs on the active elements excluding reinitialized elements
819  std::set<dof_id_type> existing_dofs;
820  for (const auto * elem : *_mesh.getActiveLocalElementRange())
821  {
822  if (_reinitialized_elems.count(elem->id()))
823  continue; // Skip reinitialized elements
824  std::vector<dof_id_type> elem_dofs;
825  dof_map.dof_indices(elem, elem_dofs, var_num);
826  existing_dofs.insert(elem_dofs.begin(), elem_dofs.end());
827  }
828 
829  // Get the DOFs on the nodes that are overridden on reinitialized elements
830  std::vector<dof_id_type> overridden_dofs;
831  std::set_intersection(reinitialized_dofs.begin(),
832  reinitialized_dofs.end(),
833  existing_dofs.begin(),
834  existing_dofs.end(),
835  std::back_inserter(overridden_dofs));
836 
837  // Values before overriding (to be restored later)
838  std::vector<Number> values;
839  for (auto dof : overridden_dofs)
840  values.push_back(current_solution(dof));
841 
842  _overridden_values_on_reinit_elems[var_name] = {overridden_dofs, values};
843 }
844 
845 void
847 {
848  const auto sn = _fe_problem.systemNumForVariable(var_name);
849  auto & sys = _fe_problem.getSystemBase(sn);
850  auto & sol = sys.solution();
851  const auto & dof_map = sys.dofMap();
852  const auto & [dof_ids, values] = _overridden_values_on_reinit_elems[var_name];
853 
854  std::unordered_map<processor_id_type, std::vector<std::pair<dof_id_type, Number>>> push_data;
855 
856  for (const int i : index_range(dof_ids))
857  {
858  if (dof_map.dof_owner(dof_ids[i]) == processor_id())
859  sol.set(dof_ids[i], values[i]);
860  else
861  push_data[dof_map.dof_owner(dof_ids[i])].emplace_back(dof_ids[i], values[i]);
862  }
863 
864  auto push_receiver = [&](const processor_id_type,
865  const std::vector<std::pair<dof_id_type, Number>> & received_data)
866  {
867  for (const auto & [id, value] : received_data)
868  sol.set(id, value);
869  };
870 
871  Parallel::push_parallel_vector_data(_mesh.comm(), push_data, push_receiver);
872 
873  sol.close();
874  sol.localize(*sys.system().current_local_solution, sys.dofMap().get_send_list());
875 }
876 
877 void
879 {
881 }
882 
885 {
887  return *_reinitialized_elem_range.get();
888 
889  // Create a vector of the newly reinitialized elements
890  std::vector<Elem *> elems;
891  for (auto elem_id : _reinitialized_elems)
892  elems.push_back(_mesh.elemPtr(elem_id));
893 
894  // Make some fake element iterators defining this vector of elements
895  Elem * const * elem_itr_begin = const_cast<Elem * const *>(elems.data());
896  Elem * const * elem_itr_end = elem_itr_begin + elems.size();
897 
898  const auto elems_begin = MeshBase::const_element_iterator(
899  elem_itr_begin, elem_itr_end, Predicates::NotNull<Elem * const *>());
900  const auto elems_end = MeshBase::const_element_iterator(
901  elem_itr_end, elem_itr_end, Predicates::NotNull<Elem * const *>());
902 
903  _reinitialized_elem_range = std::make_unique<ConstElemRange>(elems_begin, elems_end);
904 
905  return reinitializedElemRange();
906 }
907 
910 {
912  return *_reinitialized_bnd_node_range.get();
913 
914  // Create a vector of the newly reinitialized boundary nodes
915  std::vector<const BndNode *> nodes;
916  auto bnd_nodes = _mesh.getBoundaryNodeRange();
917  for (auto bnd_node : *bnd_nodes)
918  if (bnd_node->_node)
919  if (_reinitialized_nodes.count(bnd_node->_node->id()))
920  nodes.push_back(bnd_node);
921 
922  BndNode * const * bnd_node_itr_begin = const_cast<BndNode * const *>(nodes.data());
923  BndNode * const * bnd_node_itr_end = bnd_node_itr_begin + nodes.size();
924 
925  const auto bnd_nodes_begin = MooseMesh::const_bnd_node_iterator(
926  bnd_node_itr_begin, bnd_node_itr_end, Predicates::NotNull<const BndNode * const *>());
927  const auto bnd_nodes_end = MooseMesh::const_bnd_node_iterator(
928  bnd_node_itr_end, bnd_node_itr_end, Predicates::NotNull<const BndNode * const *>());
929 
931  std::make_unique<ConstBndNodeRange>(bnd_nodes_begin, bnd_nodes_end);
932 
933  return reinitializedBndNodeRange();
934 }
935 
938 {
940  return *_reinitialized_node_range.get();
941 
942  // Create a vector of the newly reinitialized nodes
943  std::vector<const Node *> nodes;
944 
945  for (auto node_id : _reinitialized_nodes)
946  nodes.push_back(_mesh.nodePtr(node_id)); // displaced mesh shares the same node object
947 
948  Node * const * node_itr_begin = const_cast<Node * const *>(nodes.data());
949  Node * const * node_itr_end = node_itr_begin + nodes.size();
950 
951  const auto nodes_begin = MeshBase::const_node_iterator(
952  node_itr_begin, node_itr_end, Predicates::NotNull<const Node * const *>());
953  const auto nodes_end = MeshBase::const_node_iterator(
954  node_itr_end, node_itr_end, Predicates::NotNull<const Node * const *>());
955 
956  _reinitialized_node_range = std::make_unique<ConstNodeRange>(nodes_begin, nodes_end);
957 
958  return *_reinitialized_node_range.get();
959 }
960 
961 void
963  ConstElemRange & elem_range,
964  ConstBndNodeRange & bnd_node_range)
965 {
966  for (auto bnd : bnd_node_range)
967  {
968  const Node * bnode = bnd->_node;
969  if (!bnode)
970  continue;
971  }
972 
973  // Don't do anything if this is a steady simulation
974  if (!sys.hasSolutionState(1))
975  return;
976 
977  NumericVector<Number> & current_solution = *sys.system().current_local_solution;
978  NumericVector<Number> & old_solution = sys.solutionOld();
979  NumericVector<Number> * older_solution = sys.hasSolutionState(2) ? &sys.solutionOlder() : nullptr;
980 
981  // Get dofs for the reinitialized elements and nodes
982  DofMap & dof_map = sys.dofMap();
983  std::vector<dof_id_type> dofs;
984 
985  for (auto & elem : elem_range)
986  {
987  std::vector<dof_id_type> elem_dofs;
988  dof_map.dof_indices(elem, elem_dofs);
989  dofs.insert(dofs.end(), elem_dofs.begin(), elem_dofs.end());
990  }
991 
992  for (auto & bnd_node : bnd_node_range)
993  {
994  std::vector<dof_id_type> bnd_node_dofs;
995  dof_map.dof_indices(bnd_node->_node, bnd_node_dofs);
996  dofs.insert(dofs.end(), bnd_node_dofs.begin(), bnd_node_dofs.end());
997  }
998 
999  // Set the old and older solutions to match the reinitialization
1000  for (auto dof : dofs)
1001  {
1002  old_solution.set(dof, current_solution(dof));
1003  if (older_solution)
1004  older_solution->set(dof, current_solution(dof));
1005  }
1006 
1007  old_solution.close();
1008  if (older_solution)
1009  older_solution->close();
1010 }
1011 
1012 void
1014  ReinitStrategy reinit_strategy)
1015 {
1016  _patch_elem_ids[var_name].clear();
1017 
1018  // First collect all elements who own dofs in the current dofmap
1019  auto & sys = _fe_problem.getSystem(var_name);
1020 
1021  // Cache evaluable elements for the system if not already done
1022  if (!_evaluable_elems.count(sys.number()))
1023  {
1024  auto & [candidate_elems, candidate_elem_ids] = _evaluable_elems[sys.number()];
1025  const auto & dof_map = sys.get_dof_map();
1026  std::vector<dof_id_type> elem_dofs;
1027  auto vn = sys.variable_number(static_cast<std::string>(var_name));
1028  for (const auto elem : *_mesh.getActiveLocalElementRange())
1029  {
1030  if (std::find(_reinitialized_elems.begin(), _reinitialized_elems.end(), elem->id()) !=
1031  _reinitialized_elems.end())
1032  continue; // Skip elements that were reinitialized
1033 
1034  dof_map.dof_indices(elem, elem_dofs, vn);
1035  if (!elem_dofs.empty())
1036  {
1037  candidate_elems.insert(elem);
1038  candidate_elem_ids.push_back(elem->id());
1039  }
1040  }
1041  }
1042  auto & [candidate_elems, candidate_elem_ids] = _evaluable_elems[sys.number()];
1043 
1044  // Now we gather patch elements based on the reinit strategy
1045  auto & patch_elems = _patch_elem_ids[var_name];
1046 
1047  switch (reinit_strategy)
1048  {
1050  {
1051  auto has_neighbor_in_reinit_elems = [&](const Elem * elem) -> bool
1052  {
1053  for (const auto & node : elem->node_ref_range())
1054  for (const auto & neigh_id : _mesh.nodeToElemMap().at(node.id()))
1055  // here we need to use _global_reinitialized_elems gathering from all processors
1056  if (_semi_local_reinitialized_elems.count(neigh_id))
1057  return true;
1058  return false;
1059  };
1060  // Loop over all candidate elements, for each element, if any of its point neighbor belongs
1061  // to the reinitialized elements, we will include that element in the patch element set.
1062  for (const auto * elem : candidate_elems)
1063  if (has_neighbor_in_reinit_elems(elem))
1064  patch_elems.push_back(elem->id());
1065  break;
1066  }
1068  {
1069  // This is simple: all candidate elements are patch elements
1070  patch_elems = candidate_elem_ids;
1071  break;
1072  }
1074  {
1075  std::vector<Point> kd_points;
1076  std::vector<dof_id_type> global_candidate_elem_ids;
1077 
1078  if (_mesh.isDistributedMesh())
1079  {
1080  std::vector<std::pair<Point, dof_id_type>> pts_ids(candidate_elem_ids.size());
1081  for (std::size_t i = 0; i < candidate_elem_ids.size(); ++i)
1082  pts_ids[i] = {_mesh.elemPtr(candidate_elem_ids[i])->vertex_average(),
1083  candidate_elem_ids[i]};
1084  _mesh.comm().allgather(pts_ids);
1085  for (const auto & [pt, id] : pts_ids)
1086  {
1087  kd_points.push_back(pt);
1088  global_candidate_elem_ids.push_back(id);
1089  }
1090  }
1091  else
1092  {
1093  _mesh.comm().allgather(candidate_elem_ids);
1094  global_candidate_elem_ids = candidate_elem_ids;
1095  for (const auto & id : candidate_elem_ids)
1096  kd_points.push_back(_mesh.elemPtr(id)->vertex_average());
1097  }
1098 
1099  const auto kd_tree = std::make_unique<KDTree>(kd_points, _leaf_max_size);
1100 
1101  std::vector<nanoflann::ResultItem<std::size_t, Real>> query_result;
1102  for (const auto & elem_id : _reinitialized_elems)
1103  {
1104  const Point & centroid = _mesh.elemPtr(elem_id)->vertex_average();
1105  kd_tree->radiusSearch(centroid, _nearby_distance_threshold, query_result);
1106  for (const auto & [qid, dist] : query_result)
1107  patch_elems.push_back(global_candidate_elem_ids[qid]);
1108  }
1109  break;
1110  }
1111  default:
1112  mooseError("Unknown reinitialization strategy");
1113  break;
1114  }
1115 
1116  // every processor should have the same patch elements to do the polynomial extrapolation,
1117  // so we gather them across all processors
1118  _mesh.comm().allgather(patch_elems);
1119 
1120  // Remove duplicates from the patch elements (espcially important for POLYNOMIAL_NEARBY)
1121  std::sort(patch_elems.begin(), patch_elems.end());
1122  patch_elems.erase(std::unique(patch_elems.begin(), patch_elems.end()), patch_elems.end());
1123 }
1124 
1125 void
1127 {
1128  const auto & coef =
1129  _pr[_var_name_to_pr_idx[var_name]]->getCachedCoefficients(_patch_elem_ids[var_name]);
1130 
1131  const unsigned dim = _mesh.dimension();
1132 
1133  libMesh::Parameters function_parameters;
1134 
1135  const auto & multi_index = _pr[_var_name_to_pr_idx[var_name]]->multiIndex();
1136 
1137  function_parameters.set<std::vector<std::vector<unsigned int>>>("multi_index") = multi_index;
1138 
1139  std::vector<Real> coef_vec(coef.size());
1140  for (auto i = 0; i < coef.size(); ++i)
1141  coef_vec[i] = coef(i);
1142 
1143  function_parameters.set<std::vector<Real>>("multi_index_coefficients") = coef_vec;
1144  function_parameters.set<unsigned int>("dimension_for_projection") = dim;
1145 
1146  // Define projection function
1147  auto poly_func = [](const Point & p,
1149  const std::string &,
1150  const std::string &) -> libMesh::Number
1151  {
1152  const auto & multi_index =
1153  parameters.get<std::vector<std::vector<unsigned int>>>("multi_index");
1154  const auto & coeffs = parameters.get<std::vector<Real>>("multi_index_coefficients");
1155 
1156  Real val = 0.0;
1157 
1158  for (unsigned int r = 0; r < multi_index.size(); r++)
1159  {
1160  Real monomial = 1.0;
1161  for (unsigned int d = 0; d < multi_index[r].size(); d++)
1162  {
1163  const auto power = multi_index[r][d];
1164  if (power == 0)
1165  continue;
1166 
1167  monomial *= std::pow(p(d), power);
1168  }
1169  val += coeffs[r] * monomial;
1170  }
1171 
1172  return val;
1173  };
1174 
1175  // Define gradient
1176  auto poly_func_grad = [](const Point & p,
1178  const std::string &,
1179  const std::string &) -> libMesh::Gradient
1180  {
1181  const unsigned int dim = parameters.get<unsigned int>("dimension_for_projection");
1182 
1183  const auto & multi_index =
1184  parameters.get<std::vector<std::vector<unsigned int>>>("multi_index");
1185  const auto & coeffs = parameters.get<std::vector<Real>>("multi_index_coefficients");
1186 
1187  libMesh::Gradient grad; // Zero-initialized
1188 
1189  for (unsigned int r = 0; r < multi_index.size(); ++r)
1190  {
1191  const auto & powers = multi_index[r];
1192  const Real coef = coeffs[r];
1193 
1194  for (unsigned int d = 0; d < dim; ++d) // Loop over dimension
1195  {
1196  const auto power_d = powers[d];
1197  if (power_d == 0)
1198  continue;
1199 
1200  // Compute partial derivative in direction d
1201  Real partial = coef * power_d;
1202 
1203  for (unsigned int i = 0; i < powers.size(); ++i)
1204  {
1205  if (i == d)
1206  {
1207  if (powers[i] > 1)
1208  partial *= std::pow(p(i), powers[i] - 1); // reduce power by 1
1209  }
1210  else
1211  {
1212  if (powers[i] > 0)
1213  partial *= std::pow(p(i), powers[i]); // full power
1214  }
1215  }
1216 
1217  grad(d) += partial;
1218  }
1219  }
1220 
1221  return grad;
1222  };
1223 
1225  reinitializedElemRange(), poly_func, poly_func_grad, function_parameters, var_name);
1226 }
ConstNodeRange & reinitializedNodeRange()
Range of reinitialized nodes.
const std::vector< UserObjectName > _pr_names
Names of the NodalPatchRecoveryVariable user objects.
virtual void meshChanged(bool intermediate_change, bool contract_mesh, bool clean_refinement_flags)
Update data after a mesh change.
bool _restep
Whether this is a re-step.
int & _t_step_old
Previous time step number.
void allgather(const T &send_data, std::vector< T, A > &recv_data) const
libMesh::ConstElemRange * getActiveLocalElementRange()
Return pointers to range objects for various types of ranges (local nodes, boundary elems...
Definition: MooseMesh.C:1276
KOKKOS_INLINE_FUNCTION const T * find(const T &target, const T *const begin, const T *const end)
Find a value in an array.
Definition: KokkosUtils.h:30
virtual libMesh::System & getSystem(const std::string &var_name) override
Returns the equation system containing the variable provided.
void projectFunctionOnCustomRange(ConstElemRange &elem_range, Number(*func)(const Point &, const libMesh::Parameters &, const std::string &, const std::string &), Gradient(*func_grad)(const Point &, const libMesh::Parameters &, const std::string &, const std::string &), const libMesh::Parameters &params, const VariableName &target_var)
Project a function onto a range of elements for a given variable.
T & getUserObject(const std::string &name, unsigned int tid=0) const
Get the user object by its name.
std::map< VariableName, unsigned int > _var_name_to_pr_idx
map from variable name to the index of the nodal patch recovery user object in _pr ...
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:439
NonlinearSystemBase & _nl_sys
Nonlinear system.
virtual Elem * elemPtr(const dof_id_type i)
Definition: MooseMesh.C:3153
void prepareVariableForReinitialization(const VariableName &var_name, ReinitStrategy reinit_strategy)
NumericVector< Number > & solution()
Definition: SystemBase.h:196
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
unsigned int number() const
Get variable number coming from libMesh.
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
static InputParameters validParams()
std::unordered_set< dof_id_type > _reinitialized_elems
Reinitialized elements.
void swap(std::vector< T > &data, const std::size_t idx0, const std::size_t idx1, const libMesh::Parallel::Communicator &comm)
Swap function for serial or distributed vector of data.
Definition: Shuffle.h:494
std::pair< SubdomainID, SubdomainID > SubdomainPair
Moving boundaries associated with each subdomain pair.
std::vector< VariableName > _vars_to_restore_overridden_dofs
List of variable names for which overridden DOF values should be restored.
void initElementStatefulProps()
Reinitialize stateful material properties on range of elements and nodes to be reinitialized.
void gatherPatchElements(const VariableName &var_name, ReinitStrategy reinit_strategy)
Gather patch elements for reinitialized elements based on the reinitialization strategy.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
std::unordered_map< BoundaryID, BoundaryName > _moving_boundary_names
Boundary names associated with each moving boundary ID.
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
Definition: MooseUtils.h:1133
std::unordered_map< dof_id_type, std::pair< SubdomainID, SubdomainID > > _moved_elems
Cached moved elements for potential restore.
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
virtual libMesh::System & system()=0
Get the reference to the libMesh system.
MeshBase & mesh
void gatherMovingBoundaryChanges(const std::unordered_map< dof_id_type, std::pair< SubdomainID, SubdomainID >> &moved_elems)
void mooseInfoRepeated(Args &&... args)
Emit an informational message with the given stringified, concatenated args.
Definition: MooseError.h:398
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:159
const Parallel::Communicator & comm() const
The definition of the const_bnd_node_iterator struct.
Definition: MooseMesh.h:2054
std::map< unsigned int, std::pair< std::unordered_set< const Elem * >, std::vector< dof_id_type > > > _evaluable_elems
local evaluable elements before reinitializing the equation systems Key of the map is the system numb...
std::vector< VariableName > _reinit_vars
List of variable names to be initialized for IC.
NumericVector< Number > & solutionOlder()
Definition: SystemBase.h:198
const MaterialWarehouse & getMaterialWarehouse() const
void initElementStatefulProps(const libMesh::ConstElemRange &elem_range, const bool threaded)
Initialize stateful properties for elements in a specific elem_range This is needed when elements/bou...
Base class for a system (of equations)
Definition: SystemBase.h:84
void gatherMovingBoundaryChangesHelper(const Elem *elem, unsigned short side, const Elem *neigh, unsigned short neigh_side)
std::vector< SubdomainID > getSubdomainIDs(const std::vector< SubdomainName > &subdomain_names) const
Get the associated subdomainIDs for the subdomain names that are passed in.
Definition: MooseMesh.C:1775
void initialSetup() override
Gets called at the beginning of the simulation before this object is asked to do its job...
bool _skip_restore_subdomain_changes
Skipping restoring the subdomain changes if the timestep is not advanced.
StoredRange< MeshBase::const_element_iterator, const Elem *> ConstElemRange
void suppressParameter(const std::string &name)
This method suppresses an inherited parameter so that it isn&#39;t required or valid in the derived class...
void updateAMRMovingBoundary(MooseMesh &mesh)
Update boundaries for adaptive mesh from the parent to children elements.
const SubdomainID INVALID_BLOCK_ID
Definition: MooseTypes.C:20
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
virtual void modify(const std::unordered_map< dof_id_type, std::pair< SubdomainID, SubdomainID >> &moved_elems)
Modify the element subdomains.
uint8_t processor_id_type
virtual libMesh::DofMap & dofMap()
Gets writeable reference to the dof map.
Definition: SystemBase.C:1163
void createMovingBoundaries(MooseMesh &mesh)
Create moving boundaries.
virtual std::vector< VariableName > getVariableNames()
Returns a list of all the variables in the problem (both from the NL and Aux systems.
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
void storeOverriddenDofValues(const VariableName &var_name)
Store values from non-reinitialized DoFs on reinitialized elements Stores the value before re-initial...
int & _t_step
The number of the time step.
ReinitStrategy
Strategies for (re)initializing the solution:
void gatherSum(T &value)
Gather the parallel sum of the variable passed in.
Definition: UserObject.h:126
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3488
ConstElemRange & reinitializedElemRange()
Range of reinitialized elements.
std::unordered_set< dof_id_type > _semi_local_reinitialized_elems
Semi-local reinitialized elements: ghosted and local reinitialized elements.
char ** sides
virtual unsigned int dimension() const
Returns MeshBase::mesh_dimension(), (not MeshBase::spatial_dimension()!) of the underlying libMesh me...
Definition: MooseMesh.C:2968
void meshChanged() override
Called on this object when the mesh changes.
void setAncestorsSubdomainIDs(Elem *elem, const SubdomainID subdomain_id)
Change the subdomain ID of all ancestor elements.
void applySubdomainChanges(const std::unordered_map< dof_id_type, std::pair< SubdomainID, SubdomainID >> &moved_elems, MooseMesh &mesh)
std::unordered_set< dof_id_type > _reinitialized_nodes
Reinitialized nodes.
virtual const Node * nodePtr(const dof_id_type i) const
Definition: MooseMesh.C:863
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:92
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
bool subdomainIsReinitialized(SubdomainID id) const
Determine if a subdomain is to be reinitialized.
int _leaf_max_size
KD-tree related members.
std::vector< SubdomainID > _subdomain_ids_to_reinitialize
Reinitialize moved elements whose new subdomain is in this list.
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
SystemBase & _sys
Reference to the system object for this user object.
Definition: UserObject.h:215
virtual bool hasSolutionState(const unsigned int state, Moose::SolutionIterationType iteration_type=Moose::SolutionIterationType::Time) const
Whether or not the system has the solution state (0 = current, 1 = old, 2 = older, etc).
Definition: SystemBase.h:1090
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1157
std::map< VariableName, std::pair< std::vector< dof_id_type >, std::vector< Number > > > _overridden_values_on_reinit_elems
A map from variable name to a pair of: (1) a vector of DOF IDs associated with non-reinitialized node...
AuxiliarySystem & getAuxiliarySystem()
unsigned int systemNumForVariable(const VariableName &variable_name) const
virtual bool hasVariable(const std::string &var_name) const
Query a system for a variable.
Definition: SystemBase.C:851
std::vector< NodalPatchRecoveryVariable * > _pr
Apply initial conditions using polynomial extrapolation.
virtual void close()=0
AuxiliarySystem & _aux_sys
Auxiliary system.
virtual MooseVariable & getStandardVariable(const THREAD_ID tid, const std::string &var_name) override
Returns the variable reference for requested MooseVariable which may be in any system.
std::unordered_map< dof_id_type, std::unordered_map< unsigned short, BoundaryID > > _add_element_sides
Element sides to be added.
void restoreOverriddenDofValues(const VariableName &var_name)
Restore values to non-reinitialized DoFs on reinitialized elements.
std::set< std::string > _depend_uo
Depend UserObjects that to be used both for determining user object sorting and by AuxKernel for find...
Definition: UserObject.h:229
virtual const SystemBase & getSystemBase(const unsigned int sys_num) const
Get constant reference to a system in this problem.
const SubdomainID ANY_BLOCK_ID
Definition: MooseTypes.C:19
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unordered_map< dof_id_type, std::unordered_map< unsigned short, BoundaryID > > _remove_element_sides
Element sides to be removed.
bool hasActiveObjects(THREAD_ID tid=0) const
void applyMovingBoundaryChanges(MooseMesh &mesh)
FEProblemBase & _fe_problem
Reference to the FEProblemBase for this user object.
Definition: UserObject.h:211
bool nodeIsNewlyReinitialized(dof_id_type node_id) const
Determine if a node is newly reinitialized.
T & set(const std::string &)
std::unique_ptr< ConstBndNodeRange > _reinitialized_bnd_node_range
Range of reinitialized boundary nodes.
std::map< VariableName, std::vector< dof_id_type > > _patch_elem_ids
A map from variable names to their corresponding patch element IDs.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
MooseMesh * _displaced_mesh
Displaced mesh.
std::vector< ReinitStrategy > _reinit_strategy
The strategies used to apply IC on newly activated elements, for each variable.
void applyIC()
Reinitialize variables on range of elements and nodes to be reinitialized.
std::unique_ptr< NumericVector< Number > > current_local_solution
std::unique_ptr< ConstElemRange > _reinitialized_elem_range
Range of reinitialized elements.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual void set(const numeric_index_type i, const Number value)=0
void setOldAndOlderSolutions(SystemBase &sys, ConstElemRange &elem_range, ConstBndNodeRange &bnd_node_range)
Set old and older solutions to reinitialized elements and nodes.
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const
Returns a vector of boundary IDs for the requested element on the requested side. ...
Patch recovery from a coupled variable.
void findReinitializedElemsAndNodes(const std::unordered_map< dof_id_type, std::pair< SubdomainID, SubdomainID >> &moved_elems)
Real Number
NumericVector< Number > & solutionOld()
Definition: SystemBase.h:197
processor_id_type processor_id() const
virtual bool isDistributedMesh() const
Returns the final Mesh distribution type.
Definition: MooseMesh.h:1012
virtual std::size_t numSolverSystems() const override
void timestepSetup() override
Gets called at the beginning of the timestep before this object is asked to do its job...
std::unique_ptr< ConstNodeRange > _reinitialized_node_range
Range of reinitialized nodes.
bool isParamSetByUser(const std::string &name) const
Test if the supplied parameter is set by a user, as opposed to not set or set to default.
Definition: MooseBase.h:205
MooseUnits pow(const MooseUnits &, int)
Definition: Units.C:537
std::unordered_map< SubdomainPair, BoundaryID > _moving_boundaries
ConstBndNodeRange & reinitializedBndNodeRange()
Range of reinitialized boundary nodes.
double _nearby_distance_threshold
Radius threshold for the k-d tree neighbor search.
libMesh::StoredRange< MooseMesh::const_bnd_node_iterator, const BndNode * > * getBoundaryNodeRange()
Definition: MooseMesh.C:1327
void ErrorVector unsigned int
auto index_range(const T &sizable)
std::unordered_map< dof_id_type, std::unordered_map< unsigned short, BoundaryID > > _remove_neighbor_sides
Neighbor sides to be removed.
std::unordered_map< dof_id_type, std::unordered_map< unsigned short, BoundaryID > > _add_neighbor_sides
Neighbor sides to be added.
ElementSubdomainModifierBase(const InputParameters &parameters)
void projectInitialConditionOnCustomRange(libMesh::ConstElemRange &elem_range, ConstBndNodeRange &bnd_node_range, const std::optional< std::set< VariableName >> &target_vars=std::nullopt)
Project initial conditions for custom elem_range and bnd_node_range This is needed when elements/boun...
const bool _old_subdomain_reinitialized
Whether to reinitialize moved elements whose old subdomain was in _reinitialize_subdomains.
uint8_t dof_id_type
const std::map< dof_id_type, std::vector< dof_id_type > > & nodeToElemMap()
If not already created, creates a map from every node to all elements to which they are connected...
Definition: MooseMesh.C:1216
const std::set< SubdomainID > & meshSubdomains() const
Returns a read-only reference to the set of subdomains currently present in the Mesh.
Definition: MooseMesh.C:3211
SubdomainID getSubdomainID(const SubdomainName &subdomain_name) const
Get the associated subdomain ID for the subdomain name.
Definition: MooseMesh.C:1769
void extrapolatePolynomial(const VariableName &var_name)
Extrapolate polynomial for the given variable onto the reinitialized elements.
const RemoteElem * remote_elem