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