Loading...
Searching...
No Matches
assemblyFaceAreaPairGAMGAgglomeration.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "lduMatrix.H"
32#include "surfaceFields.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39
41 (
45 );
46
48 (
51 geometry
52 );
53}
54
55
56// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
60{}
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
66(
67 const lduMatrix& matrix,
69)
70:
72{
73 const lduMesh& ldumesh = matrix.mesh();
74
76 {
79
80 vectorField faceAreas(mesh.lduAddr().upperAddr().size(), Zero);
81
82 const labelListList& faceMap = mesh.faceMap();
83
84 for (label i=0; i < mesh.meshes().size(); ++i)
85 {
86 const fvMesh& m = refCast<const fvMesh>(mesh.meshes()[i]);
87 const labelList& subFaceMap = faceMap[i];
88 const vectorField& areas = m.Sf().internalField();
89
90 forAll(subFaceMap, facei)
91 {
92 faceAreas[subFaceMap[facei]] = areas[facei];
93 }
94
95 const polyBoundaryMesh& patches = m.boundaryMesh();
96
97 // Fill faceAreas for new faces
98 forAll(patches, patchI)
99 {
100 const polyPatch& pp = patches[patchI];
101 label globalPatchID = mesh.patchMap()[i][patchI];
102
103 if (globalPatchID == -1)
104 {
105 if (pp.masterImplicit())
106 {
107 const vectorField& sf = m.boundary()[patchI].Sf();
108
110 {
111 const cyclicAMIPolyPatch& mpp =
113
114 const scalarListList& srcWeight =
115 mpp.AMI().srcWeights();
116
117 label subFaceI = 0;
118 forAll(pp.faceCells(), faceI)
119 {
120 const scalarList& w = srcWeight[faceI];
121
122 for(label j=0; j<w.size(); j++)
123 {
124 const label globalFaceI =
125 mesh.faceBoundMap()[i][patchI][subFaceI];
126
127 if (globalFaceI != -1)
128 {
129 faceAreas[globalFaceI] = w[j]*sf[faceI];
130 }
131 subFaceI++;
132 }
133 }
134 }
136 {
137 const cyclicACMIPolyPatch& mpp =
139
140 const scalarListList& srcWeight =
141 mpp.AMI().srcWeights();
142
143 const scalarField& mask = mpp.mask();
144 const scalar tol = mpp.tolerance();
145 label subFaceI = 0;
146 forAll(mask, faceI)
147 {
148 const scalarList& w = srcWeight[faceI];
149
150 for(label j=0; j<w.size(); j++)
151 {
152 if (mask[faceI] > tol)
153 {
154 const label globalFaceI =
155 mesh.faceBoundMap()[i]
156 [patchI][subFaceI];
157
158 faceAreas[globalFaceI] = w[j]*sf[faceI];
159 }
160 subFaceI++;
161 }
162 }
163 }
164 else
165 {
166 forAll(pp.faceCells(), faceI)
167 {
168 const label globalFaceI =
169 mesh.faceBoundMap()[i][patchI][faceI];
170
171 if (globalFaceI != -1)
172 {
173 faceAreas[globalFaceI] = sf[faceI];
174 }
175 }
176 }
177 }
178 }
179 }
180 }
181
183 (
185 0, //mesh,
186 mag
187 (
189 (
190 faceAreas/sqrt(mag(faceAreas)),
191 vector(1, 1.01, 1.02)
192 )
193 ),
194 true
195 );
196 }
197 else
198 {
199 // Revert to faceAreaPairGAMGAgglomeration
200 const fvMesh& fvmesh = refCast<const fvMesh>(matrix.mesh());
201
203 (
205 0, //matrix.mesh(),
206 mag
207 (
209 (
210 fvmesh.Sf().primitiveField()
211 /sqrt(fvmesh.magSf().primitiveField()),
212 vector(1, 1.01, 1.02)
213 //vector::one
214 )
215 ),
216 true
217 );
218 }
219}
220
221
224(
225 const lduMatrix& matrix,
226 const scalarField& cellVolumes,
227 const vectorField& faceAreas,
229)
230:
231 pairGAMGAgglomeration(matrix.mesh(), controlDict)
232{
234 (
236 0, //matrix.mesh(),
237 mag
238 (
240 (
241 faceAreas/sqrt(mag(faceAreas)),
242 vector(1, 1.01, 1.02)
243 )
244 ),
245 true
246 );
247}
248
249
250// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const scalarListList & srcWeights() const
Return const access to source patch weights.
Geometric agglomerated algebraic multigrid agglomeration class.
label nCellsInCoarsestLevel_
Number of cells in coarsest level.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
assemblyFaceAreaPairGAMGAgglomeration(const lduMatrix &matrix, const dictionary &controlDict)
Construct given mesh and controls.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
static scalar tolerance()
Overlap tolerance.
Cyclic patch for Arbitrary Mesh Interface (AMI).
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const surfaceVectorField & Sf() const
Return cell face area vectors.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition lduMatrix.H:753
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
Agglomerate using the pair algorithm.
static tmp< labelField > agglomerate(label &nCoarseCells, const lduAddressing &fineMatrixAddressing, const scalarField &faceWeights)
Calculate and return agglomeration.
pairGAMGAgglomeration(const pairGAMGAgglomeration &)=delete
No copy construct.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
runTime controlDict().readEntry("adjustTimeStep"
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
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.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.