Go to the documentation of this file.
43 assemblyFaceAreaPairGAMGAgglomeration,
50 assemblyFaceAreaPairGAMGAgglomeration,
75 if (isA<lduPrimitiveMeshAssembly>(ldumesh))
78 refCast<const lduPrimitiveMeshAssembly>(ldumesh);
84 for (label i=0; i <
mesh.meshes().
size(); ++i)
86 const fvMesh& m = refCast<const fvMesh>(
mesh.meshes()[i]);
92 faceAreas[subFaceMap[facei]] = areas[facei];
100 const polyPatch& pp =
patches[patchI];
101 label globalPatchID =
mesh.patchMap()[i][patchI];
103 if (globalPatchID == -1)
105 if (pp.masterImplicit())
109 if (isA<cyclicAMIPolyPatch>(pp))
111 const cyclicAMIPolyPatch& mpp =
112 refCast<const cyclicAMIPolyPatch>(pp);
115 mpp.AMI().srcWeights();
118 forAll(pp.faceCells(), faceI)
122 for(label j=0; j<w.size(); j++)
124 const label globalFaceI =
125 mesh.faceBoundMap()[i][patchI][subFaceI];
127 if (globalFaceI != -1)
129 faceAreas[globalFaceI] = w[j]*sf[faceI];
135 else if (isA<cyclicACMIPolyPatch>(pp))
137 const cyclicACMIPolyPatch& mpp =
138 refCast<const cyclicACMIPolyPatch>(pp);
141 mpp.AMI().srcWeights();
144 const scalar tol = mpp.tolerance();
150 for(label j=0; j<w.size(); j++)
152 if (mask[faceI] > tol)
154 const label globalFaceI =
155 mesh.faceBoundMap()[i]
158 faceAreas[globalFaceI] = w[j]*sf[faceI];
166 forAll(pp.faceCells(), faceI)
168 const label globalFaceI =
169 mesh.faceBoundMap()[i][patchI][faceI];
171 if (globalFaceI != -1)
173 faceAreas[globalFaceI] = sf[faceI];
189 faceAreas/
sqrt(
mag(faceAreas)),
198 const fvMesh& fvmesh = refCast<const fvMesh>(matrix.
mesh());
207 fvmesh.Sf().primitiveField()
208 /
sqrt(fvmesh.magSf().primitiveField()),
221 const lduMatrix& matrix,
236 faceAreas/
sqrt(
mag(faceAreas)),
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
virtual bool masterImplicit() const
List< scalarList > scalarListList
A List of scalarList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< scalar > scalarList
A List of scalars.
assemblyFaceAreaPairGAMGAgglomeration(const lduMatrix &matrix, const dictionary &controlDict)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
static constexpr const zero Zero
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
const Internal::FieldType & primitiveField() const
static scalar tolerance()
const labelListList & faceMap(const label fineLeveli) const
void agglomerate(const lduMesh &mesh, const scalarField &faceWeights)
const polyBoundaryMesh & boundaryMesh() const
const Internal & internalField() const
virtual const labelUList & upperAddr() const =0
Field< vector > vectorField
Specialisation of Field<T> for vector.
Generic templated field type.
const surfaceScalarField & magSf() const
const lduMesh & mesh() const
A patch is a list of labels that address the faces in the global face list.
virtual const lduAddressing & lduAddr() const
runTime controlDict().readEntry("adjustTimeStep"
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
const labelUList & faceCells() const
Vector< scalar > vector
A scalar version of the templated Vector.
const fvBoundaryMesh & boundary() const
const AMIPatchToPatchInterpolation & AMI() const
~assemblyFaceAreaPairGAMGAgglomeration()
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const lduMesh & mesh() const
const polyBoundaryMesh & patches
Agglomerate using the pair algorithm.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
const scalarField & mask() const
defineTypeNameAndDebug(combustionModel, 0)
virtual const lduAddressing & lduAddr() const =0
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
const scalarListList & srcWeights() const
Cyclic patch for Arbitrary Mesh Interface (AMI)
const surfaceVectorField & Sf() const