Go to the documentation of this file.
47 for (
label levelI = 0; levelI <= agglom.
size(); levelI++)
55 os <<
"Level " << levelI <<
" has no fine mesh:" <<
endl;
69 os <<
"Level " << levelI <<
" agglomeration:" <<
nl
70 <<
" nCoarseCells:" << agglom.
nCells(levelI) <<
nl
71 <<
" nCoarseFaces:" << agglom.
nFaces(levelI) <<
nl
72 <<
" cellRestriction:"
73 <<
" size:" << cellRestrict.
size()
74 <<
" max:" <<
max(cellRestrict)
76 <<
" faceRestriction:"
77 <<
" size:" << faceRestrict.
size()
78 <<
" max:" <<
max(faceRestrict)
83 forAll(patchFaceRestrict, i)
85 if (patchFaceRestrict[i].size())
90 <<
" size:" << faceRestrict.
size()
91 <<
" max:" <<
max(faceRestrict)
126 const label myProcID = Pstream::myProcNo(
mesh.comm());
137 forAll(globalIndices, cellI)
139 globalIndices[cellI] = globalNumbering.
toGlobal(myProcID, cellI);
149 if (interfaces.
set(inti))
151 interfaces[inti].initInternalFieldTransfer
153 Pstream::nonBlocking,
159 if (Pstream::parRun())
161 Pstream::waitRequests();
166 if (interfaces.
set(inti))
173 interfaces[inti].internalFieldTransfer
175 Pstream::nonBlocking,
202 if (interfaces.
set(inti))
204 const labelUList& faceCells = interfaces[inti].faceCells();
208 nNbrs[faceCells[i]]++;
220 cellCells[cellI].
setSize(nNbrs[cellI], -1);
229 label c0 = own[faceI];
232 cellCells[c0][nNbrs[c0]++] = globalIndices[
c1];
233 cellCells[
c1][nNbrs[
c1]++] = globalIndices[c0];
237 if (interfaces.
set(inti))
239 const labelUList& faceCells = interfaces[inti].faceCells();
243 label c0 = faceCells[i];
244 cellCells[c0][nNbrs[c0]++] = nbrGlobalCells[inti][i];
257 cellCells[cellI][0] = globalIndices[cellI];
266 const label fineLevelIndex,
270 const label procAgglomComm
273 const lduMesh& levelMesh = agglom_.meshLevels_[fineLevelIndex];
276 if (Pstream::myProcNo(levelComm) != -1)
281 agglom_.procAgglomerateLduAddressing
294 label levelI = fineLevelIndex+1;
295 levelI < agglom_.meshLevels_.size();
299 agglom_.procAgglomerateRestrictAddressing
307 if (Pstream::myProcNo(levelComm) == agglomProcIDs[0])
312 label levelI = fineLevelIndex;
313 levelI < agglom_.meshLevels_.size();
317 agglom_.agglomerateLduAddressing(levelI);
325 label levelI = fineLevelIndex+1;
326 levelI <= agglom_.size();
330 agglom_.clearLevel(levelI);
361 Info<<
"GAMGProcAgglomeration::New(const word&, GAMGAgglomeration&"
362 ", const dictionary&) : "
363 "constructing GAMGProcAgglomeration"
367 GAMGAgglomerationConstructorTable::iterator cstrIter =
368 GAMGAgglomerationConstructorTablePtr_->find(
type);
370 if (cstrIter == GAMGAgglomerationConstructorTablePtr_->end())
373 <<
"Unknown GAMGProcAgglomeration type "
374 <<
type <<
" for GAMGAgglomeration " << agglom.type() <<
nl <<
nl
375 <<
"Valid GAMGProcAgglomeration types are :" <<
endl
376 << GAMGAgglomerationConstructorTablePtr_->sortedToc()
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const labelList & faceRestrictAddressing(const label leveli) const
Return face restrict addressing of given level.
A class for handling words, derived from string.
Geometric agglomerated algebraic multigrid agglomeration class.
label nFaces(const label leveli) const
Return number of coarse faces (before processor agglomeration)
#define forAll(list, i)
Loop across all elements in list.
GAMGProcAgglomeration(const GAMGProcAgglomeration &)
Disallow default bitwise copy construct.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
const labelListList & patchFaceRestrictAddressing(const label leveli) const
Ostream & endl(Ostream &os)
Add newline and flush stream.
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
label size() const
Return number of equations.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
const labelField & restrictAddressing(const label leveli) const
Return cell restrict addressing of given level.
void stableSort(UList< T > &)
bool set(const label) const
Is element set.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
runTime controlDict().lookup("adjustTimeStep") >> adjustTimeStep
virtual bool agglomerate()=0
Modify agglomeration. Return true if modified.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
static autoPtr< GAMGProcAgglomeration > New(const word &type, GAMGAgglomeration &agglom, const dictionary &controlDict)
Return the selected agglomerator.
void printStats(Ostream &os, GAMGAgglomeration &agglom) const
Debug: write agglomeration info.
A list of keyword definitions, which are a keyword followed by any number of values (e....
PtrList< labelList > procAgglomMap_
Per level, per processor the processor it agglomerates into.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
errorManipArg< error, int > exit(error &err, const int errNo=1)
InfoProxy< lduMesh > info() const
Return info proxy.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
bool set(const label) const
Is element set.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
virtual label comm() const =0
Return communicator used for parallel communication.
static labelListList globalCellCells(const lduMesh &)
Debug: calculate global cell-cells.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
label nCells(const label leveli) const
Return number of coarse cells (before processor agglomeration)
defineTypeNameAndDebug(combustionModel, 0)
PtrList< labelList > procCellOffsets_
Mapping from processor to procMeshLevel cells.
PtrList< labelList > agglomProcIDs_
Per level the set of processors to agglomerate. Element 0 is.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
virtual ~GAMGProcAgglomeration()
Destructor.
label toGlobal(const label i) const
From local to global.
label size() const
Return the number of elements in the UPtrList.
labelList procCommunicator_
Communicator for given level.