37 #ifndef domainDecomposition_H
38 #define domainDecomposition_H
62 const word regionName_;
118 void addInterProcFace
121 const label ownerProc,
123 const label subPatchID,
135 template<
class BinaryOp>
136 inline void processInterCyclics
151 void decomposePoints();
155 void reconstructPoints();
158 bool completeConformal()
const;
161 bool procsConformal()
const;
174 void validateComplete()
const;
177 void validateProcs()
const;
180 void readComplete(
const bool stitch =
true);
186 void readCompleteAddressing();
189 void readProcsAddressing();
192 void readAddressing();
206 void writeCompleteAddressing()
const;
209 void writeProcsAddressing()
const;
212 void writeAddressing()
const;
217 void writeProcPoints(
const fileName& inst);
222 void writeCompletePoints(
const fileName& inst);
251 return completeMesh_();
264 return runTimes_.
nProcs();
293 return procPointAddressing_;
300 return procFaceAddressing_;
307 return procCellAddressing_;
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
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.
TypeName("domainDecomposition")
Runtime type information.
label nProcs() const
Return the number of processors in the decomposition.
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
bool readDecompose(const bool doSets)
Read in the complete mesh. Read the processor meshes if they are.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
bool readReconstruct(const bool doSets)
Read in the processor meshes. Read the complete mesh if it is.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
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.
const fvMesh & completeMesh() const
Access the global mesh.
domainDecomposition(const processorRunTimes &procRunTimes, const word ®ionName)
Construct from processor run times and region name.
fvMesh::readUpdateState readUpdateReconstruct()
Read-update for reconstruction.
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.
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.