50 for (label levelI = 0; levelI <= agglom.
size(); levelI++)
54 os <<
"Level " << levelI <<
" mesh:"
59 os <<
"Level " << levelI <<
" has no fine mesh:" <<
endl;
73 os <<
"Level " << levelI <<
" agglomeration:" <<
nl
74 <<
" nCoarseCells:" << agglom.
nCells(levelI) <<
nl
75 <<
" nCoarseFaces:" << agglom.
nFaces(levelI) <<
nl
76 <<
" cellRestriction:"
77 <<
" size:" << cellRestrict.size()
78 <<
" max:" <<
max(cellRestrict)
80 <<
" faceRestriction:"
81 <<
" size:" << faceRestrict.size()
82 <<
" max:" <<
max(faceRestrict)
87 forAll(patchFaceRestrict, i)
89 if (patchFaceRestrict[i].size())
94 <<
" size:" << faceRestrict.
size()
95 <<
" max:" <<
max(faceRestrict)
127 const lduAddressing& addr =
mesh.lduAddr();
132 const globalIndex globalNumbering
145 PtrList<labelList> nbrGlobalCells(interfaces.size());
152 if (interfaces.set(inti))
154 interfaces[inti].initInternalFieldTransfer
166 if (interfaces.set(inti))
173 interfaces[inti].internalFieldTransfer
202 if (interfaces.set(inti))
220 cellCells[celli].setSize(nNbrs[celli], -1);
229 label c0 = own[facei];
230 label
c1 = nbr[facei];
232 cellCells[c0][nNbrs[c0]++] = globalIndices[
c1];
233 cellCells[
c1][nNbrs[
c1]++] = globalIndices[c0];
237 if (interfaces.set(inti))
244 cellCells[c0][nNbrs[c0]++] = nbrGlobalCells[inti][i];
257 cellCells[celli][0] = globalIndices[celli];
266 const label fineLevelIndex,
269 const List<label>& agglomProcIDs,
270 const label procAgglomComm
273 const lduMesh& levelMesh = agglom_.meshLevels_[fineLevelIndex];
274 label levelComm = levelMesh.comm();
281 agglom_.procAgglomerateLduAddressing
294 label levelI = fineLevelIndex+1;
295 levelI < agglom_.meshLevels_.size();
299 agglom_.procAgglomerateRestrictAddressing
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 auto* ctorPtr = GAMGAgglomerationConstructorTable(
type);
366 <<
"Unknown GAMGProcAgglomeration type "
367 <<
type <<
" for GAMGAgglomeration " << agglom.type() <<
nl <<
nl
368 <<
"Valid GAMGProcAgglomeration types :" <<
endl
369 << GAMGAgglomerationConstructorTablePtr_->sortedToc()
Geometric agglomerated algebraic multigrid agglomeration class.
const labelList & faceRestrictAddressing(const label leveli) const
Return face restrict addressing of given level.
const labelListList & patchFaceRestrictAddressing(const label leveli) const
labelList procCommunicator_
Communicator for given level.
PtrList< labelList > agglomProcIDs_
Per level the set of processors to agglomerate. Element 0 is.
PtrList< labelList > procAgglomMap_
Per level, per processor the processor it agglomerates into.
PtrList< labelList > procCellOffsets_
Mapping from processor to procMeshLevel cells.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
const labelField & restrictAddressing(const label leveli) const
Return cell restrict addressing of given level.
bool hasMeshLevel(const label leveli) const
Do we have mesh for given level?
label nCells(const label leveli) const
Return number of coarse cells (before processor agglomeration).
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
label nFaces(const label leveli) const
Return number of coarse faces (before processor agglomeration).
Processor agglomeration of GAMGAgglomerations.
virtual bool agglomerate()=0
Modify agglomeration.
GAMGAgglomeration & agglom_
Reference to agglomeration.
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.
void clearCommunicators()
Clear/free allocated communicators.
DynamicList< label > comms_
Allocated communicators.
static labelListList globalCellCells(const lduMesh &)
Debug: calculate global cell-cells.
virtual ~GAMGProcAgglomeration()
Destructor. Clears allocated communicators.
GAMGProcAgglomeration(const GAMGProcAgglomeration &)=delete
No copy construct.
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().
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
void size(const label n)
Older name for setAddressableSize.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static void waitRequests()
Wait for all requests to finish.
static void freeCommunicator(const label communicator, const bool withComponents=true)
Free a previously allocated communicator.
static bool & parRun() noexcept
Test if this a parallel run.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
label size() const noexcept
The number of entries in the list.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Smooth ATC in cells next to a set of patches supplied by type.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
labelRange range(label proci) const noexcept
Return start/size range of proci data.
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
label size() const noexcept
Return number of equations.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
virtual label comm() const =0
Return communicator used for parallel communication.
InfoProxy< lduMesh > info() const noexcept
Return info proxy, used to print mesh information to a stream.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
runTime controlDict().readEntry("adjustTimeStep"
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
#define DebugInFunction
Report an information message using Foam::Info.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
Ostream & endl(Ostream &os)
Add newline and flush stream.
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...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
UList< label > labelUList
A UList of labels.
void stableSort(UList< T > &list)
Stable sort the list.
constexpr char nl
The newline '\n' character (0x0a).
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.