Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #pragma once
11 :
12 : #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;
43 : class RelationshipManager;
44 : class MooseVariableBase;
45 : class MooseAppCoordTransform;
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
59 : typedef libMesh::StoredRange<std::set<Node *>::iterator, Node *> SemiLocalNodeRange;
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 :
70 : /**
71 : * Helper object for holding qp mapping info.
72 : */
73 : class QpMap
74 : {
75 : public:
76 65270 : QpMap() : _distance(std::numeric_limits<Real>::max()) {}
77 :
78 : /// The qp to map from
79 : unsigned int _from;
80 :
81 : /// The qp to map to
82 : unsigned int _to;
83 :
84 : /// The distance between them
85 : Real _distance;
86 : };
87 :
88 : /**
89 : * MooseMesh wraps a libMesh::Mesh object and enhances its capabilities
90 : * by caching additional data and storing more state.
91 : */
92 : class MooseMesh : public MooseObject, public Restartable, public PerfGraphInterface
93 : {
94 : public:
95 : /**
96 : * Typical "Moose-style" constructor and copy constructor.
97 : */
98 : static InputParameters validParams();
99 :
100 : MooseMesh(const InputParameters & parameters);
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,
112 : DISTRIBUTED
113 : };
114 :
115 : /**
116 : * Clone method. Allocates memory you are responsible to clean up.
117 : */
118 : virtual MooseMesh & clone() const;
119 :
120 : /**
121 : * A safer version of the clone() method that hands back an
122 : * allocated object wrapped in a smart pointer. This makes it much
123 : * less likely that the caller will leak the memory in question.
124 : */
125 : virtual std::unique_ptr<MooseMesh> safeClone() const = 0;
126 :
127 : /**
128 : * Determine whether to use a distributed mesh. Should be called during construction
129 : */
130 : void determineUseDistributedMesh();
131 :
132 : /**
133 : * Method to construct a libMesh::MeshBase object that is normally set and used by the MooseMesh
134 : * object during the "init()" phase. If the parameter \p dim is not
135 : * provided, then its value will be taken from the input file mesh block.
136 : */
137 : std::unique_ptr<MeshBase> buildMeshBaseObject(unsigned int dim = libMesh::invalid_uint);
138 :
139 : /**
140 : * Shortcut method to construct a unique pointer to a libMesh mesh instance. The created
141 : * derived-from-MeshBase object will have its \p allow_remote_element_removal flag set to whatever
142 : * our value is. We will also attach any geometric \p RelationshipManagers that have been
143 : * requested by our simulation objects to the \p MeshBase object. If the parameter \p dim is not
144 : * provided, then its value will be taken from the input file mesh block.
145 : */
146 : template <typename T>
147 : std::unique_ptr<T> buildTypedMesh(unsigned int dim = libMesh::invalid_uint);
148 :
149 : /**
150 : * Method to set the mesh_base object. If this method is NOT called prior to calling init(), a
151 : * MeshBase object will be automatically constructed and set.
152 : */
153 : void setMeshBase(std::unique_ptr<MeshBase> mesh_base);
154 :
155 : /// returns MooseMesh partitioning options so other classes can use it
156 : static MooseEnum partitioning();
157 :
158 : /// returns MooseMesh element type options
159 : static MooseEnum elemTypes();
160 :
161 : /**
162 : * Initialize the Mesh object. Most of the time this will turn around
163 : * and call build_mesh so the child class can build the Mesh object.
164 : *
165 : * However, during Recovery this will read the CPA file...
166 : */
167 : virtual void init();
168 :
169 : /**
170 : * Must be overridden by child classes.
171 : *
172 : * This is where the Mesh object is actually created and filled in.
173 : */
174 : virtual void buildMesh() = 0;
175 :
176 : /**
177 : * Returns MeshBase::mesh_dimension(), (not
178 : * MeshBase::spatial_dimension()!) of the underlying libMesh mesh
179 : * object.
180 : */
181 : virtual unsigned int dimension() const;
182 :
183 : /**
184 : * Returns MeshBase::spatial_dimension
185 : */
186 56272 : virtual unsigned int spatialDimension() const { return _mesh->spatial_dimension(); }
187 :
188 : /**
189 : * Returns the effective spatial dimension determined by the coordinates actually used by the
190 : * mesh. This means that a 1D mesh that has non-zero z or y coordinates is actually a 2D or 3D
191 : * mesh, respectively. Likewise a 2D mesh that has non-zero z coordinates is actually 3D mesh.
192 : */
193 : virtual unsigned int effectiveSpatialDimension() const;
194 :
195 : /**
196 : * Returns the maximum element dimension on the given blocks
197 : */
198 : unsigned int getBlocksMaxDimension(const std::vector<SubdomainName> & blocks) const;
199 :
200 : /**
201 : * Returns a vector of boundary IDs for the requested element on the
202 : * requested side.
203 : */
204 : std::vector<BoundaryID> getBoundaryIDs(const Elem * const elem,
205 : const unsigned short int side) const;
206 :
207 : /**
208 : * Returns a const pointer to a lower dimensional element that
209 : * corresponds to a side of a higher dimensional element. This
210 : * relationship is established through an internal_parent; if there is
211 : * no lowerDElem, nullptr is returned.
212 : */
213 : const Elem * getLowerDElem(const Elem *, unsigned short int) const;
214 :
215 : /**
216 : * Returns the local side ID of the interior parent aligned with the lower dimensional element.
217 : */
218 : unsigned int getHigherDSide(const Elem * elem) const;
219 :
220 : /**
221 : * Returns a const reference to a set of all user-specified
222 : * boundary IDs. On a distributed mesh this will *only* include
223 : * boundary IDs which exist on local or ghosted elements; a copy and
224 : * a call to _communicator.set_union() will be necessary to get the
225 : * global ID set.
226 : */
227 : const std::set<BoundaryID> & getBoundaryIDs() const;
228 :
229 : /**
230 : * Calls BoundaryInfo::build_node_list()/build_side_list() and *makes separate copies* of
231 : * Nodes/Elems in those lists.
232 : *
233 : * Allocates memory which is cleaned up in the freeBndNodes()/freeBndElems() functions.
234 : */
235 : void buildNodeList();
236 : void buildBndElemList();
237 :
238 : /**
239 : * If not already created, creates a map from every node to all
240 : * elements to which they are connected.
241 : */
242 : const std::map<dof_id_type, std::vector<dof_id_type>> & nodeToElemMap();
243 :
244 : /**
245 : * If not already created, creates a map from every node to all
246 : * _active_ _semilocal_ elements to which they are connected.
247 : * Semilocal elements include local elements and elements that share at least
248 : * one node with a local element.
249 : * \note Extra ghosted elements are not included in this map!
250 : */
251 : const std::map<dof_id_type, std::vector<dof_id_type>> & nodeToActiveSemilocalElemMap();
252 :
253 : /**
254 : * These structs are required so that the bndNodes{Begin,End} and
255 : * bndElems{Begin,End} functions work...
256 : */
257 : struct bnd_node_iterator;
258 : struct const_bnd_node_iterator;
259 :
260 : struct bnd_elem_iterator;
261 : struct const_bnd_elem_iterator;
262 :
263 : /**
264 : * Return iterators to the beginning/end of the boundary nodes list.
265 : */
266 : virtual bnd_node_iterator bndNodesBegin();
267 : virtual bnd_node_iterator bndNodesEnd();
268 :
269 : /**
270 : * Return iterators to the beginning/end of the boundary elements list.
271 : */
272 : virtual bnd_elem_iterator bndElemsBegin();
273 : virtual bnd_elem_iterator bndElemsEnd();
274 :
275 : /**
276 : * Calls BoundaryInfo::build_node_list_from_side_list().
277 : */
278 : void buildNodeListFromSideList();
279 :
280 : /**
281 : * Calls BoundaryInfo::build_side_list(), returns a std::vector of
282 : * (elem-id, side-id, bc-id) tuples.
283 : */
284 : std::vector<std::tuple<dof_id_type, unsigned short int, boundary_id_type>> buildSideList();
285 :
286 : /**
287 : * Calls BoundaryInfo::build_active_side_list
288 : * @return A container of active (element, side, id) tuples.
289 : */
290 : std::vector<std::tuple<dof_id_type, unsigned short int, boundary_id_type>>
291 : buildActiveSideList() const;
292 :
293 : /**
294 : * Calls BoundaryInfo::side_with_boundary_id().
295 : */
296 : unsigned int sideWithBoundaryID(const Elem * const elem, const BoundaryID boundary_id) const;
297 :
298 : /**
299 : * Calls local_nodes_begin/end() on the underlying libMesh mesh object.
300 : */
301 : MeshBase::node_iterator localNodesBegin();
302 : MeshBase::node_iterator localNodesEnd();
303 : MeshBase::const_node_iterator localNodesBegin() const;
304 : MeshBase::const_node_iterator localNodesEnd() const;
305 :
306 : /**
307 : * Calls active_local_nodes_begin/end() on the underlying libMesh mesh object.
308 : */
309 : MeshBase::element_iterator activeLocalElementsBegin();
310 : const MeshBase::element_iterator activeLocalElementsEnd();
311 : MeshBase::const_element_iterator activeLocalElementsBegin() const;
312 : const MeshBase::const_element_iterator activeLocalElementsEnd() const;
313 :
314 : /**
315 : * Calls n_nodes/elem() on the underlying libMesh mesh object.
316 : */
317 : virtual dof_id_type nNodes() const;
318 : virtual dof_id_type nElem() const;
319 :
320 25876 : virtual dof_id_type nLocalNodes() const { return _mesh->n_local_nodes(); }
321 56635 : virtual dof_id_type nActiveElem() const { return _mesh->n_active_elem(); }
322 27358 : virtual dof_id_type nActiveLocalElem() const { return _mesh->n_active_local_elem(); }
323 58495 : virtual SubdomainID nSubdomains() const { return _mesh->n_subdomains(); }
324 25702 : virtual unsigned int nPartitions() const { return _mesh->n_partitions(); }
325 25702 : virtual bool skipPartitioning() const { return _mesh->skip_partitioning(); }
326 : virtual bool skipNoncriticalPartitioning() const;
327 :
328 : /**
329 : * Calls max_node/elem_id() on the underlying libMesh mesh object.
330 : * This may be larger than n_nodes/elem() in cases where the id
331 : * numbering is not contiguous.
332 : */
333 : virtual dof_id_type maxNodeId() const;
334 : virtual dof_id_type maxElemId() const;
335 :
336 : /**
337 : * Various accessors (pointers/references) for Node "i".
338 : *
339 : * If the requested node is a remote node on a distributed mesh,
340 : * only the query accessors are valid to call, and they return NULL.
341 : */
342 : virtual const Node & node(const dof_id_type i) const;
343 : virtual Node & node(const dof_id_type i);
344 : virtual const Node & nodeRef(const dof_id_type i) const;
345 : virtual Node & nodeRef(const dof_id_type i);
346 : virtual const Node * nodePtr(const dof_id_type i) const;
347 : virtual Node * nodePtr(const dof_id_type i);
348 : virtual const Node * queryNodePtr(const dof_id_type i) const;
349 : virtual Node * queryNodePtr(const dof_id_type i);
350 :
351 : /**
352 : * Various accessors (pointers/references) for Elem "i".
353 : *
354 : * If the requested elem is a remote element on a distributed mesh,
355 : * only the query accessors are valid to call, and they return NULL.
356 : */
357 : virtual Elem * elem(const dof_id_type i);
358 : virtual const Elem * elem(const dof_id_type i) const;
359 : virtual Elem * elemPtr(const dof_id_type i);
360 : virtual const Elem * elemPtr(const dof_id_type i) const;
361 : virtual Elem * queryElemPtr(const dof_id_type i);
362 : virtual const Elem * queryElemPtr(const dof_id_type i) const;
363 :
364 : /**
365 : * Setter/getter for whether the mesh is prepared
366 : */
367 : bool prepared() const;
368 : virtual void prepared(bool state);
369 :
370 : /**
371 : * If this method is called, we will call libMesh's prepare_for_use method when we
372 : * call Moose's prepare method. This should only be set when the mesh structure is changed
373 : * by MeshGenerators (i.e. Element deletion).
374 : */
375 : void needsPrepareForUse();
376 :
377 : /**
378 : * Declares that the MooseMesh has changed, invalidates cached data
379 : * and rebuilds caches. Sets a flag so that clients of the
380 : * MooseMesh also know when it has changed.
381 : */
382 : void meshChanged();
383 :
384 : /**
385 : * Declares a callback function that is executed at the conclusion
386 : * of meshChanged(). Ther user can implement actions required after
387 : * changing the mesh here.
388 : **/
389 : virtual void onMeshChanged();
390 :
391 : /**
392 : * Cache information about what elements were refined and coarsened in the previous step.
393 : */
394 : void cacheChangedLists();
395 :
396 : /**
397 : * Return a range that is suitable for threaded execution over elements that were just refined.
398 : *
399 : * @return The _Parent_ elements that are now set to be INACTIVE. Their _children_ are the new
400 : * elements.
401 : */
402 : ConstElemPointerRange * refinedElementRange() const;
403 :
404 : /**
405 : * Return a range that is suitable for threaded execution over elements that were just coarsened.
406 : * Note that these are the _Parent_ elements that are now set to be INACTIVE. Their _children_
407 : * are the elements that were just removed. Use coarsenedElementChildren() to get the element
408 : * IDs for the children that were just removed for a particular parent element.
409 : */
410 : ConstElemPointerRange * coarsenedElementRange() const;
411 :
412 : /**
413 : * Get the newly removed children element ids for an element that was just coarsened.
414 : *
415 : * @param elem Pointer to the parent element that was coarsened to.
416 : * @return The child element ids in Elem::child() order.
417 : */
418 : const std::vector<const Elem *> & coarsenedElementChildren(const Elem * elem) const;
419 :
420 : /**
421 : * Clears the "semi-local" node list and rebuilds it. Semi-local nodes
422 : * consist of all nodes that belong to local and ghost elements.
423 : */
424 : void updateActiveSemiLocalNodeRange(std::set<dof_id_type> & ghosted_elems);
425 :
426 : /**
427 : * Returns true if the node is semi-local
428 : * @param node Node pointer
429 : * @return true is the node is semi-local, false otherwise
430 : */
431 : bool isSemiLocal(Node * const node) const;
432 :
433 : ///@{
434 : /**
435 : * Return pointers to range objects for various types of ranges
436 : * (local nodes, boundary elems, etc.).
437 : */
438 : libMesh::ConstElemRange * getActiveLocalElementRange();
439 : libMesh::NodeRange * getActiveNodeRange();
440 : SemiLocalNodeRange * getActiveSemiLocalNodeRange() const;
441 : libMesh::ConstNodeRange * getLocalNodeRange();
442 : libMesh::StoredRange<MooseMesh::const_bnd_node_iterator, const BndNode *> *
443 : getBoundaryNodeRange();
444 : libMesh::StoredRange<MooseMesh::const_bnd_elem_iterator, const BndElement *> *
445 : getBoundaryElementRange();
446 : ///@}
447 :
448 : /**
449 : * Returns a map of boundaries to ids of elements on the boundary.
450 : */
451 : const std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> &
452 : getBoundariesToElems() const;
453 :
454 : /**
455 : * Returns a map of boundaries to ids of elements on the boundary.
456 : */
457 : const std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> &
458 : getBoundariesToActiveSemiLocalElemIds() const;
459 :
460 : /**
461 : * Return all ids of elements which have a side which is part of a sideset.
462 : * Note that boundaries are sided.
463 : * @param bid the id of the sideset of interest
464 : */
465 : std::unordered_set<dof_id_type> getBoundaryActiveSemiLocalElemIds(BoundaryID bid) const;
466 :
467 : /**
468 : * Return all ids of neighbors of elements which have a side which is part of a sideset.
469 : * Note that boundaries are sided, this is on the neighbor side. For the sideset side, use
470 : * getBoundariesActiveLocalElemIds.
471 : * Note that while the element is local and active, the neighbor is not guaranteed to be local,
472 : * it could be ghosted.
473 : * Note that if the neighbor is not ghosted, is a remote_elem, then it will not be included
474 : * @param bid the id of the sideset of interest
475 : */
476 : std::unordered_set<dof_id_type> getBoundaryActiveNeighborElemIds(BoundaryID bid) const;
477 :
478 : /**
479 : * Returns whether a boundary (given by its id) is not crossing through a group of blocks,
480 : * by which we mean that elements on both sides of the boundary are in those blocks
481 : * @param bid the id of the boundary of interest
482 : * @param blk_group the group of blocks potentially traversed
483 : * @return whether the boundary does not cross between the subdomains in the group
484 : */
485 : bool isBoundaryFullyExternalToSubdomains(BoundaryID bid,
486 : const std::set<SubdomainID> & blk_group) const;
487 :
488 : /**
489 : * Returns a read-only reference to the set of subdomains currently
490 : * present in the Mesh.
491 : */
492 : const std::set<SubdomainID> & meshSubdomains() const;
493 :
494 : /**
495 : * Returns a read-only reference to the set of boundary IDs currently
496 : * present in the Mesh.
497 : */
498 : const std::set<BoundaryID> & meshBoundaryIds() const;
499 :
500 : /**
501 : * Returns a read-only reference to the set of sidesets currently
502 : * present in the Mesh.
503 : */
504 : const std::set<BoundaryID> & meshSidesetIds() const;
505 :
506 : /**
507 : * Returns a read-only reference to the set of nodesets currently
508 : * present in the Mesh.
509 : */
510 : const std::set<BoundaryID> & meshNodesetIds() const;
511 :
512 : /**
513 : * Sets the mapping between BoundaryID and normal vector
514 : * Is called by AddAllSideSetsByNormals
515 : */
516 : void setBoundaryToNormalMap(std::unique_ptr<std::map<BoundaryID, RealVectorValue>> boundary_map);
517 :
518 : // DEPRECATED METHOD
519 : void setBoundaryToNormalMap(std::map<BoundaryID, RealVectorValue> * boundary_map);
520 :
521 : /**
522 : * Sets the set of BoundaryIDs
523 : * Is called by AddAllSideSetsByNormals
524 : */
525 : void setMeshBoundaryIDs(std::set<BoundaryID> boundary_IDs);
526 :
527 : /**
528 : * Returns the normal vector associated with a given BoundaryID.
529 : * It's only valid to call this when AddAllSideSetsByNormals is active.
530 : */
531 : const RealVectorValue & getNormalByBoundaryID(BoundaryID id) const;
532 :
533 : /**
534 : * Calls prepare_for_use() if the underlying MeshBase object isn't prepared, then communicates
535 : * various boundary information on parallel meshes. Also calls update() internally. Instead of
536 : * calling \p prepare_for_use on the currently held \p MeshBase object, a \p mesh_to_clone can be
537 : * provided. If it is provided (e.g. this method is given a non-null argument), then \p _mesh will
538 : * be assigned a clone of the \p mesh_to_clone. The provided \p mesh_to_clone must already be
539 : * prepared
540 : * @param mesh_to_clone If nonnull, we will clone this mesh instead of preparing our current one
541 : * @return Whether the libMesh mesh was prepared. This should really only be relevant in MOOSE
542 : * framework contexts where we need to make a decision about what to do with the displaced mesh.
543 : * If the reference mesh base object has \p prepare_for_use called (e.g. this method returns \p
544 : * true when called for the reference mesh), then we must pass the reference mesh base object into
545 : * this method when we call this for the displaced mesh. This is because the displaced mesh \emph
546 : * must be an exact clone of the reference mesh. We have seen that \p prepare_for_use called on
547 : * two previously identical meshes can result in two different meshes even with Metis partitioning
548 : */
549 : bool prepare(const MeshBase * mesh_to_clone);
550 :
551 : /**
552 : * Calls buildNodeListFromSideList(), buildNodeList(), and buildBndElemList().
553 : */
554 : void update();
555 :
556 : /**
557 : * Returns the level of uniform refinement requested (zero if AMR is disabled).
558 : */
559 : unsigned int uniformRefineLevel() const;
560 :
561 : /**
562 : * Set uniform refinement level
563 : */
564 : void setUniformRefineLevel(unsigned int, bool deletion = true);
565 :
566 : /**
567 : * Return a flag indicating whether or not we should skip remote deletion
568 : * and repartition after uniform refinements. If the flag is true, uniform
569 : * refinements will run more efficiently, but at the same time, there might
570 : * be extra ghosting elements. The number of layers of additional ghosting
571 : * elements depends on the number of uniform refinement levels. This flag
572 : * should be used only when you have a "fine enough" coarse mesh and want
573 : * to refine the mesh by a few levels. Otherwise, it might introduce an
574 : * unbalanced workload and too large ghosting domain.
575 : */
576 : bool skipDeletionRepartitionAfterRefine() const;
577 :
578 : /**
579 : * Whether or not skip uniform refinements when using a pre-split mesh
580 : */
581 294 : bool skipRefineWhenUseSplit() const { return _skip_refine_when_use_split; }
582 :
583 : /**
584 : * This will add the boundary ids to be ghosted to this processor
585 : */
586 : void addGhostedBoundary(BoundaryID boundary_id);
587 :
588 : /**
589 : * This sets the inflation amount for the bounding box for each partition for use in
590 : * ghosting boundaries
591 : */
592 : void setGhostedBoundaryInflation(const std::vector<Real> & inflation);
593 :
594 : /**
595 : * Return a writable reference to the set of ghosted boundary IDs.
596 : */
597 : const std::set<unsigned int> & getGhostedBoundaries() const;
598 :
599 : /**
600 : * Return a writable reference to the _ghosted_boundaries_inflation vector.
601 : */
602 : const std::vector<Real> & getGhostedBoundaryInflation() const;
603 :
604 : /**
605 : * Actually do the ghosting of boundaries that need to be ghosted to this processor.
606 : */
607 : void ghostGhostedBoundaries();
608 :
609 : /**
610 : * Whether or not we want to ghost ghosted boundaries
611 : */
612 879 : void needGhostGhostedBoundaries(bool needghost) { _need_ghost_ghosted_boundaries = needghost; }
613 :
614 : /**
615 : * Getter for the patch_size parameter.
616 : */
617 : unsigned int getPatchSize() const;
618 :
619 : /**
620 : * Getter for the ghosting_patch_size parameter.
621 : */
622 9442 : unsigned int getGhostingPatchSize() const { return _ghosting_patch_size; }
623 :
624 : /**
625 : * Getter for the maximum leaf size parameter.
626 : */
627 68972 : unsigned int getMaxLeafSize() const { return _max_leaf_size; }
628 :
629 : /**
630 : * Set the patch size update strategy
631 : */
632 : void setPatchUpdateStrategy(Moose::PatchUpdateType patch_update_strategy);
633 :
634 : /**
635 : * Get the current patch update strategy.
636 : */
637 : const Moose::PatchUpdateType & getPatchUpdateStrategy() const;
638 :
639 : /**
640 : * Get a (slightly inflated) processor bounding box.
641 : *
642 : * @param inflation_multiplier This amount will be multiplied by the length of the diagonal of the
643 : * bounding box to find the amount to inflate the bounding box by in all directions.
644 : */
645 : libMesh::BoundingBox getInflatedProcessorBoundingBox(Real inflation_multiplier = 0.01) const;
646 :
647 : /**
648 : * Implicit conversion operator from MooseMesh -> libMesh::MeshBase.
649 : */
650 : operator libMesh::MeshBase &();
651 : operator const libMesh::MeshBase &() const;
652 :
653 : /**
654 : * Accessor for the underlying libMesh Mesh object.
655 : */
656 : MeshBase & getMesh();
657 : MeshBase & getMesh(const std::string & name);
658 : const MeshBase & getMesh() const;
659 : const MeshBase & getMesh(const std::string & name) const;
660 : const MeshBase * getMeshPtr() const;
661 :
662 : /**
663 : * Accessor for Kokkos mesh object.
664 : */
665 : #ifdef MOOSE_KOKKOS_ENABLED
666 57048 : const Moose::Kokkos::Mesh * getKokkosMesh() const { return _kokkos_mesh.get(); }
667 : #endif
668 :
669 : /**
670 : * Calls print_info() on the underlying Mesh.
671 : */
672 : void printInfo(std::ostream & os = libMesh::out, const unsigned int verbosity = 0) const;
673 :
674 : /**
675 : * Return list of blocks to which the given node belongs.
676 : */
677 : const std::set<SubdomainID> & getNodeBlockIds(const Node & node) const;
678 :
679 : /**
680 : * Return a writable reference to a vector of node IDs that belong
681 : * to nodeset_id.
682 : */
683 : const std::vector<dof_id_type> & getNodeList(boundary_id_type nodeset_id) const;
684 :
685 : /**
686 : * Add a new node to the mesh. If there is already a node located at the point passed
687 : * then the node will not be added. In either case a reference to the node at that location
688 : * will be returned
689 : */
690 : const Node * addUniqueNode(const Point & p, Real tol = 1e-6);
691 :
692 : /**
693 : * Adds a fictitious "QuadratureNode". This doesn't actually add it to the libMesh mesh...
694 : * we just keep track of these here in MooseMesh.
695 : *
696 : * QuadratureNodes are fictitious "Nodes" that are located at quadrature points. This is useful
697 : * for using the geometric search system to do searches based on quadrature point locations....
698 : *
699 : * @param elem The element
700 : * @param side The side number on which we want to add a quadrature node
701 : * @param qp The number of the quadrature point
702 : * @param bid The boundary ID for the point to be added with
703 : * @param point The physical location of the point
704 : */
705 : Node * addQuadratureNode(const Elem * elem,
706 : const unsigned short int side,
707 : const unsigned int qp,
708 : BoundaryID bid,
709 : const Point & point);
710 :
711 : /**
712 : * Get a specified quadrature node.
713 : *
714 : * @param elem The element the quadrature point is on
715 : * @param side The side the quadrature point is on
716 : * @param qp The quadrature point number associated with the point
717 : */
718 : Node * getQuadratureNode(const Elem * elem, const unsigned short int side, const unsigned int qp);
719 :
720 : /**
721 : * Clear out any existing quadrature nodes.
722 : * Most likely called before re-adding them.
723 : */
724 : void clearQuadratureNodes();
725 :
726 : /**
727 : * Get the associated BoundaryID for the boundary name.
728 : *
729 : * @param boundary_name The name of the boundary.
730 : * @return the boundary id from the passed boundary name.
731 : */
732 : BoundaryID getBoundaryID(const BoundaryName & boundary_name) const;
733 :
734 : /**
735 : * Get the associated BoundaryID for the boundary names that are passed in.
736 : *
737 : * @param boundary_name The names of the boundaries.
738 : * @return the boundary ids from the passed boundary names.
739 : */
740 : std::vector<BoundaryID> getBoundaryIDs(const std::vector<BoundaryName> & boundary_name,
741 : bool generate_unknown = false) const;
742 :
743 : /**
744 : * Get the associated subdomain ID for the subdomain name.
745 : *
746 : * @param subdomain_name The name of the subdomain
747 : * @return The subdomain id from the passed subdomain name.
748 : */
749 : SubdomainID getSubdomainID(const SubdomainName & subdomain_name) const;
750 :
751 : /**
752 : * Get the associated subdomainIDs for the subdomain names that are passed in.
753 : *
754 : * @param subdomain_names The names of the subdomains
755 : * @return The subdomain ids from the passed subdomain names.
756 : */
757 : std::vector<SubdomainID>
758 : getSubdomainIDs(const std::vector<SubdomainName> & subdomain_names) const;
759 : std::set<SubdomainID> getSubdomainIDs(const std::set<SubdomainName> & subdomain_names) const;
760 :
761 : /**
762 : * This method sets the name for \p subdomain_id to \p name
763 : */
764 : void setSubdomainName(SubdomainID subdomain_id, const SubdomainName & name);
765 :
766 : /**
767 : * This method sets the name for \p subdomain_id on the provided \p mesh to \p name
768 : */
769 : static void
770 : setSubdomainName(MeshBase & mesh, SubdomainID subdomain_id, const SubdomainName & name);
771 :
772 : /**
773 : * Return the name of a block given an id.
774 : */
775 : const std::string & getSubdomainName(SubdomainID subdomain_id) const;
776 :
777 : /**
778 : * Get the associated subdomainNames for the subdomain ids that are passed in.
779 : *
780 : * @param subdomain_ids The ids of the subdomains
781 : * @return The subdomain names from the passed subdomain ids.
782 : */
783 : std::vector<SubdomainName>
784 : getSubdomainNames(const std::vector<SubdomainID> & subdomain_ids) const;
785 :
786 : /**
787 : * This method sets the boundary name of the boundary based on the id parameter
788 : */
789 : void setBoundaryName(BoundaryID boundary_id, BoundaryName name);
790 :
791 : /**
792 : * Return the name of the boundary given the id.
793 : */
794 : const std::string & getBoundaryName(BoundaryID boundary_id) const;
795 :
796 : /**
797 : * This routine builds a multimap of boundary ids to matching boundary ids across all periodic
798 : * boundaries
799 : * in the system.
800 : */
801 : void buildPeriodicNodeMap(std::multimap<dof_id_type, dof_id_type> & periodic_node_map,
802 : unsigned int var_number,
803 : libMesh::PeriodicBoundaries * pbs) const;
804 :
805 : /**
806 : * This routine builds a datastructure of node ids organized by periodic boundary ids
807 : */
808 : void buildPeriodicNodeSets(std::map<BoundaryID, std::set<dof_id_type>> & periodic_node_sets,
809 : unsigned int var_number,
810 : libMesh::PeriodicBoundaries * pbs) const;
811 :
812 : /**
813 : * Returns the width of the requested dimension
814 : */
815 : Real dimensionWidth(unsigned int component) const;
816 :
817 : ///@{
818 : /**
819 : * Returns the min or max of the requested dimension respectively
820 : */
821 : virtual Real getMinInDimension(unsigned int component) const;
822 : virtual Real getMaxInDimension(unsigned int component) const;
823 : ///@}
824 :
825 : /**
826 : * This routine determines whether the Mesh is a regular orthogonal mesh (i.e. square in 2D, cubic
827 : * in 3D). If it is, then we can use a number of convenience functions when periodic boundary
828 : * conditions are applied. This routine populates the _range vector which is necessary for these
829 : * convenience functions.
830 : *
831 : * Note: This routine can potentially identify meshes with concave faces that still "fit" in the
832 : * convex hull of the corresponding regular orthogonal mesh. This case is highly unlikely in
833 : * practice and if a user does this, well.... release the kicker!
834 : */
835 : bool detectOrthogonalDimRanges(Real tol = 1e-6);
836 :
837 : /**
838 : * For "regular orthogonal" meshes, determine if variable var_num is periodic with respect to the
839 : * primary and secondary BoundaryIDs, record this fact in the _periodic_dim data structure.
840 : */
841 : void addPeriodicVariable(unsigned int var_num, BoundaryID primary, BoundaryID secondary);
842 :
843 : /**
844 : * Returns whether this generated mesh is periodic in the given dimension for the given variable.
845 : * @param nonlinear_var_num - The nonlinear variable number
846 : * @param component - An integer representing the desired component (dimension)
847 : */
848 : bool isTranslatedPeriodic(unsigned int nonlinear_var_num, unsigned int component) const;
849 :
850 : /**
851 : * This function returns the minimum vector between two points on the mesh taking into account
852 : * periodicity for the given variable number.
853 : * @param nonlinear_var_num - The nonlinear variable number
854 : * @param p, q - The points between which to compute a minimum vector
855 : * @return RealVectorValue - The vector pointing from p to q
856 : */
857 : RealVectorValue minPeriodicVector(unsigned int nonlinear_var_num, Point p, Point q) const;
858 :
859 : /**
860 : * This function returns the distance between two points on the mesh taking into account
861 : * periodicity for the given variable number.
862 : * @param nonlinear_var_num - The nonlinear variable number
863 : * @param p, q - The points for which to compute a minimum distance
864 : * @return Real - The L2 distance between p and q
865 : */
866 : Real minPeriodicDistance(unsigned int nonlinear_var_num, Point p, Point q) const;
867 :
868 : /**
869 : * This function attempts to return the paired boundary ids for the given component. For example,
870 : * in a generated 2D mesh, passing 0 for the "x" component will return (3, 1).
871 : * @param component - An integer representing the desired component (dimension)
872 : * @return std::pair pointer - The matching boundary pairs for the passed component
873 : */
874 : const std::pair<BoundaryID, BoundaryID> * getPairedBoundaryMapping(unsigned int component);
875 :
876 : /**
877 : * Create the refinement and coarsening maps necessary for projection of stateful material
878 : * properties when using adaptivity.
879 : *
880 : * @param assembly Pointer to the Assembly object for this Mesh.
881 : */
882 : void buildRefinementAndCoarseningMaps(Assembly * assembly);
883 :
884 : /**
885 : * Get the refinement map for a given element type. This will tell you what quadrature points
886 : * to copy from and to for stateful material properties on newly created elements from Adaptivity.
887 : *
888 : * @param elem The element that represents the element type you need the refinement map for.
889 : * @param parent_side The side of the parent to map (-1 if not mapping parent sides)
890 : * @param child The child number (-1 if not mapping child internal sides)
891 : * @param child_side The side number of the child (-1 if not mapping sides)
892 : */
893 : const std::vector<std::vector<QpMap>> &
894 : getRefinementMap(const Elem & elem, int parent_side, int child, int child_side);
895 :
896 : /**
897 : * Get the coarsening map for a given element type. This will tell you what quadrature points
898 : * to copy from and to for stateful material properties on newly created elements from Adaptivity.
899 : *
900 : * @param elem The element that represents the element type you need the coarsening map for.
901 : * @param input_side The side to map
902 : */
903 : const std::vector<std::pair<unsigned int, QpMap>> & getCoarseningMap(const Elem & elem,
904 : int input_side);
905 :
906 : /**
907 : * Change all the boundary IDs for a given side from old_id to new_id. If delete_prev is true,
908 : * also actually remove the side with old_id from the BoundaryInfo object.
909 : */
910 : void
911 : changeBoundaryId(const boundary_id_type old_id, const boundary_id_type new_id, bool delete_prev);
912 :
913 : /**
914 : * Change all the boundary IDs for a given side from old_id to new_id for the given \p mesh. If
915 : * delete_prev is true, also actually remove the side with old_id from the BoundaryInfo object.
916 : */
917 : static void changeBoundaryId(MeshBase & mesh,
918 : const boundary_id_type old_id,
919 : const boundary_id_type new_id,
920 : bool delete_prev);
921 :
922 : /**
923 : * Get the list of boundary ids associated with the given subdomain id.
924 : *
925 : * @param subdomain_id The subdomain ID you want to get the boundary ids for.
926 : * @return All boundary IDs connected to elements in the give
927 : */
928 : const std::set<BoundaryID> & getSubdomainBoundaryIds(const SubdomainID subdomain_id) const;
929 :
930 : /**
931 : * Get the list of boundaries that contact the given subdomain.
932 : *
933 : * @param subdomain_id The subdomain ID you want to get the boundary ids for.
934 : * @return All boundary IDs connected to elements in the given subdomain
935 : */
936 : std::set<BoundaryID> getSubdomainInterfaceBoundaryIds(const SubdomainID subdomain_id) const;
937 :
938 : /**
939 : * Get the list of subdomains associated with the given boundary.
940 : *
941 : * @param bid The boundary ID you want to get the subdomain IDs for.
942 : * @return All subdomain IDs associated with given boundary ID
943 : */
944 : std::set<SubdomainID> getBoundaryConnectedBlocks(const BoundaryID bid) const;
945 :
946 : /**
947 : * Get the list of subdomains associated with the given boundary of its secondary side.
948 : *
949 : * @param bid The boundary ID you want to get the subdomain IDs for.
950 : * @return All subdomain IDs associated with given boundary ID
951 : */
952 : std::set<SubdomainID> getBoundaryConnectedSecondaryBlocks(const BoundaryID bid) const;
953 :
954 : /**
955 : * Get the list of subdomains contacting the given boundary.
956 : *
957 : * @param bid The boundary ID you want to get the subdomain IDs for.
958 : * @return All subdomain IDs contacting given boundary ID
959 : */
960 : std::set<SubdomainID> getInterfaceConnectedBlocks(const BoundaryID bid) const;
961 :
962 : /**
963 : * Get the list of subdomains neighboring a given subdomain.
964 : *
965 : * @param subdomain_id The boundary ID you want to get the subdomain IDs for.
966 : * @return All subdomain IDs neighboring a given subdomain
967 : */
968 : const std::set<SubdomainID> & getBlockConnectedBlocks(const SubdomainID subdomain_id) const;
969 :
970 : /**
971 : * Returns true if the requested node is in the list of boundary nodes, false otherwise.
972 : */
973 : bool isBoundaryNode(dof_id_type node_id) const;
974 :
975 : /**
976 : * Returns true if the requested node is in the list of boundary nodes for the specified boundary,
977 : * false otherwise.
978 : */
979 : bool isBoundaryNode(dof_id_type node_id, BoundaryID bnd_id) const;
980 :
981 : /**
982 : * Returns true if the requested element is in the list of boundary elements, false otherwise.
983 : */
984 : bool isBoundaryElem(dof_id_type elem_id) const;
985 :
986 : /**
987 : * Returns true if the requested element is in the list of boundary elements for the specified
988 : * boundary, false otherwise.
989 : */
990 : bool isBoundaryElem(dof_id_type elem_id, BoundaryID bnd_id) const;
991 :
992 : /**
993 : * Generate a unified error message if the underlying libMesh mesh is a DistributedMesh. Clients
994 : * of MooseMesh can use this function to throw an error if they know they don't work with
995 : * DistributedMesh.
996 : *
997 : * See, for example, the NodalVariableValue class.
998 : */
999 : void errorIfDistributedMesh(std::string name) const;
1000 :
1001 : /**
1002 : * Returns the final Mesh distribution type.
1003 : */
1004 67253 : virtual bool isDistributedMesh() const { return _use_distributed_mesh; }
1005 :
1006 : /**
1007 : * Tell the user if the distribution was overriden for any reason
1008 : */
1009 56272 : bool isParallelTypeForced() const { return _parallel_type_overridden; }
1010 :
1011 : /**
1012 : * Allow to change parallel type
1013 : */
1014 : void setParallelType(ParallelType parallel_type);
1015 :
1016 : /**
1017 : * @return The parallel type
1018 : */
1019 1242 : ParallelType getParallelType() const { return _parallel_type; }
1020 :
1021 : /*
1022 : * Set/Get the partitioner name
1023 : */
1024 25702 : const MooseEnum & partitionerName() const { return _partitioner_name; }
1025 :
1026 : /**
1027 : * Tell the user if the partitioner was overriden for any reason
1028 : */
1029 25702 : bool isPartitionerForced() const { return _partitioner_overridden; }
1030 :
1031 : /**
1032 : * Set whether or not this mesh is allowed to read a recovery file.
1033 : */
1034 11 : void allowRecovery(bool allow) { _allow_recovery = allow; }
1035 :
1036 : /**
1037 : * Method for setting the partitioner on the passed in mesh_base object.
1038 : */
1039 : static void setPartitioner(MeshBase & mesh_base,
1040 : MooseEnum & partitioner,
1041 : bool use_distributed_mesh,
1042 : const InputParameters & params,
1043 : MooseObject & context_obj);
1044 :
1045 : /**
1046 : * Setter for custom partitioner
1047 : */
1048 : void setCustomPartitioner(libMesh::Partitioner * partitioner);
1049 :
1050 : ///@{
1051 : /**
1052 : * Setter and getter for _custom_partitioner_requested
1053 : */
1054 : bool isCustomPartitionerRequested() const;
1055 : void setIsCustomPartitionerRequested(bool cpr);
1056 : ///@}
1057 :
1058 : /// Getter to query if the mesh was detected to be regular and orthogonal
1059 1212 : bool isRegularOrthogonal() { return _regular_orthogonal_mesh; }
1060 :
1061 : /// check if the mesh has SECOND order elements
1062 : bool hasSecondOrderElements();
1063 :
1064 : /**
1065 : * Proxy function to get a (sub)PointLocator from either the underlying libMesh mesh (default), or
1066 : * to allow derived meshes to return a custom point locator.
1067 : */
1068 : virtual std::unique_ptr<libMesh::PointLocatorBase> getPointLocator() const;
1069 :
1070 : /**
1071 : * Returns the name of the mesh file read to produce this mesh if any or an empty string
1072 : * otherwise.
1073 : */
1074 522 : virtual std::string getFileName() const { return ""; }
1075 :
1076 : /// Helper type for building periodic node maps
1077 : using PeriodicNodeInfo = std::pair<const Node *, BoundaryID>;
1078 :
1079 : /**
1080 : * Set whether we need to delete remote elements
1081 : */
1082 27 : void needsRemoteElemDeletion(bool need_delete) { _need_delete = need_delete; }
1083 :
1084 : /**
1085 : * Whether we need to delete remote elements
1086 : */
1087 67079 : bool needsRemoteElemDeletion() const { return _need_delete; }
1088 :
1089 : /**
1090 : * Set whether to allow remote element removal
1091 : */
1092 : void allowRemoteElementRemoval(bool allow_removal);
1093 :
1094 : /**
1095 : * Whether we are allow remote element removal
1096 : */
1097 29145 : bool allowRemoteElementRemoval() const { return _allow_remote_element_removal; }
1098 :
1099 : /**
1100 : * Delete remote elements
1101 : */
1102 : void deleteRemoteElements();
1103 :
1104 : /**
1105 : * Whether mesh base object was constructed or not
1106 : */
1107 91682 : bool hasMeshBase() const { return _mesh.get() != nullptr; }
1108 :
1109 : /**
1110 : * Whether mesh has an extra element integer with a given name
1111 : */
1112 : bool hasElementID(const std::string & id_name) const;
1113 :
1114 : /**
1115 : * Return the accessing integer for an extra element integer with its name
1116 : */
1117 : unsigned int getElementIDIndex(const std::string & id_name) const;
1118 :
1119 : /**
1120 : * Return the maximum element ID for an extra element integer with its accessing index
1121 : */
1122 : dof_id_type maxElementID(unsigned int elem_id_index) const { return _max_ids[elem_id_index]; }
1123 :
1124 : /**
1125 : * Return the minimum element ID for an extra element integer with its accessing index
1126 : */
1127 : dof_id_type minElementID(unsigned int elem_id_index) const { return _min_ids[elem_id_index]; }
1128 :
1129 : /**
1130 : * Whether or not two extra element integers are identical
1131 : */
1132 : bool areElemIDsIdentical(const std::string & id_name1, const std::string & id_name2) const;
1133 :
1134 : /**
1135 : * Return all the unique element IDs for an extra element integer with its index
1136 : */
1137 : std::set<dof_id_type> getAllElemIDs(unsigned int elem_id_index) const;
1138 :
1139 : /**
1140 : * Return all the unique element IDs for an extra element integer with its index on a set of
1141 : * subdomains
1142 : */
1143 : std::set<dof_id_type> getElemIDsOnBlocks(unsigned int elem_id_index,
1144 : const std::set<SubdomainID> & blks) const;
1145 :
1146 : /**
1147 : * Get the maximum number of sides per element
1148 : */
1149 4458 : unsigned int getMaxSidesPerElem() const { return _max_sides_per_elem; }
1150 :
1151 : /**
1152 : * Get the maximum number of nodes per element
1153 : */
1154 1488 : unsigned int getMaxNodesPerElem() const { return _max_nodes_per_elem; }
1155 :
1156 : /**
1157 : * Get the maximum number of nodes per side
1158 : */
1159 1488 : unsigned int getMaxNodesPerSide() const { return _max_nodes_per_side; }
1160 :
1161 : std::unordered_map<dof_id_type, std::set<dof_id_type>>
1162 : getElemIDMapping(const std::string & from_id_name, const std::string & to_id_name) const;
1163 :
1164 : ///@{ accessors for the FaceInfo objects
1165 : unsigned int nFace() const { return _face_info.size(); }
1166 :
1167 : /// Accessor for local \p FaceInfo objects.
1168 : const std::vector<const FaceInfo *> & faceInfo() const;
1169 :
1170 : /// Need to declare these iterators here to make sure the iterators below work
1171 : struct face_info_iterator;
1172 : struct const_face_info_iterator;
1173 :
1174 : /// Iterators to owned faceInfo objects. These faceInfo-s are required for the
1175 : /// face loops and to filter out the faceInfo-s that are not owned by this processor
1176 : /// in case we have a distributed mesh and we included FaceInfo objects that
1177 : /// are on processor boundaries
1178 : face_info_iterator ownedFaceInfoBegin();
1179 : face_info_iterator ownedFaceInfoEnd();
1180 :
1181 : /// Need to declare these iterators here to make sure the iterators below work
1182 : struct elem_info_iterator;
1183 : struct const_elem_info_iterator;
1184 :
1185 : /// Iterators to owned faceInfo objects. These faceInfo-s are required for the
1186 : /// face loops and to filter out the faceInfo-s that are not owned by this processor
1187 : /// in case we have a distributed mesh and we included FaceInfo objects that
1188 : /// are on processor boundaries
1189 : elem_info_iterator ownedElemInfoBegin();
1190 : elem_info_iterator ownedElemInfoEnd();
1191 :
1192 : /// Accessor for the local FaceInfo object on the side of one element. Returns null if ghosted.
1193 : const FaceInfo * faceInfo(const Elem * elem, unsigned int side) const;
1194 :
1195 : /// Accessor for the elemInfo object for a given element ID
1196 : const ElemInfo & elemInfo(const dof_id_type id) const;
1197 :
1198 : /// Accessor for the element info objects owned by this process
1199 24 : const std::vector<const ElemInfo *> & elemInfoVector() const { return _elem_info; }
1200 :
1201 : /// Accessor for all \p FaceInfo objects.
1202 : const std::vector<FaceInfo> & allFaceInfo() const;
1203 : ///@}
1204 :
1205 : /**
1206 : * Cache if variables live on the elements connected by the FaceInfo objects
1207 : */
1208 : void cacheFaceInfoVariableOwnership() const;
1209 :
1210 : /**
1211 : * Cache the DoF indices for FV variables on each element. These indices are used to speed up the
1212 : * setup loops of finite volume systems.
1213 : */
1214 : void cacheFVElementalDoFs() const;
1215 :
1216 : /**
1217 : * Compute the face coordinate value for all \p FaceInfo and \p ElemInfo objects. 'Coordinate'
1218 : * here means a coordinate value associated with the coordinate system. For Cartesian coordinate
1219 : * systems, 'coordinate' is simply '1'; in RZ, '2*pi*r', and in spherical, '4*pi*r^2'
1220 : */
1221 : void computeFiniteVolumeCoords() const;
1222 :
1223 : /**
1224 : * Set whether this mesh is a displaced mesh
1225 : */
1226 2257 : void isDisplaced(bool is_displaced) { _is_displaced = is_displaced; }
1227 :
1228 : /**
1229 : * whether this mesh is a displaced mesh
1230 : */
1231 : bool isDisplaced() const { return _is_displaced; }
1232 :
1233 : /**
1234 : * @return A map from nodeset ids to the vector of node ids in the nodeset
1235 : */
1236 : const std::map<boundary_id_type, std::vector<dof_id_type>> & nodeSetNodes() const;
1237 :
1238 : /**
1239 : * Get the coordinate system type, e.g. xyz, rz, or r-spherical, for the provided subdomain id \p
1240 : * sid
1241 : */
1242 : Moose::CoordinateSystemType getCoordSystem(SubdomainID sid) const;
1243 :
1244 : /**
1245 : * Get the coordinate system from the mesh, it must be the same in all subdomains otherwise this
1246 : * will error
1247 : */
1248 : Moose::CoordinateSystemType getUniqueCoordSystem() const;
1249 :
1250 : /**
1251 : * Get the map from subdomain ID to coordinate system type, e.g. xyz, rz, or r-spherical
1252 : */
1253 : const std::map<SubdomainID, Moose::CoordinateSystemType> & getCoordSystem() const;
1254 :
1255 : /**
1256 : * Set the coordinate system for the provided blocks to \p coord_sys
1257 : */
1258 : void setCoordSystem(const std::vector<SubdomainName> & blocks, const MultiMooseEnum & coord_sys);
1259 :
1260 : /**
1261 : * For axisymmetric simulations, set the symmetry coordinate axis. For r in the x-direction, z in
1262 : * the y-direction the coordinate axis would be y
1263 : */
1264 : void setAxisymmetricCoordAxis(const MooseEnum & rz_coord_axis);
1265 :
1266 : /**
1267 : * Sets the general coordinate axes for axisymmetric blocks.
1268 : *
1269 : * This method must be used if any of the following are true:
1270 : * - There are multiple axisymmetric coordinate systems
1271 : * - Any axisymmetric coordinate system axis/direction is not the +X or +Y axis
1272 : * - Any axisymmetric coordinate system does not start at (0,0,0)
1273 : *
1274 : * @param[in] blocks Subdomain names
1275 : * @param[in] axes Pair of values defining the axisymmetric coordinate axis
1276 : * for each subdomain. The first value is the point on the axis
1277 : * corresponding to the origin. The second value is the direction
1278 : * vector of the axis (normalization not necessary).
1279 : */
1280 : void setGeneralAxisymmetricCoordAxes(const std::vector<SubdomainName> & blocks,
1281 : const std::vector<std::pair<Point, RealVectorValue>> & axes);
1282 :
1283 : /**
1284 : * Gets the general axisymmetric coordinate axis for a block.
1285 : *
1286 : * @param[in] subdomain_id Subdomain ID for which to get axisymmetric coordinate axis
1287 : */
1288 : const std::pair<Point, RealVectorValue> &
1289 : getGeneralAxisymmetricCoordAxis(SubdomainID subdomain_id) const;
1290 :
1291 : /**
1292 : * Returns true if general axisymmetric coordinate axes are being used
1293 : */
1294 : bool usingGeneralAxisymmetricCoordAxes() const;
1295 :
1296 : /**
1297 : * Returns the desired radial direction for RZ coordinate transformation
1298 : * @return The coordinate direction for the radial direction
1299 : */
1300 : unsigned int getAxisymmetricRadialCoord() const;
1301 :
1302 : /**
1303 : * Performs a sanity check for every element in the mesh. If an element dimension is 3 and the
1304 : * corresponding coordinate system is RZ, then this will error. If an element dimension is greater
1305 : * than 1 and the corresponding system is RPSHERICAL then this will error
1306 : */
1307 : void checkCoordinateSystems();
1308 :
1309 : /**
1310 : * Set the coordinate system data to that of \p other_mesh
1311 : */
1312 : void setCoordData(const MooseMesh & other_mesh);
1313 :
1314 : /**
1315 : * Mark the finite volume information as dirty
1316 : */
1317 3973 : void markFiniteVolumeInfoDirty() { _finite_volume_info_dirty = true; }
1318 :
1319 : /**
1320 : * @return whether the finite volume information is dirty
1321 : */
1322 1256 : bool isFiniteVolumeInfoDirty() const { return _finite_volume_info_dirty; }
1323 :
1324 : /**
1325 : * @return the coordinate transformation object that describes how to transform this problem's
1326 : * coordinate system into the canonical/reference coordinate system
1327 : */
1328 : MooseAppCoordTransform & coordTransform();
1329 :
1330 : /**
1331 : * @return the length unit of this mesh provided through the coordinate transformation object
1332 : */
1333 : const MooseUnits & lengthUnit() const;
1334 :
1335 : /**
1336 : * This function attempts to return the map from a high-order element side to its corresponding
1337 : * lower-d element
1338 : */
1339 : const std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *> &
1340 : getLowerDElemMap() const;
1341 :
1342 : /**
1343 : * @return Whether or not this mesh comes from a split mesh
1344 : */
1345 189063 : bool isSplit() const { return _is_split; }
1346 :
1347 : /**
1348 : * Builds the face and elem info vectors that store meta-data needed for looping over and doing
1349 : * calculations based on mesh faces and elements in a finite volume setting. This should only
1350 : * be called when finite volume variables are used in the problem or when the face and elem info
1351 : * objects are necessary for functor-based evaluations.
1352 : */
1353 : void buildFiniteVolumeInfo() const;
1354 :
1355 : /**
1356 : * Sets up the additional data needed for finite volume computations.
1357 : * This involves building FaceInfo and ElemInfo objects, caching variable associations
1358 : * and elemental DoF indices for FV variables.
1359 : */
1360 : void setupFiniteVolumeMeshData() const;
1361 :
1362 : /**
1363 : * Indicate whether the kind of adaptivity we're doing is p-refinement
1364 : */
1365 221 : void doingPRefinement(bool doing_p_refinement) { _doing_p_refinement = doing_p_refinement; }
1366 :
1367 : /**
1368 : * Query whether we have p-refinement
1369 : */
1370 406992162 : [[nodiscard]] bool doingPRefinement() const { return _doing_p_refinement; }
1371 :
1372 : /**
1373 : * Returns the maximum p-refinement level of all elements
1374 : */
1375 56677 : unsigned int maxPLevel() const { return _max_p_level; }
1376 :
1377 : /**
1378 : * Returns the maximum h-refinement level of all elements
1379 : */
1380 63028 : unsigned int maxHLevel() const { return _max_h_level; }
1381 :
1382 : /**
1383 : * Get the map describing for each volumetric quadrature point (qp) on the refined level which qp
1384 : * on the previous coarser level the fine qp is closest to
1385 : */
1386 : const std::vector<QpMap> & getPRefinementMap(const Elem & elem) const;
1387 : /**
1388 : * Get the map describing for each side quadrature point (qp) on the refined level which qp
1389 : * on the previous coarser level the fine qp is closest to
1390 : */
1391 : const std::vector<QpMap> & getPRefinementSideMap(const Elem & elem) const;
1392 : /**
1393 : * Get the map describing for each volumetric quadrature point (qp) on the coarse level which qp
1394 : * on the previous finer level the coarse qp is closest to
1395 : */
1396 : const std::vector<QpMap> & getPCoarseningMap(const Elem & elem) const;
1397 : /**
1398 : * Get the map describing for each side quadrature point (qp) on the coarse level which qp
1399 : * on the previous finer level the coarse qp is closest to
1400 : */
1401 : const std::vector<QpMap> & getPCoarseningSideMap(const Elem & elem) const;
1402 :
1403 : void buildPRefinementAndCoarseningMaps(Assembly * assembly);
1404 :
1405 : /**
1406 : * @return Whether the subdomain indicated by \p subdomain_id is a lower-dimensional manifold of
1407 : * some higher-dimensional subdomain, or in implementation speak, whether the elements of this
1408 : * subdomain have non-null interior parents
1409 : */
1410 : bool isLowerD(const SubdomainID subdomain_id) const;
1411 :
1412 : /**
1413 : * @return Whether there are any lower-dimensional blocks that are manifolds of higher-dimensional
1414 : * block faces
1415 : */
1416 34538 : bool hasLowerD() const { return _has_lower_d; }
1417 :
1418 : /**
1419 : * @return The set of lower-dimensional blocks for interior sides
1420 : */
1421 437231813 : const std::set<SubdomainID> & interiorLowerDBlocks() const { return _lower_d_interior_blocks; }
1422 : /**
1423 : * @return The set of lower-dimensional blocks for boundary sides
1424 : */
1425 436024267 : const std::set<SubdomainID> & boundaryLowerDBlocks() const { return _lower_d_boundary_blocks; }
1426 : /// Return construct node list from side list boolean
1427 141 : bool getConstructNodeListFromSideList() { return _construct_node_list_from_side_list; }
1428 :
1429 : protected:
1430 : /// Deprecated (DO NOT USE)
1431 : std::vector<std::unique_ptr<libMesh::GhostingFunctor>> _ghosting_functors;
1432 :
1433 : /// The list of active geometric relationship managers (bound to the underlying MeshBase object).
1434 : std::vector<std::shared_ptr<RelationshipManager>> _relationship_managers;
1435 :
1436 : /// Whether or not this mesh was built from another mesh
1437 : bool _built_from_other_mesh = false;
1438 :
1439 : /// Can be set to DISTRIBUTED, REPLICATED, or DEFAULT. Determines whether
1440 : /// the underlying libMesh mesh is a ReplicatedMesh or DistributedMesh.
1441 : ParallelType _parallel_type;
1442 :
1443 : /// False by default. Final value is determined by several factors
1444 : /// including the 'distribution' setting in the input file, and whether
1445 : /// or not the Mesh file is a Nemesis file.
1446 : bool _use_distributed_mesh;
1447 : bool _distribution_overridden;
1448 : bool _parallel_type_overridden;
1449 :
1450 : /// Pointer to underlying libMesh mesh object
1451 : std::unique_ptr<libMesh::MeshBase> _mesh;
1452 :
1453 : /// Pointer to Kokkos mesh object
1454 : #ifdef MOOSE_KOKKOS_ENABLED
1455 : std::unique_ptr<Moose::Kokkos::Mesh> _kokkos_mesh;
1456 : #endif
1457 :
1458 : /// The partitioner used on this mesh
1459 : MooseEnum _partitioner_name;
1460 : bool _partitioner_overridden;
1461 :
1462 : /// The custom partitioner
1463 : std::unique_ptr<libMesh::Partitioner> _custom_partitioner;
1464 : bool _custom_partitioner_requested;
1465 :
1466 : /// Convenience enums
1467 : enum
1468 : {
1469 : X = 0,
1470 : Y,
1471 : Z
1472 : };
1473 : enum
1474 : {
1475 : MIN = 0,
1476 : MAX
1477 : };
1478 :
1479 : /// The level of uniform refinement requested (set to zero if AMR is disabled)
1480 : unsigned int _uniform_refine_level;
1481 :
1482 : /// Whether or not to skip uniform refinements when using a pre-split mesh
1483 : bool _skip_refine_when_use_split;
1484 :
1485 : /// Whether or not skip remote deletion and repartition after uniform refinements
1486 : bool _skip_deletion_repartition_after_refine;
1487 :
1488 : /// true if mesh is changed (i.e. after adaptivity step)
1489 : bool _is_changed;
1490 :
1491 : /// True if a Nemesis Mesh was read in
1492 : bool _is_nemesis;
1493 :
1494 : /// True if prepare has been called on the mesh
1495 : bool _moose_mesh_prepared = false;
1496 :
1497 : /// The elements that were just refined.
1498 : std::unique_ptr<ConstElemPointerRange> _refined_elements;
1499 :
1500 : /// The elements that were just coarsened.
1501 : std::unique_ptr<ConstElemPointerRange> _coarsened_elements;
1502 :
1503 : /**
1504 : * Map of Parent elements to child elements for elements that were just coarsened.
1505 : *
1506 : * NOTE: the child element pointers ARE PROBABLY INVALID. Only use them for indexing!
1507 : */
1508 : std::map<const Elem *, std::vector<const Elem *>> _coarsened_element_children;
1509 :
1510 : /// Used for generating the semilocal node range
1511 : std::set<Node *> _semilocal_node_list;
1512 :
1513 : /**
1514 : * A range for use with threading. We do this so that it doesn't have
1515 : * to get rebuilt all the time (which takes time).
1516 : */
1517 : std::unique_ptr<libMesh::ConstElemRange> _active_local_elem_range;
1518 :
1519 : std::unique_ptr<SemiLocalNodeRange> _active_semilocal_node_range;
1520 : std::unique_ptr<libMesh::NodeRange> _active_node_range;
1521 : std::unique_ptr<libMesh::ConstNodeRange> _local_node_range;
1522 : std::unique_ptr<libMesh::StoredRange<MooseMesh::const_bnd_node_iterator, const BndNode *>>
1523 : _bnd_node_range;
1524 : std::unique_ptr<libMesh::StoredRange<MooseMesh::const_bnd_elem_iterator, const BndElement *>>
1525 : _bnd_elem_range;
1526 :
1527 : /// A map of all of the current nodes to the elements that they are connected to.
1528 : std::map<dof_id_type, std::vector<dof_id_type>> _node_to_elem_map;
1529 : bool _node_to_elem_map_built;
1530 :
1531 : /// A map of all of the current nodes to the active elements that they are connected to.
1532 : std::map<dof_id_type, std::vector<dof_id_type>> _node_to_active_semilocal_elem_map;
1533 : bool _node_to_active_semilocal_elem_map_built;
1534 :
1535 : /**
1536 : * A set of subdomain IDs currently present in the mesh. For parallel meshes, includes subdomains
1537 : * defined on other processors as well.
1538 : */
1539 : std::set<SubdomainID> _mesh_subdomains;
1540 :
1541 : ///@{
1542 : /**
1543 : * A set of boundary IDs currently present in the mesh. In serial, this is equivalent to the
1544 : * values returned by _mesh.get_boundary_info().get_boundary_ids(). In parallel, it will contain
1545 : * off-processor boundary IDs as well.
1546 : */
1547 : std::set<BoundaryID> _mesh_boundary_ids;
1548 : std::set<BoundaryID> _mesh_sideset_ids;
1549 : std::set<BoundaryID> _mesh_nodeset_ids;
1550 : ///@}
1551 :
1552 : /// The boundary to normal map - valid only when AddAllSideSetsByNormals is active
1553 : std::unique_ptr<std::map<BoundaryID, RealVectorValue>> _boundary_to_normal_map;
1554 :
1555 : /// array of boundary nodes
1556 : std::vector<BndNode *> _bnd_nodes;
1557 : typedef std::vector<BndNode *>::iterator bnd_node_iterator_imp;
1558 : typedef std::vector<BndNode *>::const_iterator const_bnd_node_iterator_imp;
1559 : /// Map of sets of node IDs in each boundary
1560 : std::map<boundary_id_type, std::set<dof_id_type>> _bnd_node_ids;
1561 :
1562 : /// array of boundary elems
1563 : std::vector<BndElement *> _bnd_elems;
1564 : typedef std::vector<BndElement *>::iterator bnd_elem_iterator_imp;
1565 : typedef std::vector<BndElement *>::const_iterator const_bnd_elem_iterator_imp;
1566 :
1567 : /// Map of set of elem IDs connected to each boundary
1568 : std::unordered_map<boundary_id_type, std::unordered_set<dof_id_type>> _bnd_elem_ids;
1569 :
1570 : std::map<dof_id_type, Node *> _quadrature_nodes;
1571 : std::map<dof_id_type, std::map<unsigned int, std::map<dof_id_type, Node *>>>
1572 : _elem_to_side_to_qp_to_quadrature_nodes;
1573 : std::vector<BndNode> _extra_bnd_nodes;
1574 :
1575 : /// list of nodes that belongs to a specified block (domain)
1576 : std::map<dof_id_type, std::set<SubdomainID>> _block_node_list;
1577 :
1578 : /// list of nodes that belongs to a specified nodeset: indexing [nodeset_id] -> [array of node ids]
1579 : std::map<boundary_id_type, std::vector<dof_id_type>> _node_set_nodes;
1580 :
1581 : std::set<unsigned int> _ghosted_boundaries;
1582 : std::vector<Real> _ghosted_boundaries_inflation;
1583 :
1584 : /// The number of nodes to consider in the NearestNode neighborhood.
1585 : unsigned int _patch_size;
1586 :
1587 : /// The number of nearest neighbors to consider for ghosting purposes when iteration patch update strategy is used.
1588 : unsigned int _ghosting_patch_size;
1589 :
1590 : // The maximum number of points in each leaf of the KDTree used in the nearest neighbor search.
1591 : unsigned int _max_leaf_size;
1592 :
1593 : /// The patch update strategy
1594 : Moose::PatchUpdateType _patch_update_strategy;
1595 :
1596 : /// Vector of all the Nodes in the mesh for determining when to add a new point
1597 : std::vector<Node *> _node_map;
1598 :
1599 : /// Boolean indicating whether this mesh was detected to be regular and orthogonal
1600 : bool _regular_orthogonal_mesh;
1601 :
1602 : /// The bounds in each dimension of the mesh for regular orthogonal meshes
1603 : std::vector<std::vector<Real>> _bounds;
1604 :
1605 : /// A vector holding the paired boundaries for a regular orthogonal mesh
1606 : std::vector<std::pair<BoundaryID, BoundaryID>> _paired_boundary;
1607 :
1608 : /// Whether or not we are using a (pre-)split mesh (automatically DistributedMesh)
1609 : const bool _is_split;
1610 :
1611 : void cacheInfo();
1612 : void freeBndNodes();
1613 : void freeBndElems();
1614 : void setPartitionerHelper(MeshBase * mesh = nullptr);
1615 :
1616 : private:
1617 : /// Map connecting elems with their corresponding ElemInfo, we use the element ID as
1618 : /// the key
1619 : mutable std::unordered_map<dof_id_type, ElemInfo> _elem_to_elem_info;
1620 :
1621 : /// Holds only those \p ElemInfo objects that have \p processor_id equal to this process's id,
1622 : /// e.g. the local \p ElemInfo objects
1623 : mutable std::vector<const ElemInfo *> _elem_info;
1624 :
1625 : /// FaceInfo object storing information for face based loops. This container holds all the \p
1626 : /// FaceInfo objects accessible from this process
1627 : mutable std::vector<FaceInfo> _all_face_info;
1628 :
1629 : /// Holds only those \p FaceInfo objects that have \p processor_id equal to this process's id,
1630 : /// e.g. the local \p FaceInfo objects
1631 : mutable std::vector<const FaceInfo *> _face_info;
1632 :
1633 : /// Map from elem-side pair to FaceInfo
1634 : mutable std::unordered_map<std::pair<const Elem *, unsigned int>, FaceInfo *>
1635 : _elem_side_to_face_info;
1636 :
1637 : // true if the _face_info member needs to be rebuilt/updated.
1638 : mutable bool _finite_volume_info_dirty = true;
1639 :
1640 : // True if we have cached elemental dofs ids for the linear finite volume variables.
1641 : // This happens in the first system which has a linear finite volume variable, considering
1642 : // that currently we only support one variable per linear system.
1643 : mutable bool _linear_finite_volume_dofs_cached = false;
1644 :
1645 : /**
1646 : * A map of vectors indicating which dimensions are periodic in a regular orthogonal mesh for
1647 : * the specified variable numbers. This data structure is populated by addPeriodicVariable.
1648 : */
1649 : std::map<unsigned int, std::vector<bool>> _periodic_dim;
1650 :
1651 : /**
1652 : * A convenience vector used to hold values in each dimension representing half of the range.
1653 : */
1654 : RealVectorValue _half_range;
1655 :
1656 : /// A vector containing the nodes at the corners of a regular orthogonal mesh
1657 : std::vector<Node *> _extreme_nodes;
1658 :
1659 : /**
1660 : * This routine detects paired sidesets of a regular orthogonal mesh (.i.e. parallel sidesets
1661 : * "across" from one and other).
1662 : * The _paired_boundary datastructure is populated with this information.
1663 : */
1664 : void detectPairedSidesets();
1665 :
1666 : /**
1667 : * Build the refinement map for a given element type. This will tell you what quadrature points
1668 : * to copy from and to for stateful material properties on newly created elements from Adaptivity.
1669 : *
1670 : * @param elem The element that represents the element type you need the refinement map for.
1671 : * @param qrule The quadrature rule in use.
1672 : * @param qrule_face The current face quadrature rule
1673 : * @param parent_side The side of the parent to map (-1 if not mapping parent sides)
1674 : * @param child The child number (-1 if not mapping child internal sides)
1675 : * @param child_side The side number of the child (-1 if not mapping sides)
1676 : */
1677 : void buildRefinementMap(const Elem & elem,
1678 : libMesh::QBase & qrule,
1679 : libMesh::QBase & qrule_face,
1680 : int parent_side,
1681 : int child,
1682 : int child_side);
1683 :
1684 : /**
1685 : * Build the coarsening map for a given element type. This will tell you what quadrature points
1686 : * to copy from and to for stateful material properties on newly created elements from Adaptivity.
1687 : *
1688 : * @param elem The element that represents the element type you need the coarsening map for.
1689 : * @param qrule The quadrature rule in use.
1690 : * @param qrule_face The current face quadrature rule
1691 : * @param input_side The side to map
1692 : */
1693 : void buildCoarseningMap(const Elem & elem,
1694 : libMesh::QBase & qrule,
1695 : libMesh::QBase & qrule_face,
1696 : int input_side);
1697 :
1698 : /**
1699 : * Find the closest points that map "from" to "to" and fill up "qp_map".
1700 : * Essentially, for each point in "from" find the closest point in "to".
1701 : *
1702 : * @param from The reference positions in the parent of the the points we're mapping _from_
1703 : * @param to The reference positions in the parent of the the points we're mapping _to_
1704 : * @param qp_map This will be filled with QpMap objects holding the mappings.
1705 : */
1706 : void mapPoints(const std::vector<Point> & from,
1707 : const std::vector<Point> & to,
1708 : std::vector<QpMap> & qp_map);
1709 :
1710 : /**
1711 : * Given an elem type, get maps that tell us what qp's are closest to each other between a parent
1712 : * and it's children.
1713 : * This is mainly used for mapping stateful material properties during adaptivity.
1714 : *
1715 : * There are 3 cases here:
1716 : *
1717 : * 1. Volume to volume (parent_side = -1, child = -1, child_side = -1)
1718 : * 2. Parent side to child side (parent_side = 0+, child = -1, child_side = 0+)
1719 : * 3. Child side to parent volume (parent_side = -1, child = 0+, child_side = 0+)
1720 : *
1721 : * Case 3 only happens under refinement (need to invent data at internal child sides).
1722 : *
1723 : * @param template_elem An element of the type that we need to find the maps for
1724 : * @param qrule The quadrature rule that we need to find the maps for
1725 : * @param qrule_face The face quadrature rule that we need to find the maps for
1726 : * @param refinement_map The map to use when an element gets split
1727 : * @param coarsen_map The map to use when an element is coarsened.
1728 : * @param parent_side - the id of the parent's side
1729 : * @param child - the id of the child element
1730 : * @param child_side - The id of the child's side
1731 : */
1732 : void findAdaptivityQpMaps(const Elem * template_elem,
1733 : libMesh::QBase & qrule,
1734 : libMesh::QBase & qrule_face,
1735 : std::vector<std::vector<QpMap>> & refinement_map,
1736 : std::vector<std::pair<unsigned int, QpMap>> & coarsen_map,
1737 : int parent_side,
1738 : int child,
1739 : int child_side);
1740 :
1741 : void buildHRefinementAndCoarseningMaps(Assembly * assembly);
1742 :
1743 : const std::vector<QpMap> & getPRefinementMapHelper(
1744 : const Elem & elem,
1745 : const std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>> &) const;
1746 : const std::vector<QpMap> & getPCoarseningMapHelper(
1747 : const Elem & elem,
1748 : const std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>> &) const;
1749 :
1750 : /**
1751 : * Update the coordinate transformation object based on our coordinate system data. The coordinate
1752 : * transformation will be created if it hasn't been already
1753 : */
1754 : void updateCoordTransform();
1755 :
1756 : /**
1757 : * Loop through all subdomain IDs and check if there is name duplication used for the subdomains
1758 : * with same ID. Throw out an error if any name duplication is found.
1759 : */
1760 : void checkDuplicateSubdomainNames();
1761 :
1762 : /// Holds mappings for volume to volume and parent side to child side
1763 : /// Map key:
1764 : /// - first member corresponds to element side. It's -1 for volume quadrature points
1765 : /// - second member correponds to the element type
1766 : /// Map value:
1767 : /// - Outermost index is the child element index
1768 : /// - Once we have indexed by the child element index, we have a std::vector of QpMaps. This
1769 : /// vector is sized by the number of reference points in the child element. Then for each
1770 : /// reference point in the child element we have a QpMap whose \p _from index corresponds to
1771 : /// the child element reference point, a \p _to index which corresponds to the reference point
1772 : /// on the parent element that the child element reference point is closest to, and a
1773 : /// \p _distance member which is the distance between the mapped child and parent reference
1774 : /// quadrature points
1775 : std::map<std::pair<int, libMesh::ElemType>, std::vector<std::vector<QpMap>>>
1776 : _elem_type_to_refinement_map;
1777 :
1778 : std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1779 : _elem_type_to_p_refinement_map;
1780 : std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1781 : _elem_type_to_p_refinement_side_map;
1782 :
1783 : /// Holds mappings for "internal" child sides to parent volume. The second key is (child, child_side).
1784 : std::map<libMesh::ElemType, std::map<std::pair<int, int>, std::vector<std::vector<QpMap>>>>
1785 : _elem_type_to_child_side_refinement_map;
1786 :
1787 : /// Holds mappings for volume to volume and parent side to child side
1788 : /// Map key:
1789 : /// - first member corresponds to element side. It's -1 for volume quadrature points
1790 : /// - second member correponds to the element type
1791 : /// Map value:
1792 : /// - Vector is sized based on the number of quadrature points in the parent (e.g. coarser)
1793 : /// element.
1794 : /// - For each parent quadrature point we store a pair
1795 : /// - The first member of the pair identifies which child holds the closest refined-level
1796 : /// quadrature point
1797 : /// - The second member of the pair is the QpMap. The \p _from data member will correspond to
1798 : /// the parent quadrature point index. The \p _to data member will correspond to which child
1799 : /// element quadrature point is closest to the parent quadrature point. And \p _distance is
1800 : /// the distance between the two
1801 : std::map<std::pair<int, libMesh::ElemType>, std::vector<std::pair<unsigned int, QpMap>>>
1802 : _elem_type_to_coarsening_map;
1803 :
1804 : std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1805 : _elem_type_to_p_coarsening_map;
1806 : std::map<std::pair<libMesh::ElemType, unsigned int>, std::vector<QpMap>>
1807 : _elem_type_to_p_coarsening_side_map;
1808 :
1809 : struct SubdomainData
1810 : {
1811 : /// Neighboring subdomain ids
1812 : std::set<SubdomainID> neighbor_subs;
1813 :
1814 : /// The boundary ids that are attached. This set will include any sideset boundary ID that
1815 : /// is a side of any part of the subdomain
1816 : std::set<BoundaryID> boundary_ids;
1817 :
1818 : /// Whether this subdomain is a lower-dimensional manifold of a higher-dimensional subdomain
1819 : bool is_lower_d;
1820 : };
1821 :
1822 : /// Holds a map from subdomain ids to associated data
1823 : std::unordered_map<SubdomainID, SubdomainData> _sub_to_data;
1824 :
1825 : /// Holds a map from neighbor subomdain ids to the boundary ids that are attached to it
1826 : std::unordered_map<SubdomainID, std::set<BoundaryID>> _neighbor_subdomain_boundary_ids;
1827 :
1828 : /// Mesh blocks for interior lower-d elements in different types
1829 : std::set<SubdomainID> _lower_d_interior_blocks;
1830 : /// Mesh blocks for boundary lower-d elements in different types
1831 : std::set<SubdomainID> _lower_d_boundary_blocks;
1832 : /// Holds a map from a high-order element side to its corresponding lower-d element
1833 : std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *>
1834 : _higher_d_elem_side_to_lower_d_elem;
1835 : std::unordered_map<const Elem *, unsigned short int> _lower_d_elem_to_higher_d_elem_side;
1836 :
1837 : /// Whether there are any lower-dimensional blocks that are manifolds of higher-dimensional block
1838 : /// faces
1839 : bool _has_lower_d;
1840 :
1841 : /// Whether or not this Mesh is allowed to read a recovery file
1842 : bool _allow_recovery;
1843 :
1844 : /// Whether or not to allow generation of nodesets from sidesets
1845 : bool _construct_node_list_from_side_list;
1846 :
1847 : /// Whether we need to delete remote elements after init'ing the EquationSystems
1848 : bool _need_delete;
1849 :
1850 : /// Whether to allow removal of remote elements
1851 : bool _allow_remote_element_removal;
1852 :
1853 : /// Set of elements ghosted by ghostGhostedBoundaries
1854 : std::set<Elem *> _ghost_elems_from_ghost_boundaries;
1855 :
1856 : /// A parallel mesh generator such as DistributedRectilinearMeshGenerator
1857 : /// already make everything ready. We do not need to gather all boundaries to
1858 : /// every single processor. In general, we should avoid using ghostGhostedBoundaries
1859 : /// when possible since it is not scalable
1860 : bool _need_ghost_ghosted_boundaries;
1861 :
1862 : /// Unique element integer IDs for each subdomain and each extra element integers
1863 : std::vector<std::unordered_map<SubdomainID, std::set<dof_id_type>>> _block_id_mapping;
1864 : /// Maximum integer ID for each extra element integer
1865 : std::vector<dof_id_type> _max_ids;
1866 : /// Minimum integer ID for each extra element integer
1867 : std::vector<dof_id_type> _min_ids;
1868 : /// Flags to indicate whether or not any two extra element integers are the same
1869 : std::vector<std::vector<bool>> _id_identical_flag;
1870 :
1871 : /// The maximum number of sides per element
1872 : unsigned int _max_sides_per_elem;
1873 :
1874 : /// The maximum number of nodes per element
1875 : unsigned int _max_nodes_per_elem;
1876 :
1877 : /// The maximum number of nodes per side
1878 : unsigned int _max_nodes_per_side;
1879 :
1880 : /// Compute the maximum numbers per element and side
1881 : void computeMaxPerElemAndSide();
1882 :
1883 : /// Whether this mesh is displaced
1884 : bool _is_displaced;
1885 :
1886 : /// Build extra data for faster access to the information of extra element integers
1887 : void buildElemIDInfo();
1888 :
1889 : /// Build lower-d mesh for all sides
1890 : void buildLowerDMesh();
1891 :
1892 : /// Type of coordinate system per subdomain
1893 : std::map<SubdomainID, Moose::CoordinateSystemType> & _coord_sys;
1894 :
1895 : /// Storage for RZ axis selection
1896 : unsigned int _rz_coord_axis;
1897 :
1898 : /// Map of subdomain ID to general axisymmetric axis
1899 : std::unordered_map<SubdomainID, std::pair<Point, RealVectorValue>> _subdomain_id_to_rz_coord_axis;
1900 :
1901 : /// A coordinate transformation object that describes how to transform this problem's coordinate
1902 : /// system into the canonical/reference coordinate system
1903 : std::unique_ptr<MooseAppCoordTransform> _coord_transform;
1904 :
1905 : /// Whether the coordinate system has been set
1906 : bool _coord_system_set;
1907 :
1908 : /// Set for holding user-provided coordinate system type block names
1909 : std::vector<SubdomainName> _provided_coord_blocks;
1910 :
1911 : /// Whether we have p-refinement (as opposed to h-refinement)
1912 : bool _doing_p_refinement;
1913 : /// Maximum p-refinement level of all elements
1914 : unsigned int _max_p_level;
1915 : /// Maximum h-refinement level of all elements
1916 : unsigned int _max_h_level;
1917 :
1918 : template <typename T>
1919 : struct MeshType;
1920 : };
1921 :
1922 : inline MooseAppCoordTransform &
1923 167259 : MooseMesh::coordTransform()
1924 : {
1925 : mooseAssert(_coord_transform, "The coordinate transformation object is null.");
1926 167259 : return *_coord_transform;
1927 : }
1928 :
1929 : template <>
1930 : struct MooseMesh::MeshType<libMesh::ReplicatedMesh>
1931 : {
1932 : static const ParallelType value = ParallelType::REPLICATED;
1933 : };
1934 :
1935 : template <>
1936 : struct MooseMesh::MeshType<libMesh::DistributedMesh>
1937 : {
1938 : static const ParallelType value = ParallelType::DISTRIBUTED;
1939 : };
1940 :
1941 : /**
1942 : * The definition of the face_info_iterator struct.
1943 : */
1944 : struct MooseMesh::face_info_iterator
1945 : : variant_filter_iterator<MeshBase::Predicate, const FaceInfo *>
1946 : {
1947 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
1948 : template <typename PredType, typename IterType>
1949 414090 : face_info_iterator(const IterType & d, const IterType & e, const PredType & p)
1950 414090 : : variant_filter_iterator<MeshBase::Predicate, const FaceInfo *>(d, e, p)
1951 : {
1952 414090 : }
1953 : };
1954 :
1955 : /**
1956 : * The definition of the const_face_info_iterator struct. It is similar to the
1957 : * iterator above, but also provides an additional conversion-to-const ctor.
1958 : */
1959 : struct MooseMesh::const_face_info_iterator : variant_filter_iterator<MeshBase::Predicate,
1960 : const FaceInfo * const,
1961 : const FaceInfo * const &,
1962 : const FaceInfo * const *>
1963 : {
1964 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
1965 : template <typename PredType, typename IterType>
1966 : const_face_info_iterator(const IterType & d, const IterType & e, const PredType & p)
1967 : : variant_filter_iterator<MeshBase::Predicate,
1968 : const FaceInfo * const,
1969 : const FaceInfo * const &,
1970 : const FaceInfo * const *>(d, e, p)
1971 : {
1972 : }
1973 :
1974 : // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
1975 : // variant_filter_iterator copy constructor. Note that this one is *not* templated!
1976 414090 : const_face_info_iterator(const MooseMesh::face_info_iterator & rhs)
1977 414090 : : variant_filter_iterator<MeshBase::Predicate,
1978 : const FaceInfo * const,
1979 : const FaceInfo * const &,
1980 414090 : const FaceInfo * const *>(rhs)
1981 : {
1982 414090 : }
1983 : };
1984 :
1985 : /**
1986 : * The definition of the elem_info_iterator struct.
1987 : */
1988 : struct MooseMesh::elem_info_iterator
1989 : : variant_filter_iterator<MeshBase::Predicate, const ElemInfo *>
1990 : {
1991 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
1992 : template <typename PredType, typename IterType>
1993 44656 : elem_info_iterator(const IterType & d, const IterType & e, const PredType & p)
1994 44656 : : variant_filter_iterator<MeshBase::Predicate, const ElemInfo *>(d, e, p)
1995 : {
1996 44656 : }
1997 : };
1998 :
1999 : /**
2000 : * The definition of the const_elem_info_iterator struct. It is similar to the
2001 : * iterator above, but also provides an additional conversion-to-const ctor.
2002 : */
2003 : struct MooseMesh::const_elem_info_iterator : variant_filter_iterator<MeshBase::Predicate,
2004 : const ElemInfo * const,
2005 : const ElemInfo * const &,
2006 : const ElemInfo * const *>
2007 : {
2008 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2009 : template <typename PredType, typename IterType>
2010 : const_elem_info_iterator(const IterType & d, const IterType & e, const PredType & p)
2011 : : variant_filter_iterator<MeshBase::Predicate,
2012 : const ElemInfo * const,
2013 : const ElemInfo * const &,
2014 : const ElemInfo * const *>(d, e, p)
2015 : {
2016 : }
2017 :
2018 : // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2019 : // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2020 44656 : const_elem_info_iterator(const MooseMesh::elem_info_iterator & rhs)
2021 44656 : : variant_filter_iterator<MeshBase::Predicate,
2022 : const ElemInfo * const,
2023 : const ElemInfo * const &,
2024 44656 : const ElemInfo * const *>(rhs)
2025 : {
2026 44656 : }
2027 : };
2028 :
2029 : /**
2030 : * The definition of the bnd_node_iterator struct.
2031 : */
2032 : struct MooseMesh::bnd_node_iterator : variant_filter_iterator<MeshBase::Predicate, BndNode *>
2033 : {
2034 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2035 : template <typename PredType, typename IterType>
2036 176568 : bnd_node_iterator(const IterType & d, const IterType & e, const PredType & p)
2037 176568 : : variant_filter_iterator<MeshBase::Predicate, BndNode *>(d, e, p)
2038 : {
2039 176568 : }
2040 : };
2041 :
2042 : /**
2043 : * The definition of the const_bnd_node_iterator struct. It is similar to the
2044 : * iterator above, but also provides an additional conversion-to-const ctor.
2045 : */
2046 : struct MooseMesh::const_bnd_node_iterator : variant_filter_iterator<MeshBase::Predicate,
2047 : BndNode * const,
2048 : BndNode * const &,
2049 : BndNode * const *>
2050 : {
2051 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2052 : template <typename PredType, typename IterType>
2053 4382 : const_bnd_node_iterator(const IterType & d, const IterType & e, const PredType & p)
2054 : : variant_filter_iterator<MeshBase::Predicate,
2055 : BndNode * const,
2056 : BndNode * const &,
2057 4382 : BndNode * const *>(d, e, p)
2058 : {
2059 4382 : }
2060 :
2061 : // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2062 : // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2063 176568 : const_bnd_node_iterator(const MooseMesh::bnd_node_iterator & rhs)
2064 176568 : : variant_filter_iterator<MeshBase::Predicate,
2065 : BndNode * const,
2066 : BndNode * const &,
2067 176568 : BndNode * const *>(rhs)
2068 : {
2069 176568 : }
2070 : };
2071 :
2072 : /**
2073 : * The definition of the bnd_elem_iterator struct.
2074 : */
2075 : struct MooseMesh::bnd_elem_iterator : variant_filter_iterator<MeshBase::Predicate, BndElement *>
2076 : {
2077 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2078 : template <typename PredType, typename IterType>
2079 176566 : bnd_elem_iterator(const IterType & d, const IterType & e, const PredType & p)
2080 176566 : : variant_filter_iterator<MeshBase::Predicate, BndElement *>(d, e, p)
2081 : {
2082 176566 : }
2083 : };
2084 :
2085 : /**
2086 : * The definition of the const_bnd_elem_iterator struct. It is similar to the regular
2087 : * iterator above, but also provides an additional conversion-to-const ctor.
2088 : */
2089 : struct MooseMesh::const_bnd_elem_iterator : variant_filter_iterator<MeshBase::Predicate,
2090 : BndElement * const,
2091 : BndElement * const &,
2092 : BndElement * const *>
2093 : {
2094 : // Templated forwarding ctor -- forwards to appropriate variant_filter_iterator ctor
2095 : template <typename PredType, typename IterType>
2096 : const_bnd_elem_iterator(const IterType & d, const IterType & e, const PredType & p)
2097 : : variant_filter_iterator<MeshBase::Predicate,
2098 : BndElement * const,
2099 : BndElement * const &,
2100 : BndElement * const *>(d, e, p)
2101 : {
2102 : }
2103 :
2104 : // The conversion-to-const ctor. Takes a regular iterator and calls the appropriate
2105 : // variant_filter_iterator copy constructor. Note that this one is *not* templated!
2106 176226 : const_bnd_elem_iterator(const bnd_elem_iterator & rhs)
2107 176226 : : variant_filter_iterator<MeshBase::Predicate,
2108 : BndElement * const,
2109 : BndElement * const &,
2110 176226 : BndElement * const *>(rhs)
2111 : {
2112 176226 : }
2113 : };
2114 :
2115 : /**
2116 : * Some useful StoredRange typedefs. These are defined *outside* the
2117 : * MooseMesh class to mimic the Const{Node,Elem}Range classes in libmesh.
2118 : */
2119 : typedef libMesh::StoredRange<MooseMesh::const_bnd_node_iterator, const BndNode *> ConstBndNodeRange;
2120 : typedef libMesh::StoredRange<MooseMesh::const_bnd_elem_iterator, const BndElement *>
2121 : ConstBndElemRange;
2122 :
2123 : template <typename T>
2124 : std::unique_ptr<T>
2125 72002 : MooseMesh::buildTypedMesh(unsigned int dim)
2126 : {
2127 : // If the requested mesh type to build doesn't match our current value for _use_distributed_mesh,
2128 : // then we need to make sure to make our state consistent because other objects, like the periodic
2129 : // boundary condition action, will be querying isDistributedMesh()
2130 72002 : if (_use_distributed_mesh != std::is_same<T, libMesh::DistributedMesh>::value)
2131 : {
2132 835 : if (getMeshPtr())
2133 0 : mooseError("A MooseMesh object is being asked to build a libMesh mesh that is a different "
2134 : "parallel type than the libMesh mesh that it wraps. This is not allowed. Please "
2135 : "create another MooseMesh object to wrap the new libMesh mesh");
2136 835 : setParallelType(MeshType<T>::value);
2137 : }
2138 :
2139 72002 : if (dim == libMesh::invalid_uint)
2140 : {
2141 152010 : if (isParamValid("dim"))
2142 124827 : dim = getParam<MooseEnum>("dim");
2143 : else
2144 : // Legacy selection of the default for the 'dim' parameter
2145 9061 : dim = 1;
2146 : }
2147 :
2148 72002 : auto mesh = std::make_unique<T>(_communicator, dim);
2149 :
2150 216006 : if (!getParam<bool>("allow_renumbering"))
2151 1720 : mesh->allow_renumbering(false);
2152 :
2153 72002 : mesh->allow_remote_element_removal(_allow_remote_element_removal);
2154 72002 : _app.attachRelationshipManagers(*mesh, *this);
2155 :
2156 72002 : if (_custom_partitioner_requested)
2157 : {
2158 : // Check of partitioner is supplied (not allowed if custom partitioner is used)
2159 5160 : if (!parameters().isParamSetByAddParam("partitioner"))
2160 0 : mooseError("If partitioner block is provided, partitioner keyword cannot be used!");
2161 : // Set custom partitioner
2162 1720 : if (!_custom_partitioner.get())
2163 0 : mooseError("Custom partitioner requested but not set!");
2164 1720 : mesh->partitioner() = _custom_partitioner->clone();
2165 : }
2166 : else
2167 70282 : setPartitionerHelper(mesh.get());
2168 :
2169 72002 : return mesh;
2170 0 : }
2171 :
2172 : inline bool
2173 5152 : MooseMesh::skipDeletionRepartitionAfterRefine() const
2174 : {
2175 5152 : return _skip_deletion_repartition_after_refine;
2176 : }
2177 :
2178 : inline void
2179 902 : MooseMesh::setParallelType(ParallelType parallel_type)
2180 : {
2181 902 : _parallel_type = parallel_type;
2182 902 : determineUseDistributedMesh();
2183 902 : }
2184 :
2185 : inline bool
2186 : MooseMesh::hasElementID(const std::string & id_name) const
2187 : {
2188 : return getMesh().has_elem_integer(id_name);
2189 : }
2190 :
2191 : inline unsigned int
2192 : MooseMesh::getElementIDIndex(const std::string & id_name) const
2193 : {
2194 : if (!hasElementID(id_name))
2195 : mooseError("Mesh does not have element ID for ", id_name);
2196 : return getMesh().get_elem_integer_index(id_name);
2197 : }
2198 :
2199 : inline bool
2200 : MooseMesh::areElemIDsIdentical(const std::string & id_name1, const std::string & id_name2) const
2201 : {
2202 : auto id1 = getElementIDIndex(id_name1);
2203 : auto id2 = getElementIDIndex(id_name2);
2204 : return _id_identical_flag[id1][id2];
2205 : }
2206 :
2207 : inline const std::vector<const FaceInfo *> &
2208 43 : MooseMesh::faceInfo() const
2209 : {
2210 43 : return _face_info;
2211 : }
2212 :
2213 : inline const std::vector<FaceInfo> &
2214 4 : MooseMesh::allFaceInfo() const
2215 : {
2216 4 : return _all_face_info;
2217 : }
2218 :
2219 : inline const std::map<boundary_id_type, std::vector<dof_id_type>> &
2220 : MooseMesh::nodeSetNodes() const
2221 : {
2222 : return _node_set_nodes;
2223 : }
2224 :
2225 : inline const std::unordered_map<std::pair<const Elem *, unsigned short int>, const Elem *> &
2226 : MooseMesh::getLowerDElemMap() const
2227 : {
2228 : return _higher_d_elem_side_to_lower_d_elem;
2229 : }
2230 :
2231 : inline bool
2232 506085 : MooseMesh::isLowerD(const SubdomainID subdomain_id) const
2233 : {
2234 506085 : return libmesh_map_find(_sub_to_data, subdomain_id).is_lower_d;
2235 : }
|