https://mooseframework.inl.gov
MooseMesh.h
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 
10 #pragma once
11 
12 #ifdef MOOSE_KOKKOS_ENABLED
13 #include "KokkosMesh.h"
14 #endif
15 
16 #include "MooseObject.h"
17 #include "BndNode.h"
18 #include "BndElement.h"
19 #include "Restartable.h"
20 #include "MooseEnum.h"
21 #include "PerfGraphInterface.h"
22 #include "MooseHashing.h"
23 #include "MooseApp.h"
24 #include "FaceInfo.h"
25 #include "ElemInfo.h"
26 
27 #include <memory> //std::unique_ptr
28 #include <unordered_map>
29 #include <unordered_set>
30 
31 // libMesh
32 #include "libmesh/elem_range.h"
33 #include "libmesh/mesh_base.h"
34 #include "libmesh/replicated_mesh.h"
35 #include "libmesh/distributed_mesh.h"
36 #include "libmesh/node_range.h"
37 #include "libmesh/nanoflann.hpp"
38 #include "libmesh/vector_value.h"
39 #include "libmesh/point.h"
40 #include "libmesh/partitioner.h"
41 
42 class Assembly;
44 class MooseVariableBase;
46 class MooseUnits;
47 
48 // libMesh forward declarations
49 namespace libMesh
50 {
51 class ExodusII_IO;
52 class QBase;
53 class PeriodicBoundaries;
54 class Partitioner;
55 class GhostingFunctor;
56 class BoundingBox;
57 }
58 // Useful typedefs
60 
61 // List of supported geometrical elements
62 const std::string LIST_GEOM_ELEM = "EDGE EDGE2 EDGE3 EDGE4 "
63  "QUAD QUAD4 QUAD8 QUAD9 "
64  "TRI TRI3 TRI6 TRI7 "
65  "HEX HEX8 HEX20 HEX27 "
66  "TET TET4 TET10 TET14 "
67  "PRISM PRISM6 PRISM15 PRISM18 "
68  "PYRAMID PYRAMID5 PYRAMID13 PYRAMID14";
69 
73 class QpMap
74 {
75 public:
76  QpMap() : _distance(std::numeric_limits<Real>::max()) {}
77 
79  unsigned int _from;
80 
82  unsigned int _to;
83 
86 };
87 
92 class MooseMesh : public MooseObject, public Restartable, public PerfGraphInterface
93 {
94 public:
99 
101  MooseMesh(const MooseMesh & other_mesh);
102  MooseMesh() = delete;
103  MooseMesh & operator=(const MooseMesh & other_mesh) = delete;
104 
105  virtual ~MooseMesh();
106 
107  // The type of libMesh::MeshBase that will be used
108  enum class ParallelType
109  {
110  DEFAULT,
111  REPLICATED,
113  };
114 
118  virtual MooseMesh & clone() const;
119 
125  virtual std::unique_ptr<MooseMesh> safeClone() const = 0;
126 
131 
137  std::unique_ptr<MeshBase> buildMeshBaseObject(unsigned int dim = libMesh::invalid_uint);
138 
146  template <typename T>
147  std::unique_ptr<T> buildTypedMesh(unsigned int dim = libMesh::invalid_uint);
148 
153  void setMeshBase(std::unique_ptr<MeshBase> mesh_base);
154 
156  static MooseEnum partitioning();
157 
159  static MooseEnum elemTypes();
160 
167  virtual void init();
168 
174  virtual void buildMesh() = 0;
175 
181  virtual unsigned int dimension() const;
182 
186  virtual unsigned int spatialDimension() const { return _mesh->spatial_dimension(); }
187 
193  virtual unsigned int effectiveSpatialDimension() const;
194 
198  unsigned int getBlocksMaxDimension(const std::vector<SubdomainName> & blocks) const;
199 
204  std::vector<BoundaryID> getBoundaryIDs(const Elem * const elem,
205  const unsigned short int side) const;
206 
213  const Elem * getLowerDElem(const Elem *, unsigned short int) const;
214 
218  unsigned int getHigherDSide(const Elem * elem) const;
219 
227  const std::set<BoundaryID> & getBoundaryIDs() const;
228 
235  void buildNodeList();
236  void buildBndElemList();
237 
242  const std::map<dof_id_type, std::vector<dof_id_type>> & nodeToElemMap();
243 
251  const std::map<dof_id_type, std::vector<dof_id_type>> & nodeToActiveSemilocalElemMap();
252 
257  struct bnd_node_iterator;
258  struct const_bnd_node_iterator;
259 
260  struct bnd_elem_iterator;
261  struct const_bnd_elem_iterator;
262 
266  virtual bnd_node_iterator bndNodesBegin();
267  virtual bnd_node_iterator bndNodesEnd();
268 
272  virtual bnd_elem_iterator bndElemsBegin();
273  virtual bnd_elem_iterator bndElemsEnd();
274 
279 
286  void buildSideList(std::vector<dof_id_type> & el,
287  std::vector<unsigned short int> & sl,
288  std::vector<boundary_id_type> & il);
292  std::vector<std::tuple<dof_id_type, unsigned short int, boundary_id_type>> buildSideList();
293 
298  std::vector<std::tuple<dof_id_type, unsigned short int, boundary_id_type>>
299  buildActiveSideList() const;
300 
304  unsigned int sideWithBoundaryID(const Elem * const elem, const BoundaryID boundary_id) const;
305 
309  MeshBase::node_iterator localNodesBegin();
310  MeshBase::node_iterator localNodesEnd();
311  MeshBase::const_node_iterator localNodesBegin() const;
312  MeshBase::const_node_iterator localNodesEnd() const;
313 
317  MeshBase::element_iterator activeLocalElementsBegin();
318  const MeshBase::element_iterator activeLocalElementsEnd();
319  MeshBase::const_element_iterator activeLocalElementsBegin() const;
320  const MeshBase::const_element_iterator activeLocalElementsEnd() const;
321 
325  virtual dof_id_type nNodes() const;
326  virtual dof_id_type nElem() const;
327 
328  virtual dof_id_type nLocalNodes() const { return _mesh->n_local_nodes(); }
329  virtual dof_id_type nActiveElem() const { return _mesh->n_active_elem(); }
330  virtual dof_id_type nActiveLocalElem() const { return _mesh->n_active_local_elem(); }
331  virtual SubdomainID nSubdomains() const { return _mesh->n_subdomains(); }
332  virtual unsigned int nPartitions() const { return _mesh->n_partitions(); }
333  virtual bool skipPartitioning() const { return _mesh->skip_partitioning(); }
334  virtual bool skipNoncriticalPartitioning() const;
335 
341  virtual dof_id_type maxNodeId() const;
342  virtual dof_id_type maxElemId() const;
343 
350  virtual const Node & node(const dof_id_type i) const;
351  virtual Node & node(const dof_id_type i);
352  virtual const Node & nodeRef(const dof_id_type i) const;
353  virtual Node & nodeRef(const dof_id_type i);
354  virtual const Node * nodePtr(const dof_id_type i) const;
355  virtual Node * nodePtr(const dof_id_type i);
356  virtual const Node * queryNodePtr(const dof_id_type i) const;
357  virtual Node * queryNodePtr(const dof_id_type i);
358 
365  virtual Elem * elem(const dof_id_type i);
366  virtual const Elem * elem(const dof_id_type i) const;
367  virtual Elem * elemPtr(const dof_id_type i);
368  virtual const Elem * elemPtr(const dof_id_type i) const;
369  virtual Elem * queryElemPtr(const dof_id_type i);
370  virtual const Elem * queryElemPtr(const dof_id_type i) const;
371 
375  bool prepared() const;
376  virtual void prepared(bool state);
377 
383  void needsPrepareForUse();
384 
390  void meshChanged();
391 
397  virtual void onMeshChanged();
398 
402  void cacheChangedLists();
403 
411 
419 
426  const std::vector<const Elem *> & coarsenedElementChildren(const Elem * elem) const;
427 
432  void updateActiveSemiLocalNodeRange(std::set<dof_id_type> & ghosted_elems);
433 
439  bool isSemiLocal(Node * const node) const;
440 
442 
455 
459  const std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> &
460  getBoundariesToElems() const;
461 
465  const std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> &
467 
473  std::unordered_set<dof_id_type> getBoundaryActiveSemiLocalElemIds(BoundaryID bid) const;
474 
484  std::unordered_set<dof_id_type> getBoundaryActiveNeighborElemIds(BoundaryID bid) const;
485 
494  const std::set<SubdomainID> & blk_group) const;
495 
500  const std::set<SubdomainID> & meshSubdomains() const;
501 
506  const std::set<BoundaryID> & meshBoundaryIds() const;
507 
512  const std::set<BoundaryID> & meshSidesetIds() const;
513 
518  const std::set<BoundaryID> & meshNodesetIds() const;
519 
524  void setBoundaryToNormalMap(std::unique_ptr<std::map<BoundaryID, RealVectorValue>> boundary_map);
525 
526  // DEPRECATED METHOD
527  void setBoundaryToNormalMap(std::map<BoundaryID, RealVectorValue> * boundary_map);
528 
533  void setMeshBoundaryIDs(std::set<BoundaryID> boundary_IDs);
534 
540 
557  bool prepare(const MeshBase * mesh_to_clone);
558 
562  void update();
563 
567  unsigned int uniformRefineLevel() const;
568 
572  void setUniformRefineLevel(unsigned int, bool deletion = true);
573 
585 
590 
594  void addGhostedBoundary(BoundaryID boundary_id);
595 
600  void setGhostedBoundaryInflation(const std::vector<Real> & inflation);
601 
605  const std::set<unsigned int> & getGhostedBoundaries() const;
606 
610  const std::vector<Real> & getGhostedBoundaryInflation() const;
611 
615  void ghostGhostedBoundaries();
616 
620  void needGhostGhostedBoundaries(bool needghost) { _need_ghost_ghosted_boundaries = needghost; }
621 
625  unsigned int getPatchSize() const;
626 
630  unsigned int getGhostingPatchSize() const { return _ghosting_patch_size; }
631 
635  unsigned int getMaxLeafSize() const { return _max_leaf_size; }
636 
640  void setPatchUpdateStrategy(Moose::PatchUpdateType patch_update_strategy);
641 
646 
653  libMesh::BoundingBox getInflatedProcessorBoundingBox(Real inflation_multiplier = 0.01) const;
654 
658  operator libMesh::MeshBase &();
659  operator const libMesh::MeshBase &() const;
660 
664  MeshBase & getMesh();
665  MeshBase & getMesh(const std::string & name);
666  const MeshBase & getMesh() const;
667  const MeshBase & getMesh(const std::string & name) const;
668  const MeshBase * getMeshPtr() const;
669 
673 #ifdef MOOSE_KOKKOS_ENABLED
674  const Moose::Kokkos::Mesh * getKokkosMesh() const { return _kokkos_mesh.get(); }
675 #endif
676 
680  void printInfo(std::ostream & os = libMesh::out, const unsigned int verbosity = 0) const;
681 
685  const std::set<SubdomainID> & getNodeBlockIds(const Node & node) const;
686 
691  const std::vector<dof_id_type> & getNodeList(boundary_id_type nodeset_id) const;
692 
698  const Node * addUniqueNode(const Point & p, Real tol = 1e-6);
699 
713  Node * addQuadratureNode(const Elem * elem,
714  const unsigned short int side,
715  const unsigned int qp,
716  BoundaryID bid,
717  const Point & point);
718 
726  Node * getQuadratureNode(const Elem * elem, const unsigned short int side, const unsigned int qp);
727 
732  void clearQuadratureNodes();
733 
740  BoundaryID getBoundaryID(const BoundaryName & boundary_name) const;
741 
748  std::vector<BoundaryID> getBoundaryIDs(const std::vector<BoundaryName> & boundary_name,
749  bool generate_unknown = false) const;
750 
757  SubdomainID getSubdomainID(const SubdomainName & subdomain_name) const;
758 
765  std::vector<SubdomainID>
766  getSubdomainIDs(const std::vector<SubdomainName> & subdomain_names) const;
767  std::set<SubdomainID> getSubdomainIDs(const std::set<SubdomainName> & subdomain_names) const;
768 
772  void setSubdomainName(SubdomainID subdomain_id, const SubdomainName & name);
773 
777  static void
778  setSubdomainName(MeshBase & mesh, SubdomainID subdomain_id, const SubdomainName & name);
779 
783  const std::string & getSubdomainName(SubdomainID subdomain_id) const;
784 
791  std::vector<SubdomainName>
792  getSubdomainNames(const std::vector<SubdomainID> & subdomain_ids) const;
793 
797  void setBoundaryName(BoundaryID boundary_id, BoundaryName name);
798 
802  const std::string & getBoundaryName(BoundaryID boundary_id);
803 
809  void buildPeriodicNodeMap(std::multimap<dof_id_type, dof_id_type> & periodic_node_map,
810  unsigned int var_number,
811  libMesh::PeriodicBoundaries * pbs) const;
812 
816  void buildPeriodicNodeSets(std::map<BoundaryID, std::set<dof_id_type>> & periodic_node_sets,
817  unsigned int var_number,
818  libMesh::PeriodicBoundaries * pbs) const;
819 
823  Real dimensionWidth(unsigned int component) const;
824 
826 
829  virtual Real getMinInDimension(unsigned int component) const;
830  virtual Real getMaxInDimension(unsigned int component) const;
832 
843  bool detectOrthogonalDimRanges(Real tol = 1e-6);
844 
849  void addPeriodicVariable(unsigned int var_num, BoundaryID primary, BoundaryID secondary);
850 
856  bool isTranslatedPeriodic(unsigned int nonlinear_var_num, unsigned int component) const;
857 
865  RealVectorValue minPeriodicVector(unsigned int nonlinear_var_num, Point p, Point q) const;
866 
874  Real minPeriodicDistance(unsigned int nonlinear_var_num, Point p, Point q) const;
875 
882  const std::pair<BoundaryID, BoundaryID> * getPairedBoundaryMapping(unsigned int component);
883 
891 
901  const std::vector<std::vector<QpMap>> &
902  getRefinementMap(const Elem & elem, int parent_side, int child, int child_side);
903 
911  const std::vector<std::pair<unsigned int, QpMap>> & getCoarseningMap(const Elem & elem,
912  int input_side);
913 
918  void
919  changeBoundaryId(const boundary_id_type old_id, const boundary_id_type new_id, bool delete_prev);
920 
925  static void changeBoundaryId(MeshBase & mesh,
926  const boundary_id_type old_id,
927  const boundary_id_type new_id,
928  bool delete_prev);
929 
936  const std::set<BoundaryID> & getSubdomainBoundaryIds(const SubdomainID subdomain_id) const;
937 
944  std::set<BoundaryID> getSubdomainInterfaceBoundaryIds(const SubdomainID subdomain_id) const;
945 
952  std::set<SubdomainID> getBoundaryConnectedBlocks(const BoundaryID bid) const;
953 
960  std::set<SubdomainID> getBoundaryConnectedSecondaryBlocks(const BoundaryID bid) const;
961 
968  std::set<SubdomainID> getInterfaceConnectedBlocks(const BoundaryID bid) const;
969 
976  const std::set<SubdomainID> & getBlockConnectedBlocks(const SubdomainID subdomain_id) const;
977 
981  bool isBoundaryNode(dof_id_type node_id) const;
982 
987  bool isBoundaryNode(dof_id_type node_id, BoundaryID bnd_id) const;
988 
992  bool isBoundaryElem(dof_id_type elem_id) const;
993 
998  bool isBoundaryElem(dof_id_type elem_id, BoundaryID bnd_id) const;
999 
1007  void errorIfDistributedMesh(std::string name) const;
1008 
1012  virtual bool isDistributedMesh() const { return _use_distributed_mesh; }
1013 
1018 
1022  void setParallelType(ParallelType parallel_type);
1023 
1028 
1029  /*
1030  * Set/Get the partitioner name
1031  */
1032  const MooseEnum & partitionerName() const { return _partitioner_name; }
1033 
1038 
1042  void allowRecovery(bool allow) { _allow_recovery = allow; }
1043 
1047  static void setPartitioner(MeshBase & mesh_base,
1048  MooseEnum & partitioner,
1049  bool use_distributed_mesh,
1050  const InputParameters & params,
1051  MooseObject & context_obj);
1052 
1056  void setCustomPartitioner(libMesh::Partitioner * partitioner);
1057 
1059 
1062  bool isCustomPartitionerRequested() const;
1063  void setIsCustomPartitionerRequested(bool cpr);
1065 
1068 
1070  bool hasSecondOrderElements();
1071 
1076  virtual std::unique_ptr<libMesh::PointLocatorBase> getPointLocator() const;
1077 
1082  virtual std::string getFileName() const { return ""; }
1083 
1085  using PeriodicNodeInfo = std::pair<const Node *, BoundaryID>;
1086 
1090  void needsRemoteElemDeletion(bool need_delete) { _need_delete = need_delete; }
1091 
1095  bool needsRemoteElemDeletion() const { return _need_delete; }
1096 
1100  void allowRemoteElementRemoval(bool allow_removal);
1101 
1106 
1110  void deleteRemoteElements();
1111 
1115  bool hasMeshBase() const { return _mesh.get() != nullptr; }
1116 
1120  bool hasElementID(const std::string & id_name) const;
1121 
1125  unsigned int getElementIDIndex(const std::string & id_name) const;
1126 
1130  dof_id_type maxElementID(unsigned int elem_id_index) const { return _max_ids[elem_id_index]; }
1131 
1135  dof_id_type minElementID(unsigned int elem_id_index) const { return _min_ids[elem_id_index]; }
1136 
1140  bool areElemIDsIdentical(const std::string & id_name1, const std::string & id_name2) const;
1141 
1145  std::set<dof_id_type> getAllElemIDs(unsigned int elem_id_index) const;
1146 
1151  std::set<dof_id_type> getElemIDsOnBlocks(unsigned int elem_id_index,
1152  const std::set<SubdomainID> & blks) const;
1153 
1157  unsigned int getMaxSidesPerElem() const { return _max_sides_per_elem; }
1158 
1162  unsigned int getMaxNodesPerElem() const { return _max_nodes_per_elem; }
1163 
1167  unsigned int getMaxNodesPerSide() const { return _max_nodes_per_side; }
1168 
1169  std::unordered_map<dof_id_type, std::set<dof_id_type>>
1170  getElemIDMapping(const std::string & from_id_name, const std::string & to_id_name) const;
1171 
1173  unsigned int nFace() const { return _face_info.size(); }
1174 
1176  const std::vector<const FaceInfo *> & faceInfo() const;
1177 
1179  struct face_info_iterator;
1180  struct const_face_info_iterator;
1181 
1186  face_info_iterator ownedFaceInfoBegin();
1187  face_info_iterator ownedFaceInfoEnd();
1188 
1190  struct elem_info_iterator;
1191  struct const_elem_info_iterator;
1192 
1197  elem_info_iterator ownedElemInfoBegin();
1198  elem_info_iterator ownedElemInfoEnd();
1199 
1201  const FaceInfo * faceInfo(const Elem * elem, unsigned int side) const;
1202 
1204  const ElemInfo & elemInfo(const dof_id_type id) const;
1205 
1207  const std::vector<const ElemInfo *> & elemInfoVector() const { return _elem_info; }
1208 
1210  const std::vector<FaceInfo> & allFaceInfo() const;
1212 
1216  void cacheFaceInfoVariableOwnership() const;
1217 
1222  void cacheFVElementalDoFs() const;
1223 
1229  void computeFiniteVolumeCoords() const;
1230 
1234  void isDisplaced(bool is_displaced) { _is_displaced = is_displaced; }
1235 
1239  bool isDisplaced() const { return _is_displaced; }
1240 
1244  const std::map<boundary_id_type, std::vector<dof_id_type>> & nodeSetNodes() const;
1245 
1251 
1257 
1261  const std::map<SubdomainID, Moose::CoordinateSystemType> & getCoordSystem() const;
1262 
1266  void setCoordSystem(const std::vector<SubdomainName> & blocks, const MultiMooseEnum & coord_sys);
1267 
1272  void setAxisymmetricCoordAxis(const MooseEnum & rz_coord_axis);
1273 
1288  void setGeneralAxisymmetricCoordAxes(const std::vector<SubdomainName> & blocks,
1289  const std::vector<std::pair<Point, RealVectorValue>> & axes);
1290 
1296  const std::pair<Point, RealVectorValue> &
1297  getGeneralAxisymmetricCoordAxis(SubdomainID subdomain_id) const;
1298 
1302  bool usingGeneralAxisymmetricCoordAxes() const;
1303 
1308  unsigned int getAxisymmetricRadialCoord() const;
1309 
1315  void checkCoordinateSystems();
1316 
1320  void setCoordData(const MooseMesh & other_mesh);
1321 
1326 
1331 
1337 
1341  const MooseUnits & lengthUnit() const;
1342 
1347  const std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *> &
1348  getLowerDElemMap() const;
1349 
1353  bool isSplit() const { return _is_split; }
1354 
1361  void buildFiniteVolumeInfo() const;
1362 
1368  void setupFiniteVolumeMeshData() const;
1369 
1373  void doingPRefinement(bool doing_p_refinement) { _doing_p_refinement = doing_p_refinement; }
1374 
1378  [[nodiscard]] bool doingPRefinement() const { return _doing_p_refinement; }
1379 
1383  unsigned int maxPLevel() const { return _max_p_level; }
1384 
1388  unsigned int maxHLevel() const { return _max_h_level; }
1389 
1394  const std::vector<QpMap> & getPRefinementMap(const Elem & elem) const;
1399  const std::vector<QpMap> & getPRefinementSideMap(const Elem & elem) const;
1404  const std::vector<QpMap> & getPCoarseningMap(const Elem & elem) const;
1409  const std::vector<QpMap> & getPCoarseningSideMap(const Elem & elem) const;
1410 
1412 
1418  bool isLowerD(const SubdomainID subdomain_id) const;
1419 
1424  bool hasLowerD() const { return _has_lower_d; }
1425 
1429  const std::set<SubdomainID> & interiorLowerDBlocks() const { return _lower_d_interior_blocks; }
1433  const std::set<SubdomainID> & boundaryLowerDBlocks() const { return _lower_d_boundary_blocks; }
1436 
1437 protected:
1439  std::vector<std::unique_ptr<libMesh::GhostingFunctor>> _ghosting_functors;
1440 
1442  std::vector<std::shared_ptr<RelationshipManager>> _relationship_managers;
1443 
1446 
1450 
1457 
1459  std::unique_ptr<libMesh::MeshBase> _mesh;
1460 
1462 #ifdef MOOSE_KOKKOS_ENABLED
1463  std::unique_ptr<Moose::Kokkos::Mesh> _kokkos_mesh;
1464 #endif
1465 
1469 
1471  std::unique_ptr<libMesh::Partitioner> _custom_partitioner;
1473 
1475  enum
1476  {
1477  X = 0,
1478  Y,
1480  };
1481  enum
1482  {
1483  MIN = 0,
1485  };
1486 
1489 
1492 
1495 
1498 
1501 
1503  bool _moose_mesh_prepared = false;
1504 
1506  std::unique_ptr<ConstElemPointerRange> _refined_elements;
1507 
1509  std::unique_ptr<ConstElemPointerRange> _coarsened_elements;
1510 
1516  std::map<const Elem *, std::vector<const Elem *>> _coarsened_element_children;
1517 
1519  std::set<Node *> _semilocal_node_list;
1520 
1525  std::unique_ptr<libMesh::ConstElemRange> _active_local_elem_range;
1526 
1527  std::unique_ptr<SemiLocalNodeRange> _active_semilocal_node_range;
1528  std::unique_ptr<libMesh::NodeRange> _active_node_range;
1529  std::unique_ptr<libMesh::ConstNodeRange> _local_node_range;
1530  std::unique_ptr<libMesh::StoredRange<MooseMesh::const_bnd_node_iterator, const BndNode *>>
1532  std::unique_ptr<libMesh::StoredRange<MooseMesh::const_bnd_elem_iterator, const BndElement *>>
1534 
1536  std::map<dof_id_type, std::vector<dof_id_type>> _node_to_elem_map;
1538 
1540  std::map<dof_id_type, std::vector<dof_id_type>> _node_to_active_semilocal_elem_map;
1542 
1547  std::set<SubdomainID> _mesh_subdomains;
1548 
1550 
1555  std::set<BoundaryID> _mesh_boundary_ids;
1556  std::set<BoundaryID> _mesh_sideset_ids;
1557  std::set<BoundaryID> _mesh_nodeset_ids;
1559 
1561  std::unique_ptr<std::map<BoundaryID, RealVectorValue>> _boundary_to_normal_map;
1562 
1564  std::vector<BndNode *> _bnd_nodes;
1565  typedef std::vector<BndNode *>::iterator bnd_node_iterator_imp;
1566  typedef std::vector<BndNode *>::const_iterator const_bnd_node_iterator_imp;
1568  std::map<boundary_id_type, std::set<dof_id_type>> _bnd_node_ids;
1569 
1571  std::vector<BndElement *> _bnd_elems;
1572  typedef std::vector<BndElement *>::iterator bnd_elem_iterator_imp;
1573  typedef std::vector<BndElement *>::const_iterator const_bnd_elem_iterator_imp;
1574 
1576  std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> _bnd_elem_ids;
1577 
1578  std::map<dof_id_type, Node *> _quadrature_nodes;
1579  std::map<dof_id_type, std::map<unsigned int, std::map<dof_id_type, Node *>>>
1581  std::vector<BndNode> _extra_bnd_nodes;
1582 
1584  std::map<dof_id_type, std::set<SubdomainID>> _block_node_list;
1585 
1587  std::map<boundary_id_type, std::vector<dof_id_type>> _node_set_nodes;
1588 
1589  std::set<unsigned int> _ghosted_boundaries;
1591 
1593  unsigned int _patch_size;
1594 
1596  unsigned int _ghosting_patch_size;
1597 
1598  // The maximum number of points in each leaf of the KDTree used in the nearest neighbor search.
1599  unsigned int _max_leaf_size;
1600 
1603 
1605  std::vector<Node *> _node_map;
1606 
1609 
1611  std::vector<std::vector<Real>> _bounds;
1612 
1614  std::vector<std::pair<BoundaryID, BoundaryID>> _paired_boundary;
1615 
1617  const bool _is_split;
1618 
1619  void cacheInfo();
1620  void freeBndNodes();
1621  void freeBndElems();
1622  void setPartitionerHelper(MeshBase * mesh = nullptr);
1623 
1624 private:
1627  mutable std::unordered_map<dof_id_type, ElemInfo> _elem_to_elem_info;
1628 
1631  mutable std::vector<const ElemInfo *> _elem_info;
1632 
1635  mutable std::vector<FaceInfo> _all_face_info;
1636 
1639  mutable std::vector<const FaceInfo *> _face_info;
1640 
1642  mutable std::unordered_map<std::pair<const Elem *, unsigned int>, FaceInfo *>
1644 
1645  // true if the _face_info member needs to be rebuilt/updated.
1646  mutable bool _finite_volume_info_dirty = true;
1647 
1648  // True if we have cached elemental dofs ids for the linear finite volume variables.
1649  // This happens in the first system which has a linear finite volume variable, considering
1650  // that currently we only support one variable per linear system.
1652 
1657  std::map<unsigned int, std::vector<bool>> _periodic_dim;
1658 
1663 
1665  std::vector<Node *> _extreme_nodes;
1666 
1672  void detectPairedSidesets();
1673 
1685  void buildRefinementMap(const Elem & elem,
1686  libMesh::QBase & qrule,
1687  libMesh::QBase & qrule_face,
1688  int parent_side,
1689  int child,
1690  int child_side);
1691 
1701  void buildCoarseningMap(const Elem & elem,
1702  libMesh::QBase & qrule,
1703  libMesh::QBase & qrule_face,
1704  int input_side);
1705 
1714  void mapPoints(const std::vector<Point> & from,
1715  const std::vector<Point> & to,
1716  std::vector<QpMap> & qp_map);
1717 
1740  void findAdaptivityQpMaps(const Elem * template_elem,
1741  libMesh::QBase & qrule,
1742  libMesh::QBase & qrule_face,
1743  std::vector<std::vector<QpMap>> & refinement_map,
1744  std::vector<std::pair<unsigned int, QpMap>> & coarsen_map,
1745  int parent_side,
1746  int child,
1747  int child_side);
1748 
1750 
1751  const std::vector<QpMap> & getPRefinementMapHelper(
1752  const Elem & elem,
1753  const std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>> &) const;
1754  const std::vector<QpMap> & getPCoarseningMapHelper(
1755  const Elem & elem,
1756  const std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>> &) const;
1757 
1762  void updateCoordTransform();
1763 
1769 
1783  std::map<std::pair<int, libMesh::ElemType>, std::vector<std::vector<QpMap>>>
1785 
1786  std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1788  std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1790 
1792  std::map<libMesh::ElemType, std::map<std::pair<int, int>, std::vector<std::vector<QpMap>>>>
1794 
1809  std::map<std::pair<int, libMesh::ElemType>, std::vector<std::pair<unsigned int, QpMap>>>
1811 
1812  std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1814  std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1816 
1818  {
1820  std::set<SubdomainID> neighbor_subs;
1821 
1824  std::set<BoundaryID> boundary_ids;
1825 
1828  };
1829 
1831  std::unordered_map<SubdomainID, SubdomainData> _sub_to_data;
1832 
1834  std::unordered_map<SubdomainID, std::set<BoundaryID>> _neighbor_subdomain_boundary_ids;
1835 
1837  std::set<SubdomainID> _lower_d_interior_blocks;
1839  std::set<SubdomainID> _lower_d_boundary_blocks;
1841  std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *>
1843  std::unordered_map<const Elem *, unsigned short int> _lower_d_elem_to_higher_d_elem_side;
1844 
1848 
1851 
1854 
1857 
1860 
1863 
1869 
1871  std::vector<std::unordered_map<SubdomainID, std::set<dof_id_type>>> _block_id_mapping;
1873  std::vector<dof_id_type> _max_ids;
1875  std::vector<dof_id_type> _min_ids;
1877  std::vector<std::vector<bool>> _id_identical_flag;
1878 
1880  unsigned int _max_sides_per_elem;
1881 
1883  unsigned int _max_nodes_per_elem;
1884 
1886  unsigned int _max_nodes_per_side;
1887 
1889  void computeMaxPerElemAndSide();
1890 
1893 
1895  void buildElemIDInfo();
1896 
1898  void buildLowerDMesh();
1899 
1901  std::map<SubdomainID, Moose::CoordinateSystemType> & _coord_sys;
1902 
1904  unsigned int _rz_coord_axis;
1905 
1907  std::unordered_map<SubdomainID, std::pair<Point, RealVectorValue>> _subdomain_id_to_rz_coord_axis;
1908 
1911  std::unique_ptr<MooseAppCoordTransform> _coord_transform;
1912 
1915 
1917  std::vector<SubdomainName> _provided_coord_blocks;
1918 
1922  unsigned int _max_p_level;
1924  unsigned int _max_h_level;
1925 
1926  template <typename T>
1927  struct MeshType;
1928 };
1929 
1930 inline MooseAppCoordTransform &
1932 {
1933  mooseAssert(_coord_transform, "The coordinate transformation object is null.");
1934  return *_coord_transform;
1935 }
1936 
1937 template <>
1938 struct MooseMesh::MeshType<libMesh::ReplicatedMesh>
1939 {
1941 };
1942 
1943 template <>
1944 struct MooseMesh::MeshType<libMesh::DistributedMesh>
1945 {
1947 };
1948 
1953  : variant_filter_iterator<MeshBase::Predicate, const FaceInfo *>
1954 {
1955  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
1956  template <typename PredType, typename IterType>
1957  face_info_iterator(const IterType & d, const IterType & e, const PredType & p)
1958  : variant_filter_iterator<MeshBase::Predicate, const FaceInfo *>(d, e, p)
1959  {
1960  }
1961 };
1962 
1968  const FaceInfo * const,
1969  const FaceInfo * const &,
1970  const FaceInfo * const *>
1971 {
1972  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
1973  template <typename PredType, typename IterType>
1974  const_face_info_iterator(const IterType & d, const IterType & e, const PredType & p)
1975  : variant_filter_iterator<MeshBase::Predicate,
1976  const FaceInfo * const,
1977  const FaceInfo * const &,
1978  const FaceInfo * const *>(d, e, p)
1979  {
1980  }
1981 
1982  // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
1983  // variant_filter_iterator copy constructor. Note that this one is *not* templated!
1985  : variant_filter_iterator<MeshBase::Predicate,
1986  const FaceInfo * const,
1987  const FaceInfo * const &,
1988  const FaceInfo * const *>(rhs)
1989  {
1990  }
1991 };
1992 
1997  : variant_filter_iterator<MeshBase::Predicate, const ElemInfo *>
1998 {
1999  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2000  template <typename PredType, typename IterType>
2001  elem_info_iterator(const IterType & d, const IterType & e, const PredType & p)
2002  : variant_filter_iterator<MeshBase::Predicate, const ElemInfo *>(d, e, p)
2003  {
2004  }
2005 };
2006 
2012  const ElemInfo * const,
2013  const ElemInfo * const &,
2014  const ElemInfo * const *>
2015 {
2016  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2017  template <typename PredType, typename IterType>
2018  const_elem_info_iterator(const IterType & d, const IterType & e, const PredType & p)
2019  : variant_filter_iterator<MeshBase::Predicate,
2020  const ElemInfo * const,
2021  const ElemInfo * const &,
2022  const ElemInfo * const *>(d, e, p)
2023  {
2024  }
2025 
2026  // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2027  // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2029  : variant_filter_iterator<MeshBase::Predicate,
2030  const ElemInfo * const,
2031  const ElemInfo * const &,
2032  const ElemInfo * const *>(rhs)
2033  {
2034  }
2035 };
2036 
2040 struct MooseMesh::bnd_node_iterator : variant_filter_iterator<MeshBase::Predicate, BndNode *>
2041 {
2042  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2043  template <typename PredType, typename IterType>
2044  bnd_node_iterator(const IterType & d, const IterType & e, const PredType & p)
2045  : variant_filter_iterator<MeshBase::Predicate, BndNode *>(d, e, p)
2046  {
2047  }
2048 };
2049 
2055  BndNode * const,
2056  BndNode * const &,
2057  BndNode * const *>
2058 {
2059  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2060  template <typename PredType, typename IterType>
2061  const_bnd_node_iterator(const IterType & d, const IterType & e, const PredType & p)
2062  : variant_filter_iterator<MeshBase::Predicate,
2063  BndNode * const,
2064  BndNode * const &,
2065  BndNode * const *>(d, e, p)
2066  {
2067  }
2068 
2069  // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2070  // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2072  : variant_filter_iterator<MeshBase::Predicate,
2073  BndNode * const,
2074  BndNode * const &,
2075  BndNode * const *>(rhs)
2076  {
2077  }
2078 };
2079 
2083 struct MooseMesh::bnd_elem_iterator : variant_filter_iterator<MeshBase::Predicate, BndElement *>
2084 {
2085  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2086  template <typename PredType, typename IterType>
2087  bnd_elem_iterator(const IterType & d, const IterType & e, const PredType & p)
2088  : variant_filter_iterator<MeshBase::Predicate, BndElement *>(d, e, p)
2089  {
2090  }
2091 };
2092 
2098  BndElement * const,
2099  BndElement * const &,
2100  BndElement * const *>
2101 {
2102  // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2103  template <typename PredType, typename IterType>
2104  const_bnd_elem_iterator(const IterType & d, const IterType & e, const PredType & p)
2105  : variant_filter_iterator<MeshBase::Predicate,
2106  BndElement * const,
2107  BndElement * const &,
2108  BndElement * const *>(d, e, p)
2109  {
2110  }
2111 
2112  // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2113  // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2115  : variant_filter_iterator<MeshBase::Predicate,
2116  BndElement * const,
2117  BndElement * const &,
2118  BndElement * const *>(rhs)
2119  {
2120  }
2121 };
2122 
2130 
2131 template <typename T>
2132 std::unique_ptr<T>
2134 {
2135  // If the requested mesh type to build doesn't match our current value for _use_distributed_mesh,
2136  // then we need to make sure to make our state consistent because other objects, like the periodic
2137  // boundary condition action, will be querying isDistributedMesh()
2138  if (_use_distributed_mesh != std::is_same<T, libMesh::DistributedMesh>::value)
2139  {
2140  if (getMeshPtr())
2141  mooseError("A MooseMesh object is being asked to build a libMesh mesh that is a different "
2142  "parallel type than the libMesh mesh that it wraps. This is not allowed. Please "
2143  "create another MooseMesh object to wrap the new libMesh mesh");
2145  }
2146 
2147  if (dim == libMesh::invalid_uint)
2148  {
2149  if (isParamValid("dim"))
2150  dim = getParam<MooseEnum>("dim");
2151  else
2152  // Legacy selection of the default for the 'dim' parameter
2153  dim = 1;
2154  }
2155 
2156  auto mesh = std::make_unique<T>(_communicator, dim);
2157 
2158  if (!getParam<bool>("allow_renumbering"))
2159  mesh->allow_renumbering(false);
2160 
2161  mesh->allow_remote_element_removal(_allow_remote_element_removal);
2163 
2165  {
2166  // Check of partitioner is supplied (not allowed if custom partitioner is used)
2167  if (!parameters().isParamSetByAddParam("partitioner"))
2168  mooseError("If partitioner block is provided, partitioner keyword cannot be used!");
2169  // Set custom partitioner
2170  if (!_custom_partitioner.get())
2171  mooseError("Custom partitioner requested but not set!");
2172  mesh->partitioner() = _custom_partitioner->clone();
2173  }
2174  else
2175  setPartitionerHelper(mesh.get());
2176 
2177  return mesh;
2178 }
2179 
2180 inline bool
2182 {
2184 }
2185 
2186 inline void
2188 {
2189  _parallel_type = parallel_type;
2191 }
2192 
2193 inline bool
2194 MooseMesh::hasElementID(const std::string & id_name) const
2195 {
2196  return getMesh().has_elem_integer(id_name);
2197 }
2198 
2199 inline unsigned int
2200 MooseMesh::getElementIDIndex(const std::string & id_name) const
2201 {
2202  if (!hasElementID(id_name))
2203  mooseError("Mesh does not have element ID for ", id_name);
2204  return getMesh().get_elem_integer_index(id_name);
2205 }
2206 
2207 inline bool
2208 MooseMesh::areElemIDsIdentical(const std::string & id_name1, const std::string & id_name2) const
2209 {
2210  auto id1 = getElementIDIndex(id_name1);
2211  auto id2 = getElementIDIndex(id_name2);
2212  return _id_identical_flag[id1][id2];
2213 }
2214 
2215 inline const std::vector<const FaceInfo *> &
2217 {
2218  return _face_info;
2219 }
2220 
2221 inline const std::vector<FaceInfo> &
2223 {
2224  return _all_face_info;
2225 }
2226 
2227 inline const std::map<boundary_id_type, std::vector<dof_id_type>> &
2229 {
2230  return _node_set_nodes;
2231 }
2232 
2233 inline const std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *> &
2235 {
2237 }
2238 
2239 inline bool
2240 MooseMesh::isLowerD(const SubdomainID subdomain_id) const
2241 {
2242  return libmesh_map_find(_sub_to_data, subdomain_id).is_lower_d;
2243 }
ParallelType _parallel_type
Can be set to DISTRIBUTED, REPLICATED, or DEFAULT.
Definition: MooseMesh.h:1449
static InputParameters validParams()
Typical "Moose-style" constructor and copy constructor.
Definition: MooseMesh.C:84
virtual std::string getFileName() const
Returns the name of the mesh file read to produce this mesh if any or an empty string otherwise...
Definition: MooseMesh.h:1082
virtual bnd_node_iterator bndNodesEnd()
Definition: MooseMesh.C:1591
virtual bnd_elem_iterator bndElemsEnd()
Definition: MooseMesh.C:1607
std::vector< std::vector< Real > > _bounds
The bounds in each dimension of the mesh for regular orthogonal meshes.
Definition: MooseMesh.h:1611
const std::map< boundary_id_type, std::vector< dof_id_type > > & nodeSetNodes() const
Definition: MooseMesh.h:2228
std::set< Node * > _semilocal_node_list
Used for generating the semilocal node range.
Definition: MooseMesh.h:1519
std::map< dof_id_type, Node * > _quadrature_nodes
Definition: MooseMesh.h:1578
const std::vector< QpMap > & getPCoarseningSideMap(const Elem &elem) const
Get the map describing for each side quadrature point (qp) on the coarse level which qp on the previo...
Definition: MooseMesh.C:4436
virtual Real getMaxInDimension(unsigned int component) const
Definition: MooseMesh.C:2246
virtual unsigned int spatialDimension() const
Returns MeshBase::spatial_dimension.
Definition: MooseMesh.h:186
bool _node_to_elem_map_built
Definition: MooseMesh.h:1537
std::unique_ptr< libMesh::NodeRange > _active_node_range
Definition: MooseMesh.h:1528
std::vector< Node * > _extreme_nodes
A vector containing the nodes at the corners of a regular orthogonal mesh.
Definition: MooseMesh.h:1665
Node * addQuadratureNode(const Elem *elem, const unsigned short int side, const unsigned int qp, BoundaryID bid, const Point &point)
Adds a fictitious "QuadratureNode".
Definition: MooseMesh.C:1648
std::vector< std::tuple< dof_id_type, unsigned short int, boundary_id_type > > buildSideList()
As above, but uses the non-deprecated std::tuple interface.
Definition: MooseMesh.C:3049
const_bnd_elem_iterator(const bnd_elem_iterator &rhs)
Definition: MooseMesh.h:2114
bool isFiniteVolumeInfoDirty() const
Definition: MooseMesh.h:1330
void computeMaxPerElemAndSide()
Compute the maximum numbers per element and side.
Definition: MooseMesh.C:1074
const std::set< BoundaryID > & meshNodesetIds() const
Returns a read-only reference to the set of nodesets currently present in the Mesh.
Definition: MooseMesh.C:3229
void buildElemIDInfo()
Build extra data for faster access to the information of extra element integers.
Definition: MooseMesh.C:1097
std::vector< FaceInfo > _all_face_info
FaceInfo object storing information for face based loops.
Definition: MooseMesh.h:1635
std::vector< const FaceInfo * > _face_info
Holds only those FaceInfo objects that have processor_id equal to this process&#39;s id, e.g.
Definition: MooseMesh.h:1639
bool allowRemoteElementRemoval() const
Whether we are allow remote element removal.
Definition: MooseMesh.h:1105
bool getConstructNodeListFromSideList()
Return construct node list from side list boolean.
Definition: MooseMesh.h:1435
virtual Real getMinInDimension(unsigned int component) const
Returns the min or max of the requested dimension respectively.
Definition: MooseMesh.C:2237
libMesh::ConstElemRange * getActiveLocalElementRange()
Return pointers to range objects for various types of ranges (local nodes, boundary elems...
Definition: MooseMesh.C:1276
The Kokkos mesh object.
Definition: KokkosMesh.h:54
virtual void onMeshChanged()
Declares a callback function that is executed at the conclusion of meshChanged(). ...
Definition: MooseMesh.C:923
bool prepared() const
Setter/getter for whether the mesh is prepared.
Definition: MooseMesh.C:3177
void needsPrepareForUse()
If this method is called, we will call libMesh&#39;s prepare_for_use method when we call Moose&#39;s prepare ...
Definition: MooseMesh.C:3205
const std::set< BoundaryID > & getBoundaryIDs() const
Returns a const reference to a set of all user-specified boundary IDs.
Definition: MooseMesh.C:3017
bool _is_nemesis
True if a Nemesis Mesh was read in.
Definition: MooseMesh.h:1500
std::vector< SubdomainName > _provided_coord_blocks
Set for holding user-provided coordinate system type block names.
Definition: MooseMesh.h:1917
virtual MooseMesh & clone() const
Clone method.
Definition: MooseMesh.C:2860
void allowRecovery(bool allow)
Set whether or not this mesh is allowed to read a recovery file.
Definition: MooseMesh.h:1042
const std::vector< QpMap > & getPCoarseningMapHelper(const Elem &elem, const std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap >> &) const
Definition: MooseMesh.C:4408
const std::string LIST_GEOM_ELEM
Definition: MooseMesh.h:62
bool isCustomPartitionerRequested() const
Setter and getter for _custom_partitioner_requested.
Definition: MooseMesh.C:3745
bool _need_ghost_ghosted_boundaries
A parallel mesh generator such as DistributedRectilinearMeshGenerator already make everything ready...
Definition: MooseMesh.h:1868
A class for creating restricted objects.
Definition: Restartable.h:28
bool isDisplaced() const
whether this mesh is a displaced mesh
Definition: MooseMesh.h:1239
unsigned int _uniform_refine_level
The level of uniform refinement requested (set to zero if AMR is disabled)
Definition: MooseMesh.h:1488
const std::vector< QpMap > & getPRefinementMapHelper(const Elem &elem, const std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap >> &) const
Definition: MooseMesh.C:4397
std::vector< dof_id_type > _min_ids
Minimum integer ID for each extra element integer.
Definition: MooseMesh.h:1875
std::unordered_set< dof_id_type > getBoundaryActiveSemiLocalElemIds(BoundaryID bid) const
Return all ids of elements which have a side which is part of a sideset.
Definition: MooseMesh.C:1369
const std::set< SubdomainID > & interiorLowerDBlocks() const
Definition: MooseMesh.h:1429
const MooseUnits & lengthUnit() const
Definition: MooseMesh.C:4370
void checkDuplicateSubdomainNames()
Loop through all subdomain IDs and check if there is name duplication used for the subdomains with sa...
Definition: MooseMesh.C:4377
const unsigned int invalid_uint
const std::set< BoundaryID > & getSubdomainBoundaryIds(const SubdomainID subdomain_id) const
Get the list of boundary ids associated with the given subdomain id.
Definition: MooseMesh.C:3537
unsigned int _max_leaf_size
Definition: MooseMesh.h:1599
The definition of the const_bnd_elem_iterator struct.
Definition: MooseMesh.h:2097
RealVectorValue _half_range
A convenience vector used to hold values in each dimension representing half of the range...
Definition: MooseMesh.h:1662
Keeps track of stuff related to assembling.
Definition: Assembly.h:109
const_bnd_elem_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2104
std::map< libMesh::ElemType, std::map< std::pair< int, int >, std::vector< std::vector< QpMap > > > > _elem_type_to_child_side_refinement_map
Holds mappings for "internal" child sides to parent volume. The second key is (child, child_side).
Definition: MooseMesh.h:1793
void setCoordData(const MooseMesh &other_mesh)
Set the coordinate system data to that of other_mesh.
Definition: MooseMesh.C:4362
virtual ~MooseMesh()
Definition: MooseMesh.C:378
void freeBndElems()
Definition: MooseMesh.C:405
const std::unordered_map< std::pair< const Elem *, unsigned short int >, const Elem * > & getLowerDElemMap() const
This function attempts to return the map from a high-order element side to its corresponding lower-d ...
Definition: MooseMesh.h:2234
bool _finite_volume_info_dirty
Definition: MooseMesh.h:1646
virtual Elem * elemPtr(const dof_id_type i)
Definition: MooseMesh.C:3153
elem_info_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2001
const_bnd_node_iterator(const MooseMesh::bnd_node_iterator &rhs)
Definition: MooseMesh.h:2071
char ** blocks
The definition of the bnd_elem_iterator struct.
Definition: MooseMesh.h:2083
std::map< SubdomainID, Moose::CoordinateSystemType > & _coord_sys
Type of coordinate system per subdomain.
Definition: MooseMesh.h:1901
const std::vector< const ElemInfo * > & elemInfoVector() const
Accessor for the element info objects owned by this process.
Definition: MooseMesh.h:1207
bool isBoundaryNode(dof_id_type node_id) const
Returns true if the requested node is in the list of boundary nodes, false otherwise.
Definition: MooseMesh.C:3605
face_info_iterator ownedFaceInfoBegin()
Iterators to owned faceInfo objects.
Definition: MooseMesh.C:1548
const std::map< dof_id_type, std::vector< dof_id_type > > & nodeToActiveSemilocalElemMap()
If not already created, creates a map from every node to all active semilocal elements to which they ...
Definition: MooseMesh.C:1244
std::vector< std::tuple< dof_id_type, unsigned short int, boundary_id_type > > buildActiveSideList() const
Calls BoundaryInfo::build_active_side_list.
Definition: MooseMesh.C:3055
static void setPartitioner(MeshBase &mesh_base, MooseEnum &partitioner, bool use_distributed_mesh, const InputParameters &params, MooseObject &context_obj)
Method for setting the partitioner on the passed in mesh_base object.
Definition: MooseMesh.C:3680
void buildLowerDMesh()
Build lower-d mesh for all sides.
Definition: MooseMesh.C:673
const std::set< BoundaryID > & meshSidesetIds() const
Returns a read-only reference to the set of sidesets currently present in the Mesh.
Definition: MooseMesh.C:3223
void setParallelType(ParallelType parallel_type)
Allow to change parallel type.
Definition: MooseMesh.h:2187
unsigned int getMaxNodesPerSide() const
Get the maximum number of nodes per side.
Definition: MooseMesh.h:1167
std::unordered_map< const Elem *, unsigned short int > _lower_d_elem_to_higher_d_elem_side
Definition: MooseMesh.h:1843
void cacheFVElementalDoFs() const
Cache the DoF indices for FV variables on each element.
Definition: MooseMesh.C:4059
void cacheFaceInfoVariableOwnership() const
Cache if variables live on the elements connected by the FaceInfo objects.
Definition: MooseMesh.C:3981
bool _custom_partitioner_requested
Definition: MooseMesh.h:1472
unsigned int _max_nodes_per_side
The maximum number of nodes per side.
Definition: MooseMesh.h:1886
const std::unordered_map< boundary_id_type, std::unordered_set< dof_id_type > > & getBoundariesToElems() const
Returns a map of boundaries to ids of elements on the boundary.
Definition: MooseMesh.C:1355
Moose::CoordinateSystemType getUniqueCoordSystem() const
Get the coordinate system from the mesh, it must be the same in all subdomains otherwise this will er...
Definition: MooseMesh.C:4225
std::map< dof_id_type, std::vector< dof_id_type > > _node_to_elem_map
A map of all of the current nodes to the elements that they are connected to.
Definition: MooseMesh.h:1536
MooseMesh()=delete
unsigned int _max_sides_per_elem
The maximum number of sides per element.
Definition: MooseMesh.h:1880
const_face_info_iterator(const MooseMesh::face_info_iterator &rhs)
Definition: MooseMesh.h:1984
const_bnd_node_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2061
bool isSemiLocal(Node *const node) const
Returns true if the node is semi-local.
Definition: MooseMesh.C:1007
libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement * > ConstBndElemRange
Definition: MooseMesh.h:2129
const Elem * getLowerDElem(const Elem *, unsigned short int) const
Returns a const pointer to a lower dimensional element that corresponds to a side of a higher dimensi...
Definition: MooseMesh.C:1739
std::unique_ptr< libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement * > > _bnd_elem_range
Definition: MooseMesh.h:1533
RealVectorValue minPeriodicVector(unsigned int nonlinear_var_num, Point p, Point q) const
This function returns the minimum vector between two points on the mesh taking into account periodici...
Definition: MooseMesh.C:2300
const std::vector< std::pair< unsigned int, QpMap > > & getCoarseningMap(const Elem &elem, int input_side)
Get the coarsening map for a given element type.
Definition: MooseMesh.C:2617
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
std::unordered_map< boundary_id_type, std::unordered_set< dof_id_type > > _bnd_elem_ids
Map of set of elem IDs connected to each boundary.
Definition: MooseMesh.h:1576
bool _is_changed
true if mesh is changed (i.e. after adaptivity step)
Definition: MooseMesh.h:1497
bool _doing_p_refinement
Whether we have p-refinement (as opposed to h-refinement)
Definition: MooseMesh.h:1920
void determineUseDistributedMesh()
Determine whether to use a distributed mesh.
Definition: MooseMesh.C:2866
const std::vector< std::vector< QpMap > > & getRefinementMap(const Elem &elem, int parent_side, int child, int child_side)
Get the refinement map for a given element type.
Definition: MooseMesh.C:2553
const std::string & getBoundaryName(BoundaryID boundary_id)
Return the name of the boundary given the id.
Definition: MooseMesh.C:1830
void cacheChangedLists()
Cache information about what elements were refined and coarsened in the previous step.
Definition: MooseMesh.C:928
bool isTranslatedPeriodic(unsigned int nonlinear_var_num, unsigned int component) const
Returns whether this generated mesh is periodic in the given dimension for the given variable...
Definition: MooseMesh.C:2289
std::set< SubdomainID > neighbor_subs
Neighboring subdomain ids.
Definition: MooseMesh.h:1820
unsigned int getGhostingPatchSize() const
Getter for the ghosting_patch_size parameter.
Definition: MooseMesh.h:630
std::vector< BndNode * >::iterator bnd_node_iterator_imp
Definition: MooseMesh.h:1565
const_face_info_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:1974
The definition of the bnd_node_iterator struct.
Definition: MooseMesh.h:2040
face_info_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:1957
std::vector< const ElemInfo * > _elem_info
Holds only those ElemInfo objects that have processor_id equal to this process&#39;s id, e.g.
Definition: MooseMesh.h:1631
MeshBase & mesh
void buildHRefinementAndCoarseningMaps(Assembly *assembly)
Definition: MooseMesh.C:2349
bool usingGeneralAxisymmetricCoordAxes() const
Returns true if general axisymmetric coordinate axes are being used.
Definition: MooseMesh.C:4319
std::map< std::pair< int, libMesh::ElemType >, std::vector< std::vector< QpMap > > > _elem_type_to_refinement_map
Holds mappings for volume to volume and parent side to child side Map key:
Definition: MooseMesh.h:1784
const std::set< SubdomainID > & getBlockConnectedBlocks(const SubdomainID subdomain_id) const
Get the list of subdomains neighboring a given subdomain.
Definition: MooseMesh.C:3594
virtual const Node * queryNodePtr(const dof_id_type i) const
Definition: MooseMesh.C:875
libMesh::StoredRange< std::set< Node * >::iterator, Node * > SemiLocalNodeRange
Definition: MooseMesh.h:59
std::unordered_map< std::pair< const Elem *, unsigned short int >, const Elem * > _higher_d_elem_side_to_lower_d_elem
Holds a map from a high-order element side to its corresponding lower-d element.
Definition: MooseMesh.h:1842
Helper object for holding qp mapping info.
Definition: MooseMesh.h:73
bool _has_lower_d
Whether there are any lower-dimensional blocks that are manifolds of higher-dimensional block faces...
Definition: MooseMesh.h:1847
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const std::vector< Real > & getGhostedBoundaryInflation() const
Return a writable reference to the _ghosted_boundaries_inflation vector.
Definition: MooseMesh.C:3287
std::unique_ptr< ConstElemPointerRange > _refined_elements
The elements that were just refined.
Definition: MooseMesh.h:1506
dof_id_type maxElementID(unsigned int elem_id_index) const
Return the maximum element ID for an extra element integer with its accessing index.
Definition: MooseMesh.h:1130
std::vector< std::vector< bool > > _id_identical_flag
Flags to indicate whether or not any two extra element integers are the same.
Definition: MooseMesh.h:1877
virtual dof_id_type maxElemId() const
Definition: MooseMesh.C:3133
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:159
virtual bnd_elem_iterator bndElemsBegin()
Return iterators to the beginning/end of the boundary elements list.
Definition: MooseMesh.C:1599
The definition of the const_bnd_node_iterator struct.
Definition: MooseMesh.h:2054
void setUniformRefineLevel(unsigned int, bool deletion=true)
Set uniform refinement level.
Definition: MooseMesh.C:3262
MooseEnum _partitioner_name
The partitioner used on this mesh.
Definition: MooseMesh.h:1467
std::unique_ptr< libMesh::ConstElemRange > _active_local_elem_range
A range for use with threading.
Definition: MooseMesh.h:1525
std::map< dof_id_type, std::set< SubdomainID > > _block_node_list
list of nodes that belongs to a specified block (domain)
Definition: MooseMesh.h:1584
bool areElemIDsIdentical(const std::string &id_name1, const std::string &id_name2) const
Whether or not two extra element integers are identical.
Definition: MooseMesh.h:2208
bool _node_to_active_semilocal_elem_map_built
Definition: MooseMesh.h:1541
std::map< boundary_id_type, std::set< dof_id_type > > _bnd_node_ids
Map of sets of node IDs in each boundary.
Definition: MooseMesh.h:1568
ConstElemPointerRange * refinedElementRange() const
Return a range that is suitable for threaded execution over elements that were just refined...
Definition: MooseMesh.C:946
virtual void init()
Initialize the Mesh object.
Definition: MooseMesh.C:2913
The definition of the const_face_info_iterator struct.
Definition: MooseMesh.h:1967
bool _skip_refine_when_use_split
Whether or not to skip uniform refinements when using a pre-split mesh.
Definition: MooseMesh.h:1491
unsigned int getHigherDSide(const Elem *elem) const
Returns the local side ID of the interior parent aligned with the lower dimensional element...
Definition: MooseMesh.C:1750
std::unordered_map< std::pair< const Elem *, unsigned int >, FaceInfo * > _elem_side_to_face_info
Map from elem-side pair to FaceInfo.
Definition: MooseMesh.h:1643
void detectPairedSidesets()
This routine detects paired sidesets of a regular orthogonal mesh (.i.e.
Definition: MooseMesh.C:2037
const std::set< SubdomainID > & getNodeBlockIds(const Node &node) const
Return list of blocks to which the given node belongs.
Definition: MooseMesh.C:1537
unsigned int getMaxNodesPerElem() const
Get the maximum number of nodes per element.
Definition: MooseMesh.h:1162
const Parallel::Communicator & _communicator
std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap > > _elem_type_to_p_coarsening_side_map
Definition: MooseMesh.h:1815
void setMeshBoundaryIDs(std::set< BoundaryID > boundary_IDs)
Sets the set of BoundaryIDs Is called by AddAllSideSetsByNormals.
Definition: MooseMesh.C:3235
std::map< boundary_id_type, std::vector< dof_id_type > > _node_set_nodes
list of nodes that belongs to a specified nodeset: indexing [nodeset_id] -> [array of node ids] ...
Definition: MooseMesh.h:1587
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
std::set< SubdomainID > _lower_d_boundary_blocks
Mesh blocks for boundary lower-d elements in different types.
Definition: MooseMesh.h:1839
bool _linear_finite_volume_dofs_cached
Definition: MooseMesh.h:1651
void markFiniteVolumeInfoDirty()
Mark the finite volume information as dirty.
Definition: MooseMesh.h:1325
void cacheInfo()
Definition: MooseMesh.C:1445
std::basic_ostream< charT, traits > * os
Definition: InfixIterator.h:33
void changeBoundaryId(const boundary_id_type old_id, const boundary_id_type new_id, bool delete_prev)
Change all the boundary IDs for a given side from old_id to new_id.
Definition: MooseMesh.C:2800
std::vector< std::unique_ptr< libMesh::GhostingFunctor > > _ghosting_functors
Deprecated (DO NOT USE)
Definition: MooseMesh.h:1439
std::set< Elem * > _ghost_elems_from_ghost_boundaries
Set of elements ghosted by ghostGhostedBoundaries.
Definition: MooseMesh.h:1862
virtual void buildMesh()=0
Must be overridden by child classes.
void setPartitionerHelper(MeshBase *mesh=nullptr)
Definition: MooseMesh.C:3668
void deleteRemoteElements()
Delete remote elements.
Definition: MooseMesh.C:3969
unsigned int _to
The qp to map to.
Definition: MooseMesh.h:82
libMesh::ConstNodeRange * getLocalNodeRange()
Definition: MooseMesh.C:1313
virtual const Node & nodeRef(const dof_id_type i) const
Definition: MooseMesh.C:849
bool _allow_recovery
Whether or not this Mesh is allowed to read a recovery file.
Definition: MooseMesh.h:1850
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 buildNodeListFromSideList()
Calls BoundaryInfo::build_node_list_from_side_list().
Definition: MooseMesh.C:3023
const std::string & getSubdomainName(SubdomainID subdomain_id) const
Return the name of a block given an id.
Definition: MooseMesh.C:1801
void setPatchUpdateStrategy(Moose::PatchUpdateType patch_update_strategy)
Set the patch size update strategy.
Definition: MooseMesh.C:3447
void buildFiniteVolumeInfo() const
Builds the face and elem info vectors that store meta-data needed for looping over and doing calculat...
Definition: MooseMesh.C:3779
std::unordered_map< SubdomainID, std::set< BoundaryID > > _neighbor_subdomain_boundary_ids
Holds a map from neighbor subomdain ids to the boundary ids that are attached to it.
Definition: MooseMesh.h:1834
const std::pair< Point, RealVectorValue > & getGeneralAxisymmetricCoordAxis(SubdomainID subdomain_id) const
Gets the general axisymmetric coordinate axis for a block.
Definition: MooseMesh.C:4309
auto max(const L &left, const R &right)
const std::set< BoundaryID > & meshBoundaryIds() const
Returns a read-only reference to the set of boundary IDs currently present in the Mesh...
Definition: MooseMesh.C:3217
std::set< SubdomainID > _lower_d_interior_blocks
Mesh blocks for interior lower-d elements in different types.
Definition: MooseMesh.h:1837
virtual Elem * queryElemPtr(const dof_id_type i)
Definition: MooseMesh.C:3165
bool isLowerD(const SubdomainID subdomain_id) const
Definition: MooseMesh.h:2240
elem_info_iterator ownedElemInfoBegin()
Iterators to owned faceInfo objects.
Definition: MooseMesh.C:1566
void setIsCustomPartitionerRequested(bool cpr)
Definition: MooseMesh.C:3767
std::set< BoundaryID > boundary_ids
The boundary ids that are attached.
Definition: MooseMesh.h:1824
unsigned int _max_h_level
Maximum h-refinement level of all elements.
Definition: MooseMesh.h:1924
virtual bnd_node_iterator bndNodesBegin()
Return iterators to the beginning/end of the boundary nodes list.
Definition: MooseMesh.C:1583
void mapPoints(const std::vector< Point > &from, const std::vector< Point > &to, std::vector< QpMap > &qp_map)
Find the closest points that map "from" to "to" and fill up "qp_map".
Definition: MooseMesh.C:2628
bnd_elem_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2087
bool is_lower_d
Whether this subdomain is a lower-dimensional manifold of a higher-dimensional subdomain.
Definition: MooseMesh.h:1827
void attachRelationshipManagers(Moose::RelationshipManagerType rm_type, bool attach_geometric_rm_final=false)
Attach the relationship managers of the given type Note: Geometric relationship managers that are sup...
Definition: MooseApp.C:3275
const Moose::Kokkos::Mesh * getKokkosMesh() const
Accessor for Kokkos mesh object.
Definition: MooseMesh.h:674
void errorIfDistributedMesh(std::string name) const
Generate a unified error message if the underlying libMesh mesh is a DistributedMesh.
Definition: MooseMesh.C:3657
const MeshBase * getMeshPtr() const
Definition: MooseMesh.C:3482
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
void updateCoordTransform()
Update the coordinate transformation object based on our coordinate system data.
Definition: MooseMesh.C:4325
const std::vector< const FaceInfo * > & faceInfo() const
Accessor for local FaceInfo objects.
Definition: MooseMesh.h:2216
dof_id_type minElementID(unsigned int elem_id_index) const
Return the minimum element ID for an extra element integer with its accessing index.
Definition: MooseMesh.h:1135
const std::vector< FaceInfo > & allFaceInfo() const
Accessor for all FaceInfo objects.
Definition: MooseMesh.h:2222
bool _use_distributed_mesh
False by default.
Definition: MooseMesh.h:1454
bool isSplit() const
Definition: MooseMesh.h:1353
std::unique_ptr< std::map< BoundaryID, RealVectorValue > > _boundary_to_normal_map
The boundary to normal map - valid only when AddAllSideSetsByNormals is active.
Definition: MooseMesh.h:1561
std::vector< BndElement * >::const_iterator const_bnd_elem_iterator_imp
Definition: MooseMesh.h:1573
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:103
bool hasElementID(const std::string &id_name) const
Whether mesh has an extra element integer with a given name.
Definition: MooseMesh.h:2194
bool _built_from_other_mesh
Whether or not this mesh was built from another mesh.
Definition: MooseMesh.h:1445
virtual dof_id_type nLocalNodes() const
Definition: MooseMesh.h:328
bool _allow_remote_element_removal
Whether to allow removal of remote elements.
Definition: MooseMesh.h:1859
std::unique_ptr< Moose::Kokkos::Mesh > _kokkos_mesh
Pointer to Kokkos mesh object.
Definition: MooseMesh.h:1463
bnd_node_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2044
virtual bool skipNoncriticalPartitioning() const
Definition: MooseMesh.C:4442
SemiLocalNodeRange * getActiveSemiLocalNodeRange() const
Definition: MooseMesh.C:1304
int8_t boundary_id_type
std::unordered_set< dof_id_type > getBoundaryActiveNeighborElemIds(BoundaryID bid) const
Return all ids of neighbors of elements which have a side which is part of a sideset.
Definition: MooseMesh.C:1380
void setSubdomainName(SubdomainID subdomain_id, const SubdomainName &name)
This method sets the name for subdomain_id to name.
Definition: MooseMesh.C:1787
void clearQuadratureNodes()
Clear out any existing quadrature nodes.
Definition: MooseMesh.C:1718
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3488
const std::map< SubdomainID, Moose::CoordinateSystemType > & getCoordSystem() const
Get the map from subdomain ID to coordinate system type, e.g.
Definition: MooseMesh.C:4247
std::vector< BndNode * > _bnd_nodes
array of boundary nodes
Definition: MooseMesh.h:1564
std::unique_ptr< T > buildTypedMesh(unsigned int dim=libMesh::invalid_uint)
Shortcut method to construct a unique pointer to a libMesh mesh instance.
Definition: MooseMesh.h:2133
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:27
std::vector< std::shared_ptr< RelationshipManager > > _relationship_managers
The list of active geometric relationship managers (bound to the underlying MeshBase object)...
Definition: MooseMesh.h:1442
virtual unsigned int dimension() const
Returns MeshBase::mesh_dimension(), (not MeshBase::spatial_dimension()!) of the underlying libMesh me...
Definition: MooseMesh.C:2968
unsigned int _from
The qp to map from.
Definition: MooseMesh.h:79
std::pair< const Node *, BoundaryID > PeriodicNodeInfo
Helper type for building periodic node maps.
Definition: MooseMesh.h:1085
virtual SubdomainID nSubdomains() const
Definition: MooseMesh.h:331
const_elem_info_iterator(const IterType &d, const IterType &e, const PredType &p)
Definition: MooseMesh.h:2018
std::unique_ptr< MeshBase > buildMeshBaseObject(unsigned int dim=libMesh::invalid_uint)
Method to construct a libMesh::MeshBase object that is normally set and used by the MooseMesh object ...
Definition: MooseMesh.C:2894
unsigned int getAxisymmetricRadialCoord() const
Returns the desired radial direction for RZ coordinate transformation.
Definition: MooseMesh.C:4334
bool isPartitionerForced() const
Tell the user if the partitioner was overriden for any reason.
Definition: MooseMesh.h:1037
const std::set< unsigned int > & getGhostedBoundaries() const
Return a writable reference to the set of ghosted boundary IDs.
Definition: MooseMesh.C:3281
bool _construct_node_list_from_side_list
Whether or not to allow generation of nodesets from sidesets.
Definition: MooseMesh.h:1853
boundary_id_type BoundaryID
unsigned int sideWithBoundaryID(const Elem *const elem, const BoundaryID boundary_id) const
Calls BoundaryInfo::side_with_boundary_id().
Definition: MooseMesh.C:3061
const std::vector< dof_id_type > & getNodeList(boundary_id_type nodeset_id) const
Return a writable reference to a vector of node IDs that belong to nodeset_id.
Definition: MooseMesh.C:3510
unsigned int getMaxLeafSize() const
Getter for the maximum leaf size parameter.
Definition: MooseMesh.h:635
virtual const Node * nodePtr(const dof_id_type i) const
Definition: MooseMesh.C:863
bool hasLowerD() const
Definition: MooseMesh.h:1424
std::vector< dof_id_type > _max_ids
Maximum integer ID for each extra element integer.
Definition: MooseMesh.h:1873
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:92
void update()
Calls buildNodeListFromSideList(), buildNodeList(), and buildBndElemList().
Definition: MooseMesh.C:629
const std::pair< BoundaryID, BoundaryID > * getPairedBoundaryMapping(unsigned int component)
This function attempts to return the paired boundary ids for the given component. ...
Definition: MooseMesh.C:2331
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
std::set< BoundaryID > _mesh_nodeset_ids
Definition: MooseMesh.h:1557
std::unique_ptr< MooseAppCoordTransform > _coord_transform
A coordinate transformation object that describes how to transform this problem&#39;s coordinate system i...
Definition: MooseMesh.h:1911
std::set< SubdomainID > getBoundaryConnectedBlocks(const BoundaryID bid) const
Get the list of subdomains associated with the given boundary.
Definition: MooseMesh.C:3561
void checkCoordinateSystems()
Performs a sanity check for every element in the mesh.
Definition: MooseMesh.C:4347
ParallelType getParallelType() const
Definition: MooseMesh.h:1027
std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap > > _elem_type_to_p_refinement_side_map
Definition: MooseMesh.h:1789
std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap > > _elem_type_to_p_coarsening_map
Definition: MooseMesh.h:1813
const bool _is_split
Whether or not we are using a (pre-)split mesh (automatically DistributedMesh)
Definition: MooseMesh.h:1617
void setCoordSystem(const std::vector< SubdomainName > &blocks, const MultiMooseEnum &coord_sys)
Set the coordinate system for the provided blocks to coord_sys.
Definition: MooseMesh.C:4121
std::unique_ptr< ConstElemPointerRange > _coarsened_elements
The elements that were just coarsened.
Definition: MooseMesh.h:1509
std::set< dof_id_type > getAllElemIDs(unsigned int elem_id_index) const
Return all the unique element IDs for an extra element integer with its index.
Definition: MooseMesh.C:1167
std::set< BoundaryID > _mesh_boundary_ids
A set of boundary IDs currently present in the mesh.
Definition: MooseMesh.h:1555
const Node * addUniqueNode(const Point &p, Real tol=1e-6)
Add a new node to the mesh.
Definition: MooseMesh.C:1614
std::vector< BndNode > _extra_bnd_nodes
Definition: MooseMesh.h:1581
MooseApp & _app
The MOOSE application this is associated with.
Definition: MooseBase.h:357
bool _moose_mesh_prepared
True if prepare has been called on the mesh.
Definition: MooseMesh.h:1503
void buildRefinementMap(const Elem &elem, libMesh::QBase &qrule, libMesh::QBase &qrule_face, int parent_side, int child, int child_side)
Build the refinement map for a given element type.
Definition: MooseMesh.C:2510
The definition of the face_info_iterator struct.
Definition: MooseMesh.h:1952
unsigned int uniformRefineLevel() const
Returns the level of uniform refinement requested (zero if AMR is disabled).
Definition: MooseMesh.C:3256
std::vector< BndElement * > _bnd_elems
array of boundary elems
Definition: MooseMesh.h:1571
bool _coord_system_set
Whether the coordinate system has been set.
Definition: MooseMesh.h:1914
const std::set< SubdomainID > & boundaryLowerDBlocks() const
Definition: MooseMesh.h:1433
std::vector< SubdomainName > getSubdomainNames(const std::vector< SubdomainID > &subdomain_ids) const
Get the associated subdomainNames for the subdomain ids that are passed in.
Definition: MooseMesh.C:1807
MooseMesh & operator=(const MooseMesh &other_mesh)=delete
const std::vector< QpMap > & getPRefinementMap(const Elem &elem) const
Get the map describing for each volumetric quadrature point (qp) on the refined level which qp on the...
Definition: MooseMesh.C:4418
virtual dof_id_type nActiveLocalElem() const
Definition: MooseMesh.h:330
Real minPeriodicDistance(unsigned int nonlinear_var_num, Point p, Point q) const
This function returns the distance between two points on the mesh taking into account periodicity for...
Definition: MooseMesh.C:2325
bool hasSecondOrderElements()
check if the mesh has SECOND order elements
Definition: MooseMesh.C:3751
unsigned int getPatchSize() const
Getter for the patch_size parameter.
Definition: MooseMesh.C:3441
void buildRefinementAndCoarseningMaps(Assembly *assembly)
Create the refinement and coarsening maps necessary for projection of stateful material properties wh...
Definition: MooseMesh.C:2500
bool _parallel_type_overridden
Definition: MooseMesh.h:1456
bool isRegularOrthogonal()
Getter to query if the mesh was detected to be regular and orthogonal.
Definition: MooseMesh.h:1067
unsigned int getMaxSidesPerElem() const
Get the maximum number of sides per element.
Definition: MooseMesh.h:1157
Interface for objects interacting with the PerfGraph.
MeshBase::element_iterator activeLocalElementsBegin()
Calls active_local_nodes_begin/end() on the underlying libMesh mesh object.
Definition: MooseMesh.C:3091
std::vector< Node * > _node_map
Vector of all the Nodes in the mesh for determining when to add a new point.
Definition: MooseMesh.h:1605
std::map< dof_id_type, std::map< unsigned int, std::map< dof_id_type, Node * > > > _elem_to_side_to_qp_to_quadrature_nodes
Definition: MooseMesh.h:1580
const std::vector< QpMap > & getPRefinementSideMap(const Elem &elem) const
Get the map describing for each side quadrature point (qp) on the refined level which qp on the previ...
Definition: MooseMesh.C:4424
std::unordered_map< dof_id_type, ElemInfo > _elem_to_elem_info
Map connecting elems with their corresponding ElemInfo, we use the element ID as the key...
Definition: MooseMesh.h:1627
bool skipRefineWhenUseSplit() const
Whether or not skip uniform refinements when using a pre-split mesh.
Definition: MooseMesh.h:589
std::vector< BndElement * >::iterator bnd_elem_iterator_imp
Definition: MooseMesh.h:1572
Node * getQuadratureNode(const Elem *elem, const unsigned short int side, const unsigned int qp)
Get a specified quadrature node.
Definition: MooseMesh.C:1700
void printInfo(std::ostream &os=libMesh::out, const unsigned int verbosity=0) const
Calls print_info() on the underlying Mesh.
Definition: MooseMesh.C:3502
libMesh::StoredRange< MooseMesh::const_bnd_node_iterator, const BndNode * > ConstBndNodeRange
Some useful StoredRange typedefs.
Definition: MooseMesh.h:2127
virtual dof_id_type nNodes() const
Calls n_nodes/elem() on the underlying libMesh mesh object.
Definition: MooseMesh.C:3115
static MooseEnum partitioning()
returns MooseMesh partitioning options so other classes can use it
Definition: MooseMesh.C:3938
MooseAppCoordTransform & coordTransform()
Definition: MooseMesh.h:1931
void buildCoarseningMap(const Elem &elem, libMesh::QBase &qrule, libMesh::QBase &qrule_face, int input_side)
Build the coarsening map for a given element type.
Definition: MooseMesh.C:2589
const std::vector< const Elem * > & coarsenedElementChildren(const Elem *elem) const
Get the newly removed children element ids for an element that was just coarsened.
Definition: MooseMesh.C:958
std::unordered_map< dof_id_type, std::set< dof_id_type > > getElemIDMapping(const std::string &from_id_name, const std::string &to_id_name) const
Definition: MooseMesh.C:1138
std::map< std::pair< int, libMesh::ElemType >, std::vector< std::pair< unsigned int, QpMap > > > _elem_type_to_coarsening_map
Holds mappings for volume to volume and parent side to child side Map key:
Definition: MooseMesh.h:1810
void setBoundaryName(BoundaryID boundary_id, BoundaryName name)
This method sets the boundary name of the boundary based on the id parameter.
Definition: MooseMesh.C:1818
Moose::PatchUpdateType _patch_update_strategy
The patch update strategy.
Definition: MooseMesh.h:1602
Physical unit management class with runtime unit string parsing, unit checking, unit conversion...
Definition: Units.h:32
void addPeriodicVariable(unsigned int var_num, BoundaryID primary, BoundaryID secondary)
For "regular orthogonal" meshes, determine if variable var_num is periodic with respect to the primar...
Definition: MooseMesh.C:2255
std::set< BoundaryID > _mesh_sideset_ids
Definition: MooseMesh.h:1556
void setGeneralAxisymmetricCoordAxes(const std::vector< SubdomainName > &blocks, const std::vector< std::pair< Point, RealVectorValue >> &axes)
Sets the general coordinate axes for axisymmetric blocks.
Definition: MooseMesh.C:4261
void setBoundaryToNormalMap(std::unique_ptr< std::map< BoundaryID, RealVectorValue >> boundary_map)
Sets the mapping between BoundaryID and normal vector Is called by AddAllSideSetsByNormals.
Definition: MooseMesh.C:3241
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _partitioner_overridden
Definition: MooseMesh.h:1468
bool detectOrthogonalDimRanges(Real tol=1e-6)
This routine determines whether the Mesh is a regular orthogonal mesh (i.e.
Definition: MooseMesh.C:1963
RelationshipManagers are used for describing what kinds of non-local resources are needed for an obje...
void updateActiveSemiLocalNodeRange(std::set< dof_id_type > &ghosted_elems)
Clears the "semi-local" node list and rebuilds it.
Definition: MooseMesh.C:966
void isDisplaced(bool is_displaced)
Set whether this mesh is a displaced mesh.
Definition: MooseMesh.h:1234
CoordinateSystemType
Definition: MooseTypes.h:810
std::set< dof_id_type > getElemIDsOnBlocks(unsigned int elem_id_index, const std::set< SubdomainID > &blks) const
Return all the unique element IDs for an extra element integer with its index on a set of subdomains...
Definition: MooseMesh.C:1177
unsigned int getBlocksMaxDimension(const std::vector< SubdomainName > &blocks) const
Returns the maximum element dimension on the given blocks.
Definition: MooseMesh.C:2989
virtual std::unique_ptr< MooseMesh > safeClone() const =0
A safer version of the clone() method that hands back an allocated object wrapped in a smart pointer...
bool _is_displaced
Whether this mesh is displaced.
Definition: MooseMesh.h:1892
bool hasMeshBase() const
Whether mesh base object was constructed or not.
Definition: MooseMesh.h:1115
std::unique_ptr< libMesh::MeshBase > _mesh
Pointer to underlying libMesh mesh object.
Definition: MooseMesh.h:1459
libMesh::NodeRange * getActiveNodeRange()
Definition: MooseMesh.C:1290
virtual unsigned int nPartitions() const
Definition: MooseMesh.h:332
void setGhostedBoundaryInflation(const std::vector< Real > &inflation)
This sets the inflation amount for the bounding box for each partition for use in ghosting boundaries...
Definition: MooseMesh.C:3275
void freeBndNodes()
Definition: MooseMesh.C:386
The definition of the elem_info_iterator struct.
Definition: MooseMesh.h:1996
Real dimensionWidth(unsigned int component) const
Returns the width of the requested dimension.
Definition: MooseMesh.C:2231
PatchUpdateType
Type of patch update strategy for modeling node-face constraints or contact.
Definition: MooseTypes.h:953
bool _skip_deletion_repartition_after_refine
Whether or not skip remote deletion and repartition after uniform refinements.
Definition: MooseMesh.h:1494
std::set< SubdomainID > getBoundaryConnectedSecondaryBlocks(const BoundaryID bid) const
Get the list of subdomains associated with the given boundary of its secondary side.
Definition: MooseMesh.C:3572
unsigned int _max_nodes_per_elem
The maximum number of nodes per element.
Definition: MooseMesh.h:1883
bool _need_delete
Whether we need to delete remote elements after init&#39;ing the EquationSystems.
Definition: MooseMesh.h:1856
std::map< std::pair< libMesh::ElemType, unsigned int >, std::vector< QpMap > > _elem_type_to_p_refinement_map
Definition: MooseMesh.h:1787
void setAxisymmetricCoordAxis(const MooseEnum &rz_coord_axis)
For axisymmetric simulations, set the symmetry coordinate axis.
Definition: MooseMesh.C:4253
std::set< BoundaryID > getSubdomainInterfaceBoundaryIds(const SubdomainID subdomain_id) const
Get the list of boundaries that contact the given subdomain.
Definition: MooseMesh.C:3548
virtual const Node & node(const dof_id_type i) const
Various accessors (pointers/references) for Node "i".
Definition: MooseMesh.C:835
unsigned int maxPLevel() const
Returns the maximum p-refinement level of all elements.
Definition: MooseMesh.h:1383
libMesh::BoundingBox getInflatedProcessorBoundingBox(Real inflation_multiplier=0.01) const
Get a (slightly inflated) processor bounding box.
Definition: MooseMesh.C:3459
std::map< unsigned int, std::vector< bool > > _periodic_dim
A map of vectors indicating which dimensions are periodic in a regular orthogonal mesh for the specif...
Definition: MooseMesh.h:1657
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
void setMeshBase(std::unique_ptr< MeshBase > mesh_base)
Method to set the mesh_base object.
Definition: MooseMesh.C:2906
void buildBndElemList()
Definition: MooseMesh.C:1193
std::set< SubdomainID > getInterfaceConnectedBlocks(const BoundaryID bid) const
Get the list of subdomains contacting the given boundary.
Definition: MooseMesh.C:3583
bool isParallelTypeForced() const
Tell the user if the distribution was overriden for any reason.
Definition: MooseMesh.h:1017
std::vector< Real > _ghosted_boundaries_inflation
Definition: MooseMesh.h:1590
std::vector< std::unordered_map< SubdomainID, std::set< dof_id_type > > > _block_id_mapping
Unique element integer IDs for each subdomain and each extra element integers.
Definition: MooseMesh.h:1871
void needGhostGhostedBoundaries(bool needghost)
Whether or not we want to ghost ghosted boundaries.
Definition: MooseMesh.h:620
void findAdaptivityQpMaps(const Elem *template_elem, libMesh::QBase &qrule, libMesh::QBase &qrule_face, std::vector< std::vector< QpMap >> &refinement_map, std::vector< std::pair< unsigned int, QpMap >> &coarsen_map, int parent_side, int child, int child_side)
Given an elem type, get maps that tell us what qp&#39;s are closest to each other between a parent and it...
Definition: MooseMesh.C:2659
bool needsRemoteElemDeletion() const
Whether we need to delete remote elements.
Definition: MooseMesh.h:1095
std::vector< BndNode * >::const_iterator const_bnd_node_iterator_imp
Definition: MooseMesh.h:1566
unsigned int _patch_size
The number of nodes to consider in the NearestNode neighborhood.
Definition: MooseMesh.h:1593
bool skipDeletionRepartitionAfterRefine() const
Return a flag indicating whether or not we should skip remote deletion and repartition after uniform ...
Definition: MooseMesh.h:2181
elem_info_iterator ownedElemInfoEnd()
Definition: MooseMesh.C:1574
void buildPeriodicNodeSets(std::map< BoundaryID, std::set< dof_id_type >> &periodic_node_sets, unsigned int var_number, libMesh::PeriodicBoundaries *pbs) const
This routine builds a datastructure of node ids organized by periodic boundary ids.
Definition: MooseMesh.C:1936
virtual std::unique_ptr< libMesh::PointLocatorBase > getPointLocator() const
Proxy function to get a (sub)PointLocator from either the underlying libMesh mesh (default)...
Definition: MooseMesh.C:3773
void addGhostedBoundary(BoundaryID boundary_id)
This will add the boundary ids to be ghosted to this processor.
Definition: MooseMesh.C:3269
virtual dof_id_type maxNodeId() const
Calls max_node/elem_id() on the underlying libMesh mesh object.
Definition: MooseMesh.C:3127
virtual Elem * elem(const dof_id_type i)
Various accessors (pointers/references) for Elem "i".
Definition: MooseMesh.C:3139
libMesh::StoredRange< MooseMesh::const_bnd_elem_iterator, const BndElement * > * getBoundaryElementRange()
Definition: MooseMesh.C:1341
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:199
virtual dof_id_type nActiveElem() const
Definition: MooseMesh.h:329
const_elem_info_iterator(const MooseMesh::elem_info_iterator &rhs)
Definition: MooseMesh.h:2028
face_info_iterator ownedFaceInfoEnd()
Definition: MooseMesh.C:1557
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
std::unordered_map< SubdomainID, std::pair< Point, RealVectorValue > > _subdomain_id_to_rz_coord_axis
Map of subdomain ID to general axisymmetric axis.
Definition: MooseMesh.h:1907
const MooseEnum & partitionerName() const
Definition: MooseMesh.h:1032
Class used for caching additional information for elements such as the volume and centroid...
Definition: ElemInfo.h:25
MeshBase::node_iterator localNodesEnd()
Definition: MooseMesh.C:3073
const RealVectorValue & getNormalByBoundaryID(BoundaryID id) const
Returns the normal vector associated with a given BoundaryID.
Definition: MooseMesh.C:2850
void doingPRefinement(bool doing_p_refinement)
Indicate whether the kind of adaptivity we&#39;re doing is p-refinement.
Definition: MooseMesh.h:1373
std::set< SubdomainID > _mesh_subdomains
A set of subdomain IDs currently present in the mesh.
Definition: MooseMesh.h:1547
ConstElemPointerRange * coarsenedElementRange() const
Return a range that is suitable for threaded execution over elements that were just coarsened...
Definition: MooseMesh.C:952
const Moose::PatchUpdateType & getPatchUpdateStrategy() const
Get the current patch update strategy.
Definition: MooseMesh.C:3453
bool doingPRefinement() const
Query whether we have p-refinement.
Definition: MooseMesh.h:1378
std::vector< std::pair< BoundaryID, BoundaryID > > _paired_boundary
A vector holding the paired boundaries for a regular orthogonal mesh.
Definition: MooseMesh.h:1614
virtual bool skipPartitioning() const
Definition: MooseMesh.h:333
unsigned int maxHLevel() const
Returns the maximum h-refinement level of all elements.
Definition: MooseMesh.h:1388
void ghostGhostedBoundaries()
Actually do the ghosting of boundaries that need to be ghosted to this processor. ...
Definition: MooseMesh.C:3354
std::set< unsigned int > _ghosted_boundaries
Definition: MooseMesh.h:1589
MeshBase::node_iterator localNodesBegin()
Calls local_nodes_begin/end() on the underlying libMesh mesh object.
Definition: MooseMesh.C:3067
unsigned int _rz_coord_axis
Storage for RZ axis selection.
Definition: MooseMesh.h:1904
void computeFiniteVolumeCoords() const
Compute the face coordinate value for all FaceInfo and ElemInfo objects.
Definition: MooseMesh.C:3917
bool _distribution_overridden
Definition: MooseMesh.h:1455
void buildNodeList()
Calls BoundaryInfo::build_node_list()/build_side_list() and makes separate copies of Nodes/Elems in t...
Definition: MooseMesh.C:1040
std::unordered_map< SubdomainID, SubdomainData > _sub_to_data
Holds a map from subdomain ids to associated data.
Definition: MooseMesh.h:1831
std::unique_ptr< libMesh::Partitioner > _custom_partitioner
The custom partitioner.
Definition: MooseMesh.h:1471
bool isBoundaryElem(dof_id_type elem_id) const
Returns true if the requested element is in the list of boundary elements, false otherwise.
Definition: MooseMesh.C:3631
unsigned int getElementIDIndex(const std::string &id_name) const
Return the accessing integer for an extra element integer with its name.
Definition: MooseMesh.h:2200
void needsRemoteElemDeletion(bool need_delete)
Set whether we need to delete remote elements.
Definition: MooseMesh.h:1090
std::map< const Elem *, std::vector< const Elem * > > _coarsened_element_children
Map of Parent elements to child elements for elements that were just coarsened.
Definition: MooseMesh.h:1516
virtual bool isDistributedMesh() const
Returns the final Mesh distribution type.
Definition: MooseMesh.h:1012
const std::vector< QpMap > & getPCoarseningMap(const Elem &elem) const
Get the map describing for each volumetric quadrature point (qp) on the coarse level which qp on the ...
Definition: MooseMesh.C:4430
OStreamProxy out(std::cout)
void buildPeriodicNodeMap(std::multimap< dof_id_type, dof_id_type > &periodic_node_map, unsigned int var_number, libMesh::PeriodicBoundaries *pbs) const
This routine builds a multimap of boundary ids to matching boundary ids across all periodic boundarie...
Definition: MooseMesh.C:1851
unsigned int _ghosting_patch_size
The number of nearest neighbors to consider for ghosting purposes when iteration patch update strateg...
Definition: MooseMesh.h:1596
const std::unordered_map< boundary_id_type, std::unordered_set< dof_id_type > > & getBoundariesToActiveSemiLocalElemIds() const
Returns a map of boundaries to ids of elements on the boundary.
Definition: MooseMesh.C:1363
const MeshBase::element_iterator activeLocalElementsEnd()
Definition: MooseMesh.C:3097
std::unique_ptr< libMesh::ConstNodeRange > _local_node_range
Definition: MooseMesh.h:1529
QpMap()
Definition: MooseMesh.h:76
libMesh::StoredRange< MooseMesh::const_bnd_node_iterator, const BndNode * > * getBoundaryNodeRange()
Definition: MooseMesh.C:1327
The definition of the const_elem_info_iterator struct.
Definition: MooseMesh.h:2011
bool _regular_orthogonal_mesh
Boolean indicating whether this mesh was detected to be regular and orthogonal.
Definition: MooseMesh.h:1608
std::map< dof_id_type, std::vector< dof_id_type > > _node_to_active_semilocal_elem_map
A map of all of the current nodes to the active elements that they are connected to.
Definition: MooseMesh.h:1540
bool prepare(const MeshBase *mesh_to_clone)
Calls prepare_for_use() if the underlying MeshBase object isn&#39;t prepared, then communicates various b...
Definition: MooseMesh.C:419
const ElemInfo & elemInfo(const dof_id_type id) const
Accessor for the elemInfo object for a given element ID.
Definition: MooseMesh.C:3911
void setCustomPartitioner(libMesh::Partitioner *partitioner)
Setter for custom partitioner.
Definition: MooseMesh.C:3739
Real _distance
The distance between them.
Definition: MooseMesh.h:85
unsigned int nFace() const
accessors for the FaceInfo objects
Definition: MooseMesh.h:1173
static MooseEnum elemTypes()
returns MooseMesh element type options
Definition: MooseMesh.C:3946
Base variable class.
void meshChanged()
Declares that the MooseMesh has changed, invalidates cached data and rebuilds caches.
Definition: MooseMesh.C:897
void buildPRefinementAndCoarseningMaps(Assembly *assembly)
Definition: MooseMesh.C:2406
BoundaryID getBoundaryID(const BoundaryName &boundary_name) const
Get the associated BoundaryID for the boundary name.
Definition: MooseMesh.C:1730
std::unique_ptr< libMesh::StoredRange< MooseMesh::const_bnd_node_iterator, const BndNode * > > _bnd_node_range
Definition: MooseMesh.h:1531
uint8_t dof_id_type
virtual dof_id_type nElem() const
Definition: MooseMesh.C:3121
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
bool isBoundaryFullyExternalToSubdomains(BoundaryID bid, const std::set< SubdomainID > &blk_group) const
Returns whether a boundary (given by its id) is not crossing through a group of blocks, by which we mean that elements on both sides of the boundary are in those blocks.
Definition: MooseMesh.C:1411
std::unique_ptr< SemiLocalNodeRange > _active_semilocal_node_range
Definition: MooseMesh.h:1527
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
unsigned int _max_p_level
Maximum p-refinement level of all elements.
Definition: MooseMesh.h:1922
ParallelType
void setupFiniteVolumeMeshData() const
Sets up the additional data needed for finite volume computations.
Definition: MooseMesh.C:4112
virtual unsigned int effectiveSpatialDimension() const
Returns the effective spatial dimension determined by the coordinates actually used by the mesh...
Definition: MooseMesh.C:2974