Go to the documentation of this file.
78 Info<<
"GAMGAgglomeration:" <<
nl
79 <<
" local agglomerator : " <<
type() <<
nl;
82 Info<<
" processor agglomerator : "
88 <<
setw(20) <<
"nFaces/nCells"
89 <<
setw(20) <<
"nInterfaces"
90 <<
setw(20) <<
"nIntFaces/nCells"
91 <<
setw(12) <<
"profile"
94 <<
setw(8) <<
"nProcs"
110 <<
setw(8) <<
"-----"
111 <<
setw(8) <<
"------"
129 for (
label levelI = 0; levelI <=
size(); levelI++)
133 scalar faceCellRatio = 0;
134 label nInterfaces = 0;
137 scalar profile = 0.0;
152 if (interfaces.
set(i))
155 nIntFaces += interfaces[i].faceCells().
size();
158 ratio = scalar(nIntFaces)/
nCells;
168 scalar maxFaceCellRatio =
170 scalar totFaceCellRatio =
181 int oldPrecision =
Info().precision(4);
184 <<
setw(8) << totNprocs
186 <<
setw(8) << totNCells/totNprocs
187 <<
setw(8) << maxNCells
189 <<
setw(8) << totFaceCellRatio/totNprocs
190 <<
setw(8) << maxFaceCellRatio
192 <<
setw(8) << scalar(totNInt)/totNprocs
193 <<
setw(8) << maxNInt
195 <<
setw(8) << totRatio/totNprocs
196 <<
setw(8) << maxRatio
197 <<
setw(12) << totProfile/totNprocs
200 Info().precision(oldPrecision);
209 const label nCoarseCells
213 bool contAgg = nCoarseCells >= nCellsInCoarsestLevel_;
231 nCellsInCoarsestLevel_
252 restrictAddressing_(maxLevels_),
254 faceRestrictAddressing_(maxLevels_),
255 faceFlipMap_(maxLevels_),
256 nPatchFaces_(maxLevels_),
257 patchFaceRestrictAddressing_(maxLevels_),
259 meshLevels_(maxLevels_)
261 procCommunicator_.setSize(maxLevels_ + 1, -1);
262 if (processorAgglomerate())
264 procAgglomMap_.setSize(maxLevels_);
265 agglomProcIDs_.setSize(maxLevels_);
266 procCellOffsets_.setSize(maxLevels_);
267 procFaceMap_.setSize(maxLevels_);
268 procBoundaryMap_.setSize(maxLevels_);
269 procBoundaryFaceMap_.setSize(maxLevels_);
284 GAMGAgglomeration::typeName
288 const word agglomeratorType
296 "geometricGAMGAgglomerationLibs",
297 lduMeshConstructorTablePtr_
300 lduMeshConstructorTable::iterator cstrIter =
301 lduMeshConstructorTablePtr_->
find(agglomeratorType);
303 if (cstrIter == lduMeshConstructorTablePtr_->end())
306 <<
"Unknown GAMGAgglomeration type "
307 << agglomeratorType <<
".\n"
308 <<
"Valid matrix GAMGAgglomeration types are :"
309 << lduMatrixConstructorTablePtr_->sortedToc() <<
endl
310 <<
"Valid geometric GAMGAgglomeration types are :"
311 << lduMeshConstructorTablePtr_->sortedToc()
321 GAMGAgglomeration::typeName
339 GAMGAgglomeration::typeName
343 const word agglomeratorType
351 "algebraicGAMGAgglomerationLibs",
352 lduMatrixConstructorTablePtr_
357 !lduMatrixConstructorTablePtr_
358 || !lduMatrixConstructorTablePtr_->found(agglomeratorType)
365 lduMatrixConstructorTable::iterator cstrIter =
366 lduMatrixConstructorTablePtr_->find(agglomeratorType);
368 return store(cstrIter()(matrix,
controlDict).ptr());
375 GAMGAgglomeration::typeName
389 const word agglomeratorType
397 "geometricGAMGAgglomerationLibs",
398 geometryConstructorTablePtr_
401 geometryConstructorTable::iterator cstrIter =
402 geometryConstructorTablePtr_->
find(agglomeratorType);
404 if (cstrIter == geometryConstructorTablePtr_->end())
407 <<
"Unknown GAMGAgglomeration type "
408 << agglomeratorType <<
".\n"
409 <<
"Valid geometric GAMGAgglomeration types are :"
410 << geometryConstructorTablePtr_->sortedToc()
446 return meshLevels_[i - 1];
459 return meshLevels_.set(i - 1);
471 return meshInterfaces_;
475 return meshLevels_[i - 1].rawInterfaces();
484 meshLevels_.set(i - 1, NULL);
486 if (i < nCells_.size())
489 restrictAddressing_.set(i, NULL);
491 faceRestrictAddressing_.set(i, NULL);
492 faceFlipMap_.set(i, NULL);
493 nPatchFaces_.set(i, NULL);
494 patchFaceRestrictAddressing_.set(i, NULL);
505 return procAgglomMap_[leveli];
514 return agglomProcIDs_[leveli];
520 return procCommunicator_[leveli] != -1;
526 return procCommunicator_[leveli];
535 return procCellOffsets_[leveli];
544 return procFaceMap_[leveli];
553 return procBoundaryMap_[leveli];
562 return procBoundaryFaceMap_[leveli];
575 if (fineAddressing.
size() != restrict.
size())
578 <<
"nCells:" << fineAddressing.
size()
579 <<
" agglom:" << restrict.
size()
596 label own = lower[faceI];
597 label nei = upper[faceI];
599 if (restrict[own] == restrict[nei])
603 if (master[own] < master[nei])
605 master[nei] = master[own];
608 else if (master[own] > master[nei])
610 master[own] = master[nei];
630 labelList& masters = coarseToMasters[restrict[cellI]];
632 if (
findIndex(masters, master[cellI]) == -1)
634 masters.
append(master[cellI]);
639 if (nNewCoarse > nCoarse)
650 nNewCoarse = nCoarse;
652 forAll(coarseToMasters, coarseI)
654 const labelList& masters = coarseToMasters[coarseI];
656 labelList& newCoarse = coarseToNewCoarse[coarseI];
658 newCoarse[0] = coarseI;
659 for (
label i = 1; i < newCoarse.
size(); i++)
661 newCoarse[i] = nNewCoarse++;
668 label coarseI = restrict[cellI];
671 newRestrict[cellI] = coarseToNewCoarse[coarseI][index];
The class contains the addressing required by the lduMatrix: upper, lower and losort.
static bool checkRestriction(labelList &newRestrict, label &nNewCoarse, const lduAddressing &fineAddressing, const labelUList &restrict, const label nCoarse)
Given restriction determines if coarse cells are connected.
PtrList< labelListList > procFaceMap_
Mapping from processor to procMeshLevel face.
autoPtr< GAMGProcAgglomeration > procAgglomeratorPtr_
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from string.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Geometric agglomerated algebraic multigrid agglomeration class.
#define forAll(list, i)
Loop across all elements in list.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
labelList nCells_
The number of cells in each level.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
static label nProcs(const label communicator=0)
Number of processes in parallel run.
const Time & time() const
Return time.
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
const labelListList & faceMap(const label fineLeveli) const
Mapping from processor to procMesh face.
Ostream & endl(Ostream &os)
Add newline and flush stream.
PtrList< labelList > faceRestrictAddressing_
Face restriction addressing array.
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
PtrList< lduPrimitiveMesh > meshLevels_
Hierarchy of mesh addressing.
label size() const
Return number of equations.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
const labelList & cellOffsets(const label fineLeveli) const
Mapping from processor to procMesh cells.
bool continueAgglomerating(const label nCoarseCells) const
Check the need for further agglomeration.
virtual lduInterfacePtrsList interfaces() const =0
Return a list of pointers for each patch.
const labelList & agglomProcIDs(const label fineLeveli) const
Set of processors to agglomerate. Element 0 is the.
virtual label comm() const
Return communicator used for parallel communication.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const labelListList & boundaryMap(const label fineLeveli) const
Mapping from processor to procMesh boundary.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
runTime controlDict().lookup("adjustTimeStep") >> adjustTimeStep
void append(const T &)
Append an element at the end of the list.
void clearLevel(const label leveli)
Istream and Ostream manipulators taking arguments.
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.
PtrList< labelList > nPatchFaces_
The number of (coarse) patch faces in each level.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
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.
errorManip< error > abort(error &err)
Omanip< int > setw(const int i)
static const GAMGAgglomeration & New(const lduMesh &mesh, const dictionary &controlDict)
Return the selected geometric agglomerator.
GAMGAgglomeration(const GAMGAgglomeration &)
Disallow default bitwise copy construct.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void setSize(const label)
Reset size of List.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool foundObject(const word &name) const
Is the named Type found?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
void reduce(T &Value, const BinaryOp &bop) const
Helper: reduce with current communicator.
const labelList & procAgglomMap(const label fineLeveli) const
Mapping from processor to agglomerated processor (global, all.
labelList nFaces_
The number of (coarse) faces in each level.
const lduInterfacePtrsList & interfaceLevel(const label leveli) const
Return LDU interface addressing of given level.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
PtrList< labelListList > procBoundaryMap_
Mapping from processor to procMeshLevel boundary.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void compactLevels(const label nCreatedLevels)
Shrink the number of levels to that specified.
bool set(const label) const
Is element set.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
const labelListListList & boundaryFaceMap(const label fineLeveli) const
Mapping from processor to procMesh boundary face.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
bool hasProcMesh(const label fineLeveli) const
Check that level has combined mesh.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
void size(const label)
Override size to be inconsistent with allocated storage.
label procCommunicator(const label fineLeveli) const
Communicator for current level or -1.
label size() const
Return the number of elements in the UList.
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)
PtrList< boolList > faceFlipMap_
Face flip: for faces mapped to internal faces stores whether.
PtrList< labelListListList > procBoundaryFaceMap_
Mapping from processor to procMeshLevel boundary face.
defineTypeNameAndDebug(combustionModel, 0)
PtrList< labelListList > patchFaceRestrictAddressing_
Patch-local face restriction addressing array.
PtrList< labelList > procCellOffsets_
Mapping from processor to procMeshLevel cells.
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
bool processorAgglomerate() const
Whether to agglomerate across processors.
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...
Tuple2< label, scalar > band() const
Calculate bandwidth and profile of addressing.
label size() const
Return the number of elements in the UPtrList.
~GAMGAgglomeration()
Destructor.
labelList procCommunicator_
Communicator for given level.