procFacesGAMGProcAgglomeration.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 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
31#include "GAMGAgglomeration.H"
32#include "Random.H"
33#include "lduMesh.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
43
45 (
49 );
50}
51
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
56Foam::procFacesGAMGProcAgglomeration::singleCellMesh
57(
58 const label singleCellMeshComm,
59 const lduMesh& mesh,
60 scalarField& faceWeights
61) const
62{
63 // Count number of faces per processor
64 List<Map<label>> procFaces(UPstream::nProcs(mesh.comm()));
65 Map<label>& myNeighbours = procFaces[UPstream::myProcNo(mesh.comm())];
66
67 {
68 const lduInterfacePtrsList interfaces(mesh.interfaces());
69 forAll(interfaces, intI)
70 {
71 if (interfaces.set(intI))
72 {
73 const processorLduInterface& pp =
74 refCast<const processorLduInterface>
75 (
76 interfaces[intI]
77 );
78 label size = interfaces[intI].faceCells().size();
79 myNeighbours.insert(pp.neighbProcNo(), size);
80 }
81 }
82 }
83
85
86 autoPtr<lduPrimitiveMesh> singleCellMeshPtr;
87
89 {
90 // I am master
91 label nCells = Pstream::nProcs(mesh.comm());
92
93 DynamicList<label> l(3*nCells);
94 DynamicList<label> u(3*nCells);
95 DynamicList<scalar> weight(3*nCells);
96
97 DynamicList<label> nbrs;
98 DynamicList<scalar> weights;
99
100 forAll(procFaces, proci)
101 {
102 const Map<label>& neighbours = procFaces[proci];
103
104 // Add all the higher processors
105 nbrs.clear();
106 weights.clear();
107 forAllConstIters(neighbours, iter)
108 {
109 if (iter.key() > proci)
110 {
111 nbrs.append(iter.key());
112 weights.append(iter());
113 }
114 sort(nbrs);
115 forAll(nbrs, i)
116 {
117 l.append(proci);
118 u.append(nbrs[i]);
119 weight.append(weights[i]);
120 }
121 }
122 }
123
124 faceWeights.transfer(weight);
125
126 PtrList<const lduInterface> primitiveInterfaces(0);
127 const lduSchedule ps(0);
128
129 singleCellMeshPtr.reset
130 (
131 new lduPrimitiveMesh
132 (
133 nCells,
134 l,
135 u,
136 primitiveInterfaces,
137 ps,
138 singleCellMeshComm
139 )
140 );
141 }
142 return singleCellMeshPtr;
143}
144
145
147Foam::procFacesGAMGProcAgglomeration::processorAgglomeration
148(
149 const lduMesh& mesh
150) const
151{
152 label singleCellMeshComm = UPstream::allocateCommunicator
153 (
154 mesh.comm(),
155 labelList(1, Zero) // only processor 0
156 );
157
158 scalarField faceWeights;
159 autoPtr<lduPrimitiveMesh> singleCellMeshPtr
160 (
161 singleCellMesh
162 (
163 singleCellMeshComm,
164 mesh,
165 faceWeights
166 )
167 );
168
169 tmp<labelField> tfineToCoarse(new labelField(0));
170 labelField& fineToCoarse = tfineToCoarse.ref();
171
172 if (singleCellMeshPtr)
173 {
174 // On master call the agglomerator
175 const lduPrimitiveMesh& singleCellMesh = *singleCellMeshPtr;
176
177 label nCoarseProcs;
179 (
180 nCoarseProcs,
181 singleCellMesh,
182 faceWeights
183 );
184
185 labelList coarseToMaster(nCoarseProcs, labelMax);
186 forAll(fineToCoarse, celli)
187 {
188 label coarseI = fineToCoarse[celli];
189 coarseToMaster[coarseI] = min(coarseToMaster[coarseI], celli);
190 }
191
192 // Sort according to master and redo restriction
193 labelList newToOld(sortedOrder(coarseToMaster));
194 labelList oldToNew(invert(newToOld.size(), newToOld));
195
196 fineToCoarse = labelUIndList(oldToNew, fineToCoarse)();
197 }
198
199 Pstream::scatter(fineToCoarse, Pstream::msgType(), mesh.comm());
200 UPstream::freeCommunicator(singleCellMeshComm);
201
202 return tfineToCoarse;
203}
204
205
206bool Foam::procFacesGAMGProcAgglomeration::doProcessorAgglomeration
207(
208 const lduMesh& mesh
209) const
210{
211 // Check the need for further agglomeration on all processors
212 bool doAgg = mesh.lduAddr().size() < nAgglomeratingCells_;
213 mesh.reduce(doAgg, orOp<bool>());
214 return doAgg;
215}
216
217
218// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219
221(
222 GAMGAgglomeration& agglom,
224)
225:
227 nAgglomeratingCells_(controlDict.get<label>("nAgglomeratingCells"))
228{}
229
230
231// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232
234{
235 forAllReverse(comms_, i)
236 {
237 if (comms_[i] != -1)
238 {
240 }
241 }
242}
243
244
245// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
246
248{
249 if (debug)
250 {
251 Pout<< nl << "Starting mesh overview" << endl;
252 printStats(Pout, agglom_);
253 }
254
255 if (agglom_.size() >= 1)
256 {
257 Random rndGen(0);
258
259 for
260 (
261 label fineLevelIndex = 2;
262 fineLevelIndex < agglom_.size();
263 fineLevelIndex++
264 )
265 {
266 if (agglom_.hasMeshLevel(fineLevelIndex))
267 {
268 // Get the fine mesh
269 const lduMesh& levelMesh = agglom_.meshLevel(fineLevelIndex);
270
271 label levelComm = levelMesh.comm();
272 label nProcs = UPstream::nProcs(levelComm);
273
274 if (nProcs > 1 && doProcessorAgglomeration(levelMesh))
275 {
276 tmp<labelField> tprocAgglomMap
277 (
278 processorAgglomeration(levelMesh)
279 );
280 const labelField& procAgglomMap = tprocAgglomMap();
281
282 // Master processor
283 labelList masterProcs;
284 // Local processors that agglomerate. agglomProcIDs[0] is in
285 // masterProc.
286 List<label> agglomProcIDs;
288 (
289 levelComm,
290 procAgglomMap,
291 masterProcs,
292 agglomProcIDs
293 );
294
295 // Allocate a communicator for the processor-agglomerated
296 // matrix
297 comms_.append
298 (
300 (
301 levelComm,
302 masterProcs
303 )
304 );
305
306
307 // Use processor agglomeration maps to do the actual
308 // collecting.
310 (
311 fineLevelIndex,
312 procAgglomMap,
313 masterProcs,
314 agglomProcIDs,
315 comms_.last()
316 );
317 }
318 }
319 }
320 }
321
322 // Print a bit
323 if (debug)
324 {
325 Pout<< nl << "Agglomerated mesh overview" << endl;
326 printStats(Pout, agglom_);
327 }
328
329 return true;
330}
331
332
333// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Geometric agglomerated algebraic multigrid agglomeration class.
static void calculateRegionMaster(const label comm, const labelList &procAgglomMap, labelList &masterProcs, List< label > &agglomProcIDs)
Given fine to coarse processor map determine:
Processor agglomeration of GAMGAgglomerations.
virtual bool agglomerate()=0
Modify agglomeration. Return true if modified.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void allGatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void scatter(const List< commsStruct > &comms, T &value, const int tag, const label comm)
Broadcast data: Distribute without modification.
Random number generator.
Definition: Random.H:60
T & last()
Return the last element of the list.
Definition: UListI.H:216
static void freeCommunicator(const label communicator, const bool doPstream=true)
Free a previously allocated communicator.
Definition: UPstream.C:174
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
static label allocateCommunicator(const label parent, const labelList &subRanks, const bool doPstream=true)
Allocate a new communicator.
Definition: UPstream.C:108
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:735
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:718
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:326
label size() const
Return number of equations.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:63
virtual label comm() const =0
Return communicator used for parallel communication.
void reduce(T &Value, const BinaryOp &bop) const
Helper: reduce with current communicator.
Processor agglomeration of GAMGAgglomerations. Needs nAgglomeratingCells which is when to start agglo...
virtual bool agglomerate()
Modify agglomeration. Return true if modified.
int myProcNo() const noexcept
Return processor number.
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
runTime controlDict().readEntry("adjustTimeStep"
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition: label.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:342
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:36
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:68
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Random rndGen
Definition: createFields.H:23