55void Foam::domainDecomposition::mark
57 const labelList& zoneElems,
59 labelList& elementToZone
62 for (
const label pointi : zoneElems)
64 if (elementToZone[pointi] == -1)
67 elementToZone[pointi] = zoneI;
69 else if (elementToZone[pointi] >= 0)
72 elementToZone[pointi] = -2;
87 facesInstancePointsPtr_
89 pointsInstance() != facesInstance()
105 decompDictFile_(decompDictFile),
118 cellToProc_(nCells()),
119 procPointAddressing_(nProcs_),
120 procFaceAddressing_(nProcs_),
121 procCellAddressing_(nProcs_),
122 procPatchSize_(nProcs_),
123 procPatchStartIndex_(nProcs_),
124 procNeighbourProcessors_(nProcs_),
125 procProcessorPatchSize_(nProcs_),
126 procProcessorPatchStartIndex_(nProcs_),
127 procProcessorPatchSubPatchIDs_(nProcs_),
128 procProcessorPatchSubPatchStarts_(nProcs_)
130 updateParameters(this->model());
147 params.readIfPresent(
"distributed", distributed_);
153 Info<<
"\nConstructing processor meshes" <<
endl;
166 forAll(pointZones(), zonei)
168 mark(pointZones()[zonei], zonei, pointToZone);
172 labelList faceToZone(faces().size(), -1);
174 forAll(faceZones(), zonei)
176 mark(faceZones()[zonei], zonei, faceToZone);
182 forAll(cellZones(), zonei)
184 mark(cellZones()[zonei], zonei, cellToZone);
194 IOobjectList objects(*
this, facesInstance(),
"polyMesh/sets");
197 cellSets.emplace_back(
io);
201 faceSets.emplace_back(
io);
205 pointSets.emplace_back(
io);
226 label maxProcCells = 0;
227 label maxProcFaces = 0;
228 label totProcFaces = 0;
229 label maxProcPatches = 0;
230 label totProcPatches = 0;
233 for (label proci = 0; proci < nProcs_; proci++)
236 const labelList& curPointLabels = procPointAddressing_[proci];
244 forAll(curPointLabels, pointi)
246 procPoints[pointi] = meshPoints[curPointLabels[pointi]];
248 pointLookup[curPointLabels[pointi]] = pointi;
252 const labelList& curFaceLabels = procFaceAddressing_[proci];
254 const faceList& meshFaces = faces();
258 faceList procFaces(curFaceLabels.size());
260 forAll(curFaceLabels, facei)
264 label curF =
mag(curFaceLabels[facei]) - 1;
266 faceLookup[curF] = facei;
271 if (curFaceLabels[facei] >= 0)
274 origFaceLabels = meshFaces[curF];
278 origFaceLabels = meshFaces[curF].reverseFace();
282 face& procFaceLabels = procFaces[facei];
284 procFaceLabels.
setSize(origFaceLabels.size());
286 forAll(origFaceLabels, pointi)
288 procFaceLabels[pointi] = pointLookup[origFaceLabels[pointi]];
293 const labelList& curCellLabels = procCellAddressing_[proci];
297 cellList procCells(curCellLabels.size());
299 forAll(curCellLabels, celli)
301 const labelList& origCellLabels = meshCells[curCellLabels[celli]];
303 cell& curCell = procCells[celli];
305 curCell.
setSize(origCellLabels.size());
307 forAll(origCellLabels, cellFacei)
309 curCell[cellFacei] = faceLookup[origCellLabels[cellFacei]];
320 time().caseName()/(
"processor" +
Foam::name(proci)),
324 processorDb.setTime(time());
338 if (facesInstancePointsPtr_)
343 facesInstancePointsPtr_(),
357 std::move(facesInstancePoints),
358 std::move(procFaces),
374 std::move(procPoints),
375 std::move(procFaces),
383 const labelList& curPatchSizes = procPatchSize_[proci];
385 const labelList& curPatchStarts = procPatchStartIndex_[proci];
387 const labelList& curNeighbourProcessors =
388 procNeighbourProcessors_[proci];
390 const labelList& curProcessorPatchSizes =
391 procProcessorPatchSize_[proci];
393 const labelList& curProcessorPatchStarts =
394 procProcessorPatchStartIndex_[proci];
397 procProcessorPatchSubPatchIDs_[proci];
400 procProcessorPatchSubPatchStarts_[proci];
406 label nInterProcPatches = 0;
407 forAll(curSubPatchIDs, procPatchi)
409 nInterProcPatches += curSubPatchIDs[procPatchi].size();
414 curPatchSizes.size() + nInterProcPatches
419 forAll(curPatchSizes, patchi)
423 const polyPatch& meshPatch = meshPatches[patchi];
430 curPatchSizes[patchi],
431 curPatchStarts[patchi]
442 procMesh.boundaryMesh(),
444 patchMapper.directAddressing(),
445 curPatchStarts[patchi]
452 forAll(curProcessorPatchSizes, procPatchi)
454 const labelList& subPatchID = curSubPatchIDs[procPatchi];
455 const labelList& subStarts = curSubStarts[procPatchi];
457 label curStart = curProcessorPatchStarts[procPatchi];
463 i < subPatchID.size()-1
464 ? subStarts[i+1] - subStarts[i]
465 : curProcessorPatchSizes[procPatchi] - subStarts[i]
468 if (subPatchID[i] == -1)
479 procMesh.boundaryMesh(),
481 curNeighbourProcessors[procPatchi]
501 procMesh.boundaryMesh(),
503 curNeighbourProcessors[procPatchi],
516 procMesh.addPatches(procPatches);
532 zonePoints[zonei].setCapacity(pz[zonei].size()/nProcs_);
537 forAll(curPointLabels, pointi)
539 label curPoint = curPointLabels[pointi];
541 label zonei = pointToZone[curPoint];
546 zonePoints[zonei].append(pointi);
548 else if (zonei == -2)
553 label index = pz[zonei].whichPoint(curPoint);
557 zonePoints[zonei].append(pointi);
563 procMesh.pointZones().clearAddressing();
564 procMesh.pointZones().setSize(zonePoints.size());
567 procMesh.pointZones().set
572 procMesh.pointZones(),
574 zonePoints[zonei].shrink()
599 label procSize = fz[zonei].size() / nProcs_;
601 zoneFaces[zonei].setCapacity(procSize);
602 zoneFaceFlips[zonei].setCapacity(procSize);
608 forAll(curFaceLabels, facei)
612 label curF =
mag(curFaceLabels[facei]) - 1;
614 label zonei = faceToZone[curF];
619 zoneFaces[zonei].append(facei);
621 label index = fz[zonei].whichFace(curF);
623 bool flip = fz[zonei].flipMap()[index];
625 if (curFaceLabels[facei] < 0)
630 zoneFaceFlips[zonei].append(flip);
632 else if (zonei == -2)
637 label index = fz[zonei].whichFace(curF);
641 zoneFaces[zonei].append(facei);
643 bool flip = fz[zonei].flipMap()[index];
645 if (curFaceLabels[facei] < 0)
650 zoneFaceFlips[zonei].append(flip);
656 procMesh.faceZones().clearAddressing();
657 procMesh.faceZones().setSize(zoneFaces.size());
660 procMesh.faceZones().set
665 zoneFaces[zonei].shrink(),
666 zoneFaceFlips[zonei].shrink(),
692 zoneCells[zonei].setCapacity(cz[zonei].size()/nProcs_);
695 forAll(curCellLabels, celli)
697 label curCelli = curCellLabels[celli];
699 label zonei = cellToZone[curCelli];
704 zoneCells[zonei].append(celli);
706 else if (zonei == -2)
711 label index = cz[zonei].whichCell(curCelli);
715 zoneCells[zonei].append(celli);
721 procMesh.cellZones().clearAddressing();
722 procMesh.cellZones().setSize(zoneCells.size());
725 procMesh.cellZones().set
730 zoneCells[zonei].shrink(),
750 const auto* pMeshPtr =
751 thisDb().cfindObject<
pointMesh>(pointMesh::typeName);
754 const auto& pMesh = *pMeshPtr;
755 const auto& pMeshBoundary = pMesh.boundary();
767 bool differsFromPoly =
false;
770 forAll(pMeshBoundary, patchi)
773 if (mppPtr && (procBoundary.findPatchID(mppPtr->name()) == -1))
775 const auto& mpp = *mppPtr;
779 forAll(mpp.meshPoints(), i)
781 const label pointi = mpp.meshPoints()[i];
782 const label procPointi = pointLookup[pointi];
783 if (procPointi != -1)
785 procMeshPoints.append(procPointi);
786 procNormals.append(mpp.pointNormals()[i]);
790 procBoundary.push_back
799 meshPointPatch::typeName
802 differsFromPoly =
true;
809 forAll(procBoundary, patchi)
813 oldToNew[patchi] = newPatchi;
815 if (newPatchi != patchi)
817 differsFromPoly =
true;
826 boundaryProcAddressing.setSize(procBoundary.size(), -1);
828 forAll(procBoundary, patchi)
832 oldToNew[patchi] = newPatchi++;
835 procBoundary.reorder(oldToNew,
true);
841 procBoundary.write();
846 "boundaryProcAddressing",
847 procMesh.facesInstance(),
849 procPointMesh.thisDb()
856 if (facesInstancePointsPtr_)
870 std::move(procPoints)
872 pointsInstancePoints.write();
881 const cellSet& cs = cellSets[i];
882 cellSet set(procMesh, cs.name(), cs.size()/nProcs_);
885 if (cs.found(curCellLabels[i]))
894 const faceSet& cs = faceSets[i];
895 faceSet set(procMesh, cs.name(), cs.size()/nProcs_);
898 if (cs.found(
mag(curFaceLabels[i])-1))
908 pointSet set(procMesh, cs.name(), cs.size()/nProcs_);
911 if (cs.found(curPointLabels[i]))
935 procCellAddressing_[proci],
936 procPointAddressing_[proci]
941 Info<<
nl <<
"Processor " << proci;
943 if (procMesh.nCells())
952 Info<<
"Number of cells = " << procMesh.nCells() <<
nl;
954 if (procMesh.nCells())
956 Info<<
" Number of points = " << procMesh.nPoints() <<
nl;
959 maxProcCells =
max(maxProcCells, procMesh.nCells());
961 label nBoundaryFaces = 0;
963 label nProcFaces = 0;
965 for (
const polyPatch&
pp : procMesh.boundaryMesh())
971 const auto& procPatch = *cpp;
973 Info<<
" Number of faces shared with processor "
974 << procPatch.neighbProcNo() <<
" = "
975 << procPatch.size() <<
nl;
977 nProcFaces += procPatch.size();
982 nBoundaryFaces +=
pp.
size();
986 if (procMesh.nCells() && (nBoundaryFaces || nProcFaces))
989 <<
" Number of processor faces = " << nProcFaces <<
nl
990 <<
" Number of boundary faces = " << nBoundaryFaces <<
nl;
993 totProcFaces += nProcFaces;
995 maxProcFaces =
max(maxProcFaces, nProcFaces);
1003 procMesh.facesInstance(),
1012 ioAddr.rename(
"pointProcAddressing");
1016 ioAddr.rename(
"faceProcAddressing");
1020 ioAddr.rename(
"cellProcAddressing");
1025 label nMeshPatches = curPatchSizes.size();
1027 procBoundaryAddr.resize(nMeshPatches+
nProcPatches, -1);
1030 ioAddr.rename(
"boundaryProcAddressing");
1037 <<
"Number of processor faces = " << (totProcFaces/2) <<
nl
1038 <<
"Max number of cells = " << maxProcCells;
1040 if (maxProcCells != nCells())
1042 scalar avgValue = scalar(nCells())/nProcs_;
1044 Info<<
" (" << 100.0*(maxProcCells-avgValue)/avgValue
1045 <<
"% above average " << avgValue <<
')';
1049 Info<<
"Max number of processor patches = " << maxProcPatches;
1052 scalar avgValue = scalar(totProcPatches)/nProcs_;
1054 Info<<
" (" << 100.0*(maxProcPatches-avgValue)/avgValue
1055 <<
"% above average " << avgValue <<
')';
1059 Info<<
"Max number of faces between processors = " << maxProcFaces;
1062 scalar avgValue = scalar(totProcFaces)/nProcs_;
1064 Info<<
" (" << 100.0*(maxProcFaces-avgValue)/avgValue
1065 <<
"% above average " << avgValue <<
')';
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
static void writeContents(const IOobject &io, const UList< T > &content)
Write contents. The IOobject is never registered.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
label size() const noexcept
The number of elements in the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
A dynamically resizable PtrList with allocation management.
A non-owning sub-view of a List (allocated or unallocated storage).
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static word controlDictName
The default control dictionary name (normally "controlDict").
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A collection of cell labels.
A cell is defined as a list of faces with extra functionality.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Abstract base class for domain decomposition.
MeshObject wrapper of decompositionMethod.
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="", const dictionary *fallback=nullptr)
Read and register on mesh, optionally with alternative decomposeParDict path/name or with fallback co...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
domainDecomposition(const IOobject &io, const fileName &decompDictFile="")
Construct from components.
void updateParameters(const dictionary ¶ms)
Update flags based on the decomposition model settings.
const decompositionModel & model() const
Return decomposition model used.
bool writeDecomposition(const bool decomposeSets)
Write decomposition.
A face is a list of labels corresponding to mesh vertices.
A class for handling file names.
Patch field decomposer class.
Mesh data needed to do the Finite Volume discretisation.
Various for reading/decomposing/reconstructing/distributing refinement data.
pointPatch with explicitly provided points instead of using the points of a polyPatch.
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
Mesh representing a set of points created from polyMesh.
static word meshSubDir
Return the mesh sub-directory name (usually "pointMesh").
Mesh consisting of general polyhedral cells.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
A patch is a list of labels that address the faces in the global face list.
Neighbour processor patch.
Neighbour processor patch.
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
vectorIOField pointIOField
pointIOField is a vectorIOField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
List< cell > cellList
List of cell.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.