Go to the documentation of this file.
41 Foam::hexRef8Data::hexRef8Data(
const IOobject& io)
51 cellLevelPtr_.reset(
new Type(rio));
62 pointLevelPtr_.reset(
new Type(rio));
73 level0EdgePtr_.reset(
new Type(rio));
78 IOobject rio(io,
"refinementHistory");
84 refHistoryPtr_.reset(
new Type(rio));
90 Foam::hexRef8Data::hexRef8Data
98 if (
data.cellLevelPtr_)
111 if (
data.pointLevelPtr_)
124 if (
data.level0EdgePtr_)
133 if (
data.refHistoryPtr_)
137 refHistoryPtr_ =
data.refHistoryPtr_().
clone(rio, cellMap);
142 Foam::hexRef8Data::hexRef8Data
154 if (procDatas[0].cellLevelPtr_)
159 auto& cellLevel = *cellLevelPtr_;
163 const labelList& procCellLevel = procDatas[procI].cellLevelPtr_();
171 if (procDatas[0].pointLevelPtr_)
176 auto& pointLevel = *pointLevelPtr_;
180 const labelList& procPointLevel = procDatas[procI].pointLevelPtr_();
181 labelUIndList(pointLevel, pointMaps[procI]) = procPointLevel;
188 if (procDatas[0].level0EdgePtr_)
197 procDatas[0].level0EdgePtr_()
205 if (procDatas[0].refHistoryPtr_)
212 procRefs.set(i, &procDatas[i].refHistoryPtr_());
241 if (hasCellLevel && !cellLevelPtr_)
252 if (hasPointLevel && !pointLevelPtr_)
284 if (hasHistory && !refHistoryPtr_)
286 IOobject rio(io,
"refinementHistory");
298 (cellLevelPtr_ && cellLevelPtr_().size() != map.
nOldCells())
299 || (pointLevelPtr_ && pointLevelPtr_().size() != map.
nOldPoints())
302 cellLevelPtr_.clear();
303 pointLevelPtr_.clear();
304 level0EdgePtr_.clear();
305 refHistoryPtr_.clear();
318 label oldCelli = cellMap[newCelli];
322 newCellLevel[newCelli] = 0;
326 newCellLevel[newCelli] = cellLevel[oldCelli];
335 labelList& pointLevel = pointLevelPtr_();
337 labelList newPointLevel(pointMap.size());
348 newPointLevel[
newPointi] = pointLevel[oldPointi];
356 if (refHistoryPtr_ && refHistoryPtr_().active())
358 refHistoryPtr_().updateMesh(map);
377 if (refHistoryPtr_ && refHistoryPtr_().active())
379 refHistoryPtr_().distribute(map);
389 ok = ok && cellLevelPtr_().write();
393 ok = ok && pointLevelPtr_().write();
397 ok = ok && level0EdgePtr_().write();
401 ok = ok && refHistoryPtr_().write();
List< label > labelList
A List of labels.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
label size() const noexcept
The number of elements in the list.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static constexpr const zero Zero
Global zero (0)
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
const fileName & facesInstance() const
Return the current instance directory for faces.
static bool master(const label communicator=worldComm)
Am I the master process.
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const word & name() const
IOList< label > labelIOList
Label container classes.
Mesh consisting of general polyhedral cells.
label nPoints() const noexcept
Number of mesh points.
#define forAll(list, i)
Loop across all elements in list.
const mapDistribute & cellMap() const
Cell distribute map.
All refinement history. Used in unrefinement.
label nCells() const noexcept
Number of mesh cells.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label nOldPoints() const
Number of old points.
UniformDimensionedField< scalar > uniformDimensionedScalarField
void transfer(List< T > &list)
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void distribute(const mapDistributePolyMesh &)
In-place distribute.
Various for reading/decomposing/reconstructing/distributing refinement data.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
void sync(const IOobject &io)
Parallel synchronise. This enforces valid objects on all processors.
readOption readOpt() const noexcept
The read option.
const word & name() const noexcept
Return name.
const labelList & pointMap() const
Old point map.
label nOldCells() const
Number of old cells.
autoPtr< dictionary > clone() const
Construct and return clone.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
~hexRef8Data()
Destructor.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const mapDistribute & pointMap() const
Point distribute map.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const polyMesh & mesh() const
Return polyMesh.
Database for solution data, solver performance and other reduced data.
const labelList & cellMap() const
Old cell map.
UIndirectList< label > labelUIndList
UIndirectList of labels.
const objectRegistry & db() const noexcept
Return the local objectRegistry.