38 #ifndef domainDecomposition_H
39 #define domainDecomposition_H
53 class multiDomainDecomposition;
76 static bool sortReconstructNonConformalCyclicAddressing_;
85 const word regionName_;
142 void addInterProcFace
145 const label ownerProc,
147 const label subPatchID,
159 void processInterCyclics
177 void decomposePoints();
186 const label masterMeshProcStart,
187 const label masterMeshProcEnd,
189 const label meshToAddProcStart,
190 const label meshToAddProcEnd,
200 void reconstructPoints();
206 bool completeConformal()
const;
209 bool procsConformal()
const;
214 labelList completeFaceAddressing()
const;
225 nonConformalMappedWallProcOffsets(
const bool appendSize)
const;
229 void decomposeNonConformalCyclicAddressing
231 const label nccPatchi,
238 void decomposeNonConformalMappedWallAddressing
240 const label ncmwPatchi,
247 void decomposeNonConformalErrorAddressing
249 const label ncePatchi,
255 void reconstructNonConformalCyclicAddressing
257 const label nccPatchi,
272 void sortReconstructNonConformalCyclicAddressing
274 const label nccPatchi,
281 void reconstructNonConformalMappedWallAddressing
283 const label ncmwPatchi,
289 void reconstructNonConformalErrorAddressing
291 const label ncePatchi,
298 nonConformalProcFaceAddressingBf()
const;
302 void unconformComplete();
306 void unconformProcs();
319 void validateComplete()
const;
322 void validateProcs()
const;
328 void readProcs(
const bool doPost);
331 void readCompleteAddressing();
334 void readProcsAddressing();
337 void readAddressing();
343 void writeCompleteAddressing()
const;
346 void writeProcsAddressing()
const;
349 void writeAddressing()
const;
354 void writeProcPoints(
const fileName& inst);
359 void writeCompletePoints(
const fileName& inst);
362 template<
class SetType>
381 NullObjectRef<multiDomainDecomposition>()
412 return completeMesh_.
valid();
418 return procMeshes_.
size() && procMeshes_.
set(0);
425 return completeMesh_();
438 return runTimes_.
nProcs();
503 return procPointAddressing_;
510 return procFaceAddressing_;
517 return procCellAddressing_;
531 template<
class SetType>
535 template<
class SetType>
551 label domainDecomposition::setIndex<cellSet>(
const label,
const label)
const;
554 label domainDecomposition::setIndex<faceSet>(
const label,
const label)
const;
557 label domainDecomposition::setIndex<pointSet>(
const label,
const label)
const;
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
An STL-conforming hash table.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
bool set(const label) const
Is element set.
label size() const
Return the number of elements in the UPtrList.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Automatic domain decomposition class for finite-volume meshes.
const labelListList & procPointAddressing() const
Access the labels of points for each processor.
const fileName & meshPath() const
Access the mesh path.
bool readDecompose(const bool doPost)
Read in the complete mesh. Read the processor meshes if they are.
fvMesh::readUpdateState readUpdateDecompose(const bool doPost)
Read-update for decomposition.
fvMesh::readUpdateState readUpdateReconstruct(const bool doPost)
Read-update for reconstruction.
void postReadReconstruct()
Post-read-construction steps for the meshes after read-reconstruct.
const word & regionName() const
Access the region name.
void postReadDecompose()
Post-read-construction steps for the meshes after read-decompose.
PtrList< SetType > decomposeSet(const word &name) const
Decompose a set.
void unconformReadReconstruct()
Synchronise non-conformal structure after read-reconstruct.
void readComplete()
Read the complete mesh only.
TypeName("domainDecomposition")
Runtime type information.
autoPtr< SetType > reconstructSet(const word &name) const
Reconstruct a set.
void unconformReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-reconstruct.
void postReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Complete read-update-reconstruct.
fvMesh::readUpdateState readUpdateComplete()
Read-update the complete mesh only.
label nProcs() const
Return the number of processors in the decomposition.
const processorRunTimes & procRunTimes() const
Access the run times.
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
void writeReadReconstruct(const bool reconstructed, const bool doSets)
Write following read-reconstruct.
bool haveComplete() const
Do we have the global mesh?
void writeReadDecompose(const bool decomposed, const bool doSets)
Write following read-decompose.
bool readReconstruct(const bool doPost)
Read in the processor meshes. Read the complete mesh if it is.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
void postReadUpdateDecompose(const fvMesh::readUpdateState stat)
Complete read-update-decompose.
const PtrList< surfaceLabelField::Boundary > & procFaceAddressingBf() const
Access the labels of faces for each processor (see notes above)
const labelList & cellProc() const
Return the distribution as an FV field for writing.
const labelListList & procCellAddressing() const
Access the labels of cells for each processor.
domainDecomposition(const processorRunTimes &procRunTimes, const fileName &meshPath, const word ®ionName, const multiDomainDecomposition ®ionMeshes=NullObjectRef< multiDomainDecomposition >())
Construct from processor run times and region name.
void unconformReadDecompose()
Synchronise non-conformal structure after read-decompose.
const fvMesh & completeMesh() const
Access the global mesh.
void unconformReadUpdateDecompose(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-decompose.
virtual ~domainDecomposition()
Destructor.
bool haveProcs() const
Do we have the global mesh?
A class for handling file names.
Mesh data needed to do the Finite Volume discretisation.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Mesh consisting of general polyhedral cells.
label nProcs() const
Return the number of processors.
A class for handling words, derived from string.
const fvPatchList & patches
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & regionName(const solver ®ion)
labelList first(const UList< labelPair > &p)
const word & regionName< domainDecomposition >(const domainDecomposition ®ion)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.