Loading...
Searching...
No Matches
MGridGenGAMGAgglomeration.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "fvMesh.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39
41 (
45 );
46}
47
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50
51void Foam::MGridGenGAMGAgglomeration::swap
52(
53 const lduInterfacePtrsList& interfaces,
54 const labelUList& cellValues,
55 PtrList<labelList>& nbrValues
56) const
57{
58 const label startOfRequests = UPstream::nRequests();
59
60 // Initialise transfer of restrict addressing on the interface
61 forAll(interfaces, inti)
62 {
63 if (interfaces.set(inti))
64 {
65 interfaces[inti].initInternalFieldTransfer
66 (
68 cellValues
69 );
70 }
71 }
72
73 // Wait for outstanding requests
74 // (commsType == UPstream::commsTypes::nonBlocking)
75 UPstream::waitRequests(startOfRequests);
76
77
78 // Get the interface agglomeration
79 nbrValues.setSize(interfaces.size());
80 forAll(interfaces, inti)
81 {
82 if (interfaces.set(inti))
83 {
84 nbrValues.set
85 (
86 inti,
87 new labelList
88 (
89 interfaces[inti].internalFieldTransfer
90 (
92 cellValues
93 )
94 )
95 );
96 }
97 }
98}
99
100
101void Foam::MGridGenGAMGAgglomeration::getNbrAgglom
102(
103 const lduAddressing& addr,
104 const lduInterfacePtrsList& interfaces,
105 const PtrList<labelList>& nbrGlobalAgglom,
106 labelList& cellToNbrAgglom
107) const
108{
109 cellToNbrAgglom.setSize(addr.size());
110 cellToNbrAgglom = -1;
111
112 forAll(interfaces, inti)
113 {
114 if (interfaces.set(inti))
115 {
116 if (isA<processorLduInterface>(interfaces[inti]))
117 {
118 const processorLduInterface& pldui =
119 refCast<const processorLduInterface>(interfaces[inti]);
120
121 if (pldui.myProcNo() > pldui.neighbProcNo())
122 {
123 const labelUList& faceCells =
124 interfaces[inti].faceCells();
125 const labelList& nbrData = nbrGlobalAgglom[inti];
126
127 forAll(faceCells, i)
128 {
129 cellToNbrAgglom[faceCells[i]] = nbrData[i];
130 }
131 }
132 }
133 }
134 }
135}
136
137
138void Foam::MGridGenGAMGAgglomeration::detectSharedFaces
139(
140 const lduMesh& mesh,
141 const labelList& value,
142 labelHashSet& sharedFaces
143) const
144{
145 const lduAddressing& addr = mesh.lduAddr();
146 const labelUList& lower = addr.lowerAddr();
147 const labelUList& upper = addr.upperAddr();
148
149 sharedFaces.clear();
150 sharedFaces.reserve(addr.lowerAddr().size()/100);
151
152 // Detect any faces inbetween same value
153 forAll(lower, facei)
154 {
155 label lowerData = value[lower[facei]];
156 label upperData = value[upper[facei]];
157
158 if (lowerData != -1 && lowerData == upperData)
159 {
160 sharedFaces.insert(facei);
162 }
163}
164
165
166// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167
168Foam::MGridGenGAMGAgglomeration::MGridGenGAMGAgglomeration
169(
170 const lduMesh& mesh,
172)
173:
175 fvMesh_(refCast<const fvMesh>(mesh)),
176 minSize_(controlDict.get<label>("minSize")),
177 maxSize_(controlDict.get<label>("maxSize")),
178 nProcConsistencyIter_
179 (
180 controlDict.get<label>("nProcConsistencyIter")
181 )
182{
183 agglomerate
184 (
185 1, //nCellsInCoarsestLevel
186 0,
188 true
189 );
190}
191
192
193// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194
195void Foam::MGridGenGAMGAgglomeration::agglomerate
196(
197 const label nCellsInCoarsestLevel,
198 const label startLevel,
199 const scalarField& startFaceWeights,
200 const bool doProcessorAgglomerate
201)
202{
203 // Start geometric agglomeration from the cell volumes and areas of the mesh
204 scalarField* VPtr = const_cast<scalarField*>(&fvMesh_.cellVolumes());
205
206 scalarField magFaceAreas(sqrt(3.0)*mag(fvMesh_.faceAreas()));
207 SubField<scalar> magSf(magFaceAreas, fvMesh_.nInternalFaces());
208
209 scalarField* magSfPtr = const_cast<scalarField*>
210 (
211 &magSf.operator const scalarField&()
212 );
213
214 // Create the boundary area cell field
215 scalarField* magSbPtr(new scalarField(fvMesh_.nCells(), Zero));
216
217 {
218 scalarField& magSb = *magSbPtr;
219
220 const labelList& own = fvMesh_.faceOwner();
221 const vectorField& Sf = fvMesh_.faceAreas();
222
223 forAll(Sf, facei)
224 {
225 if (!fvMesh_.isInternalFace(facei))
226 {
227 magSb[own[facei]] += mag(Sf[facei]);
228 }
229 }
230 }
231
232 // Agglomerate until the required number of cells in the coarsest level
233 // is reached
234
235 label nCreatedLevels = 0;
236
237 while (nCreatedLevels < maxLevels_ - 1)
238 {
239 label nCoarseCells = -1;
240
241 tmp<labelField> finalAgglomPtr = agglomerate
242 (
243 nCoarseCells,
244 minSize_,
245 maxSize_,
246 meshLevel(nCreatedLevels).lduAddr(),
247 *VPtr,
248 *magSfPtr,
249 *magSbPtr
250 );
251
252 // Adjust weights only
253 for (int i=0; i<nProcConsistencyIter_; i++)
254 {
255 const lduMesh& mesh = meshLevel(nCreatedLevels);
256 const lduAddressing& addr = mesh.lduAddr();
257 const lduInterfacePtrsList interfaces = mesh.interfaces();
258
259 const labelField& agglom = finalAgglomPtr();
260
261 // Global numbering
262 const globalIndex globalNumbering(nCoarseCells);
263
264 labelField globalAgglom(addr.size());
265 forAll(agglom, celli)
266 {
267 globalAgglom[celli] = globalNumbering.toGlobal(agglom[celli]);
268 }
269
270 // Get the interface agglomeration
271 PtrList<labelList> nbrGlobalAgglom;
272 swap(interfaces, globalAgglom, nbrGlobalAgglom);
273
274
275 // Get the interface agglomeration on a cell basis (-1 for all
276 // other cells)
277 labelList cellToNbrAgglom;
278 getNbrAgglom(addr, interfaces, nbrGlobalAgglom, cellToNbrAgglom);
279
280
281 // Mark all faces inbetween cells with same nbragglomeration
282 labelHashSet sharedFaces(addr.size()/100);
283 detectSharedFaces(mesh, cellToNbrAgglom, sharedFaces);
284
285
286 //- Note: in-place update of weights is more effective it seems?
287 // Should not be. fluke?
288 //scalarField weights(*faceWeightsPtr);
289 scalarField weights = *magSfPtr;
290 for (const label facei : sharedFaces)
291 {
292 weights[facei] *= 2.0;
293 }
294
295 // Redo the agglomeration using the new weights
296 finalAgglomPtr = agglomerate
297 (
298 nCoarseCells,
299 minSize_,
300 maxSize_,
301 meshLevel(nCreatedLevels).lduAddr(),
302 *VPtr,
303 weights,
304 *magSbPtr
305 );
306 }
307
308 if
309 (
310 continueAgglomerating
311 (
312 nCellsInCoarsestLevel_,
313 finalAgglomPtr().size(),
314 nCoarseCells,
315 meshLevel(nCreatedLevels).comm()
316 )
317 )
318 {
319 nCells_[nCreatedLevels] = nCoarseCells;
320 restrictAddressing_.set(nCreatedLevels, finalAgglomPtr);
321 }
322 else
323 {
324 break;
325 }
326
327 agglomerateLduAddressing(nCreatedLevels);
328
329 // Agglomerate the cell volumes field for the next level
330 {
331 scalarField* aggVPtr
332 (
333 new scalarField(meshLevels_[nCreatedLevels].size())
334 );
335
336 // Restrict but no parallel agglomeration (not supported)
337 restrictField(*aggVPtr, *VPtr, nCreatedLevels, false);
338
339 if (nCreatedLevels)
340 {
341 delete VPtr;
342 }
343
344 VPtr = aggVPtr;
345 }
346
347 // Agglomerate the face areas field for the next level
348 {
349 scalarField* aggMagSfPtr
350 (
351 new scalarField
352 (
353 meshLevels_[nCreatedLevels].upperAddr().size(),
354 0
355 )
356 );
357
358 restrictFaceField(*aggMagSfPtr, *magSfPtr, nCreatedLevels);
359
360 if (nCreatedLevels)
361 {
362 delete magSfPtr;
363 }
364
365 magSfPtr = aggMagSfPtr;
366 }
367
368 // Agglomerate the cell boundary areas field for the next level
369 {
370 scalarField* aggMagSbPtr
371 (
372 new scalarField(meshLevels_[nCreatedLevels].size())
373 );
374
375 // Restrict but no parallel agglomeration (not supported)
376 restrictField(*aggMagSbPtr, *magSbPtr, nCreatedLevels, false);
377
378 delete magSbPtr;
379 magSbPtr = aggMagSbPtr;
380 }
381
382 nCreatedLevels++;
383 }
384
385 // Shrink the storage of the levels to those created
386 compactLevels(nCreatedLevels, doProcessorAgglomerate);
387
388 // Delete temporary geometry storage
389 if (nCreatedLevels)
390 {
391 delete VPtr;
392 delete magSfPtr;
393 }
394 delete magSbPtr;
395}
396
397
398// ************************************************************************* //
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
Definition Field.H:192
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.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition SubField.H:58
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition fvMesh.C:724
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition fvMesh.C:707
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
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...
Definition lduMesh.H:54
An abstract base class for processor coupled interfaces.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
runTime controlDict().readEntry("adjustTimeStep"
dynamicFvMesh & mesh
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.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
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.
Definition typeInfo.H:87
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.
Definition labelField.H:48
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299