38 #ifndef domainDecomposition_H
39 #define domainDecomposition_H
53 class multiDomainDecomposition;
73 static bool sortReconstructNonConformalCyclicAddressing_;
79 const word regionName_;
136 void addInterProcFace
139 const label ownerProc,
141 const label subPatchID,
153 void processInterCyclics
171 void decomposePoints();
180 const label masterMeshProcStart,
181 const label masterMeshProcEnd,
183 const label meshToAddProcStart,
184 const label meshToAddProcEnd,
194 void reconstructPoints();
200 bool completeConformal()
const;
203 bool procsConformal()
const;
208 labelList completeFaceAddressing()
const;
219 nonConformalMappedWallProcOffsets(
const bool appendSize)
const;
223 void decomposeNonConformalCyclicAddressing
225 const label nccPatchi,
232 void decomposeNonConformalMappedWallAddressing
234 const label ncmwPatchi,
241 void decomposeNonConformalErrorAddressing
243 const label ncePatchi,
249 void reconstructNonConformalCyclicAddressing
251 const label nccPatchi,
266 void sortReconstructNonConformalCyclicAddressing
268 const label nccPatchi,
275 void reconstructNonConformalMappedWallAddressing
277 const label ncmwPatchi,
283 void reconstructNonConformalErrorAddressing
285 const label ncePatchi,
292 nonConformalProcFaceAddressingBf()
const;
296 void unconformComplete();
300 void unconformProcs();
313 void validateComplete()
const;
316 void validateProcs()
const;
325 void readCompleteAddressing();
328 void readProcsAddressing();
331 void readAddressing();
337 void writeCompleteAddressing()
const;
340 void writeProcsAddressing()
const;
343 void writeAddressing()
const;
348 void writeProcPoints(
const fileName& inst);
353 void writeCompletePoints(
const fileName& inst);
370 NullObjectRef<multiDomainDecomposition>()
396 return completeMesh_();
409 return runTimes_.
nProcs();
468 return procPointAddressing_;
475 return procFaceAddressing_;
482 return procCellAddressing_;
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...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Automatic domain decomposition class for finite-volume meshes.
const labelListList & procPointAddressing() const
Access the labels of points for each processor.
bool readReconstruct()
Read in the processor meshes. Read the complete mesh if it is.
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.
bool readDecompose()
Read in the complete mesh. Read the processor meshes if they are.
void unconformReadReconstruct()
Synchronise non-conformal structure after read-reconstruct.
TypeName("domainDecomposition")
Runtime type information.
void unconformReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-reconstruct.
domainDecomposition(const processorRunTimes &procRunTimes, const word ®ionName, const multiDomainDecomposition ®ionMeshes=NullObjectRef< multiDomainDecomposition >())
Construct from processor run times and region name.
void postReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Complete read-update-reconstruct.
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.
void writeReadDecompose(const bool decomposed, const bool doSets)
Write following read-decompose.
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.
void unconformReadDecompose()
Synchronise non-conformal structure after read-decompose.
const fvMesh & completeMesh() const
Access the global mesh.
fvMesh::readUpdateState readUpdateReconstruct()
Read-update for reconstruction.
void unconformReadUpdateDecompose(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-decompose.
fvMesh::readUpdateState readUpdateDecompose()
Read-update for decomposition.
virtual ~domainDecomposition()
Destructor.
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)