51void Foam::MGridGenGAMGAgglomeration::swap
53 const lduInterfacePtrsList& interfaces,
54 const labelUList& cellValues,
55 PtrList<labelList>& nbrValues
63 if (interfaces.set(inti))
65 interfaces[inti].initInternalFieldTransfer
79 nbrValues.setSize(interfaces.size());
82 if (interfaces.set(inti))
89 interfaces[inti].internalFieldTransfer
101void Foam::MGridGenGAMGAgglomeration::getNbrAgglom
109 cellToNbrAgglom.setSize(addr.size());
110 cellToNbrAgglom = -1;
114 if (interfaces.set(inti))
121 if (pldui.myProcNo() > pldui.neighbProcNo())
124 interfaces[inti].faceCells();
125 const labelList& nbrData = nbrGlobalAgglom[inti];
129 cellToNbrAgglom[
faceCells[i]] = nbrData[i];
138void Foam::MGridGenGAMGAgglomeration::detectSharedFaces
150 sharedFaces.reserve(addr.lowerAddr().size()/100);
155 label lowerData = value[
lower[facei]];
156 label upperData = value[
upper[facei]];
158 if (lowerData != -1 && lowerData == upperData)
160 sharedFaces.insert(facei);
168Foam::MGridGenGAMGAgglomeration::MGridGenGAMGAgglomeration
178 nProcConsistencyIter_
195void Foam::MGridGenGAMGAgglomeration::agglomerate
197 const label nCellsInCoarsestLevel,
198 const label startLevel,
200 const bool doProcessorAgglomerate
220 const labelList& own = fvMesh_.faceOwner();
225 if (!fvMesh_.isInternalFace(facei))
227 magSb[own[facei]] +=
mag(Sf[facei]);
235 label nCreatedLevels = 0;
237 while (nCreatedLevels < maxLevels_ - 1)
239 label nCoarseCells = -1;
246 meshLevel(nCreatedLevels).lduAddr(),
253 for (
int i=0; i<nProcConsistencyIter_; i++)
267 globalAgglom[celli] = globalNumbering.toGlobal(agglom[celli]);
272 swap(interfaces, globalAgglom, nbrGlobalAgglom);
278 getNbrAgglom(addr, interfaces, nbrGlobalAgglom, cellToNbrAgglom);
283 detectSharedFaces(
mesh, cellToNbrAgglom, sharedFaces);
290 for (
const label facei : sharedFaces)
292 weights[facei] *= 2.0;
296 finalAgglomPtr = agglomerate
301 meshLevel(nCreatedLevels).lduAddr(),
310 continueAgglomerating
312 nCellsInCoarsestLevel_,
313 finalAgglomPtr().size(),
315 meshLevel(nCreatedLevels).comm()
319 nCells_[nCreatedLevels] = nCoarseCells;
320 restrictAddressing_.set(nCreatedLevels, finalAgglomPtr);
327 agglomerateLduAddressing(nCreatedLevels);
333 new scalarField(meshLevels_[nCreatedLevels].size())
337 restrictField(*aggVPtr, *VPtr, nCreatedLevels,
false);
353 meshLevels_[nCreatedLevels].upperAddr().size(),
358 restrictFaceField(*aggMagSfPtr, *magSfPtr, nCreatedLevels);
365 magSfPtr = aggMagSfPtr;
372 new scalarField(meshLevels_[nCreatedLevels].size())
376 restrictField(*aggMagSbPtr, *magSbPtr, nCreatedLevels,
false);
379 magSbPtr = aggMagSbPtr;
386 compactLevels(nCreatedLevels, doProcessorAgglomerate);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static const Field< scalar > & null() noexcept
Geometric agglomerated algebraic multigrid agglomeration class.
void restrictField(Field< Type > &cf, const Field< Type > &ff, const label fineLevelIndex, const bool procAgglom) const
Restrict (integrate by summation) cell field.
void agglomerateLduAddressing(const label fineLevelIndex)
Assemble coarse mesh addressing.
const label maxLevels_
Max number of levels.
void compactLevels(const label nCreatedLevels, const bool doProcessorAgglomerate)
Shrink the number of levels to that specified. Optionally do.
label nCellsInCoarsestLevel_
Number of cells in coarsest level.
PtrList< lduPrimitiveMesh > meshLevels_
Hierarchy of mesh addressing.
GAMGAgglomeration(const GAMGAgglomeration &)=delete
No copy construct.
labelList nCells_
The number of cells in each level.
bool continueAgglomerating(const label nCellsInCoarsestLevel, const label nCells, const label nCoarseCells, const label comm) const
Check the need for further agglomeration.
void restrictFaceField(Field< Type > &cf, const Field< Type > &ff, const label fineLevelIndex) const
Restrict (integrate by summation) face field.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
const lduMesh & meshLevel(const label leveli) const
Return LDU mesh of given level.
Agglomerate using the MGridGen algorithm.
const lduMesh & mesh() const noexcept
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
SubField is a Field obtained as a section of another Field, without its own allocation....
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
static void waitRequests()
Wait for all requests to finish.
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.
Mesh data needed to do the Finite Volume discretisation.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
label toGlobal(const label proci, const label i) const
From local to global on proci.
The class contains the addressing required by the lduMatrix: upper, lower and losort.
label size() const noexcept
Return number of equations.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
An abstract base class for processor coupled interfaces.
A class for managing temporary objects.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
runTime controlDict().readEntry("adjustTimeStep"
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
static constexpr const zero Zero
Global zero (0).
UList< label > labelUList
A UList of labels.
#define forAll(list, i)
Loop across all elements in list.