Loading...
Searching...
No Matches
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-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
31#include "GAMGAgglomeration.H"
32#include "Random.H"
33#include "lduMesh.H"
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 =
75 (
76 interfaces[intI]
77 );
78 label size = interfaces[intI].faceCells().size();
79 myNeighbours.insert(pp.neighbProcNo(), size);
80 }
81 }
82 }
83
84 Pstream::allGatherList(procFaces, Pstream::msgType(), mesh.comm());
85
86 autoPtr<lduPrimitiveMesh> singleCellMeshPtr;
87
88 if (Pstream::master(mesh.comm()))
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
146Foam::tmp<Foam::labelField>
147Foam::procFacesGAMGProcAgglomeration::processorAgglomeration
148(
149 const lduMesh& mesh
150) const
151{
152 UPstream::communicator singleCellMeshComm
153 (
154 mesh.comm(),
155 labelRange(1) // Processor 0 only
156 );
157
158 scalarField faceWeights;
159 autoPtr<lduPrimitiveMesh> singleCellMeshPtr
160 (
161 singleCellMesh
162 (
163 singleCellMeshComm.comm(),
164 mesh,
165 faceWeights
166 )
167 );
168
169 auto tfineToCoarse = tmp<labelField>::New();
170 auto& 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::broadcast(fineToCoarse, mesh.comm());
200
201 return tfineToCoarse;
202}
203
204
205bool Foam::procFacesGAMGProcAgglomeration::doProcessorAgglomeration
206(
207 const lduMesh& mesh
208) const
209{
210 // Check the need for further agglomeration on any processors
211 bool doAgg = mesh.lduAddr().size() < nAgglomeratingCells_;
212 UPstream::reduceOr(doAgg, mesh.comm());
213 return doAgg;
214}
215
216
217// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218
220(
221 GAMGAgglomeration& agglom,
223)
224:
226 nAgglomeratingCells_(controlDict.get<label>("nAgglomeratingCells"))
227{}
228
229
230// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233{}
234
235
236// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237
239{
240 if (debug)
241 {
242 Pout<< nl << "Starting mesh overview" << endl;
243 printStats(Pout, agglom_);
244 }
245
246 if (agglom_.size() >= 1)
247 {
248 Random rndGen(0);
249
250 for
251 (
252 label fineLevelIndex = 2;
253 fineLevelIndex < agglom_.size();
254 fineLevelIndex++
255 )
256 {
257 if (agglom_.hasMeshLevel(fineLevelIndex))
258 {
259 // Get the fine mesh
260 const lduMesh& levelMesh = agglom_.meshLevel(fineLevelIndex);
261
262 label levelComm = levelMesh.comm();
263 label nProcs = UPstream::nProcs(levelComm);
264
265 if (nProcs > 1 && doProcessorAgglomeration(levelMesh))
266 {
267 tmp<labelField> tprocAgglomMap
268 (
269 processorAgglomeration(levelMesh)
270 );
271 const labelField& procAgglomMap = tprocAgglomMap();
272
273 // Master processor
274 labelList masterProcs;
275 // Local processors that agglomerate. agglomProcIDs[0] is in
276 // masterProc.
277 List<label> agglomProcIDs;
279 (
280 levelComm,
281 procAgglomMap,
282 masterProcs,
283 agglomProcIDs
284 );
285
286 // Communicator for the processor-agglomerated matrix
287 comms_.push_back
288 (
290 (
291 levelComm,
292 masterProcs
293 )
294 );
295
296 // Use processor agglomeration maps to do the actual
297 // collecting.
299 (
300 fineLevelIndex,
301 procAgglomMap,
302 masterProcs,
303 agglomProcIDs,
304 comms_.back()
305 );
306 }
307 }
308 }
309 }
310
311 // Print a bit
312 if (debug)
313 {
314 Pout<< nl << "Agglomerated mesh overview" << endl;
315 printStats(Pout, agglom_);
316 }
317
318 return true;
319}
320
321
322// ************************************************************************* //
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())
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.
GAMGAgglomeration & agglom_
Reference to agglomeration.
void printStats(Ostream &os, GAMGAgglomeration &agglom) const
Debug: write agglomeration info.
DynamicList< label > comms_
Allocated communicators.
GAMGProcAgglomeration(const GAMGProcAgglomeration &)=delete
No copy construct.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
Random number generator.
Definition Random.H:56
Wrapper for internally indexed communicator label. Always invokes UPstream::allocateCommunicatorCompo...
Definition UPstream.H:2546
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static label newCommunicator(const label parent, const labelRange &subRanks, const bool withComponents=true)
Create new communicator with sub-ranks on the parent communicator.
Definition UPstream.C:272
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition fvMesh.C:707
virtual label comm() const
Return communicator used for parallel communication.
Definition fvMesh.H:415
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
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
virtual label comm() const =0
Return communicator used for parallel communication.
Simplest concrete lduMesh that stores the addressing needed by lduMatrix.
static tmp< labelField > agglomerate(label &nCoarseCells, const lduAddressing &fineMatrixAddressing, const scalarField &faceWeights)
Calculate and return agglomeration.
Processor agglomeration of GAMGAgglomerations. Needs nAgglomeratingCells which is when to start agglo...
procFacesGAMGProcAgglomeration(const procFacesGAMGProcAgglomeration &)=delete
No copy construct.
virtual bool agglomerate()
Modify agglomeration. Return true if modified.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
runTime controlDict().readEntry("adjustTimeStep"
dynamicFvMesh & mesh
Namespace for handling debugging switches.
Definition debug.C:45
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
UIndirectList< label > labelUIndList
UIndirectList of labels.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition lduSchedule.H:46
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
void sort(UList< T > &list)
Sort the list.
Definition UList.C:283
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Random rndGen