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 -------------------------------------------------------------------------------
11 License
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"
34 #include "processorLduInterface.H"
35 #include "processorGAMGInterface.H"
36 #include "pairGAMGAgglomeration.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(procFacesGAMGProcAgglomeration, 0);
43 
45  (
46  GAMGProcAgglomeration,
47  procFacesGAMGProcAgglomeration,
48  GAMGAgglomeration
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 Foam::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 
84  Pstream::gatherList(procFaces, Pstream::msgType(), mesh.comm());
85  Pstream::scatterList(procFaces, Pstream::msgType(), mesh.comm());
86 
87  autoPtr<lduPrimitiveMesh> singleCellMeshPtr;
88 
89  if (Pstream::master(mesh.comm()))
90  {
91  // I am master
92  label nCells = Pstream::nProcs(mesh.comm());
93 
94  DynamicList<label> l(3*nCells);
95  DynamicList<label> u(3*nCells);
96  DynamicList<scalar> weight(3*nCells);
97 
98  DynamicList<label> nbrs;
99  DynamicList<scalar> weights;
100 
101  forAll(procFaces, proci)
102  {
103  const Map<label>& neighbours = procFaces[proci];
104 
105  // Add all the higher processors
106  nbrs.clear();
107  weights.clear();
108  forAllConstIters(neighbours, iter)
109  {
110  if (iter.key() > proci)
111  {
112  nbrs.append(iter.key());
113  weights.append(iter());
114  }
115  sort(nbrs);
116  forAll(nbrs, i)
117  {
118  l.append(proci);
119  u.append(nbrs[i]);
120  weight.append(weights[i]);
121  }
122  }
123  }
124 
125  faceWeights.transfer(weight);
126 
127  PtrList<const lduInterface> primitiveInterfaces(0);
128  const lduSchedule ps(0);
129 
130  singleCellMeshPtr.reset
131  (
132  new lduPrimitiveMesh
133  (
134  nCells,
135  l,
136  u,
137  primitiveInterfaces,
138  ps,
139  singleCellMeshComm
140  )
141  );
142  }
143  return singleCellMeshPtr;
144 }
145 
146 
148 Foam::procFacesGAMGProcAgglomeration::processorAgglomeration
149 (
150  const lduMesh& mesh
151 ) const
152 {
153  label singleCellMeshComm = UPstream::allocateCommunicator
154  (
155  mesh.comm(),
156  labelList(1, Zero) // only processor 0
157  );
158 
159  scalarField faceWeights;
160  autoPtr<lduPrimitiveMesh> singleCellMeshPtr
161  (
162  singleCellMesh
163  (
164  singleCellMeshComm,
165  mesh,
166  faceWeights
167  )
168  );
169 
170  tmp<labelField> tfineToCoarse(new labelField(0));
171  labelField& fineToCoarse = tfineToCoarse.ref();
172 
173  if (singleCellMeshPtr)
174  {
175  // On master call the agglomerator
176  const lduPrimitiveMesh& singleCellMesh = *singleCellMeshPtr;
177 
178  label nCoarseProcs;
180  (
181  nCoarseProcs,
182  singleCellMesh,
183  faceWeights
184  );
185 
186  labelList coarseToMaster(nCoarseProcs, labelMax);
187  forAll(fineToCoarse, celli)
188  {
189  label coarseI = fineToCoarse[celli];
190  coarseToMaster[coarseI] = min(coarseToMaster[coarseI], celli);
191  }
192 
193  // Sort according to master and redo restriction
194  labelList newToOld(sortedOrder(coarseToMaster));
195  labelList oldToNew(invert(newToOld.size(), newToOld));
196 
197  fineToCoarse = labelUIndList(oldToNew, fineToCoarse)();
198  }
199 
200  Pstream::scatter(fineToCoarse, Pstream::msgType(), mesh.comm());
201  UPstream::freeCommunicator(singleCellMeshComm);
202 
203  return tfineToCoarse;
204 }
205 
206 
207 bool Foam::procFacesGAMGProcAgglomeration::doProcessorAgglomeration
208 (
209  const lduMesh& mesh
210 ) const
211 {
212  // Check the need for further agglomeration on all processors
213  bool doAgg = mesh.lduAddr().size() < nAgglomeratingCells_;
214  mesh.reduce(doAgg, orOp<bool>());
215  return doAgg;
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
220 
221 Foam::procFacesGAMGProcAgglomeration::procFacesGAMGProcAgglomeration
222 (
223  GAMGAgglomeration& agglom,
224  const dictionary& controlDict
225 )
226 :
228  nAgglomeratingCells_(controlDict.get<label>("nAgglomeratingCells"))
229 {}
230 
231 
232 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233 
235 {
236  forAllReverse(comms_, i)
237  {
238  if (comms_[i] != -1)
239  {
240  UPstream::freeCommunicator(comms_[i]);
241  }
242  }
243 }
244 
245 
246 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247 
249 {
250  if (debug)
251  {
252  Pout<< nl << "Starting mesh overview" << endl;
253  printStats(Pout, agglom_);
254  }
255 
256  if (agglom_.size() >= 1)
257  {
258  Random rndGen(0);
259 
260  for
261  (
262  label fineLevelIndex = 2;
263  fineLevelIndex < agglom_.size();
264  fineLevelIndex++
265  )
266  {
267  if (agglom_.hasMeshLevel(fineLevelIndex))
268  {
269  // Get the fine mesh
270  const lduMesh& levelMesh = agglom_.meshLevel(fineLevelIndex);
271 
272  label levelComm = levelMesh.comm();
273  label nProcs = UPstream::nProcs(levelComm);
274 
275  if (nProcs > 1 && doProcessorAgglomeration(levelMesh))
276  {
277  tmp<labelField> tprocAgglomMap
278  (
279  processorAgglomeration(levelMesh)
280  );
281  const labelField& procAgglomMap = tprocAgglomMap();
282 
283  // Master processor
284  labelList masterProcs;
285  // Local processors that agglomerate. agglomProcIDs[0] is in
286  // masterProc.
287  List<label> agglomProcIDs;
289  (
290  levelComm,
291  procAgglomMap,
292  masterProcs,
293  agglomProcIDs
294  );
295 
296  // Allocate a communicator for the processor-agglomerated
297  // matrix
298  comms_.append
299  (
301  (
302  levelComm,
303  masterProcs
304  )
305  );
306 
307 
308  // Use processor agglomeration maps to do the actual
309  // collecting.
311  (
312  fineLevelIndex,
313  procAgglomMap,
314  masterProcs,
315  agglomProcIDs,
316  comms_.last()
317  );
318  }
319  }
320  }
321  }
322 
323  // Print a bit
324  if (debug)
325  {
326  Pout<< nl << "Agglomerated mesh overview" << endl;
327  printStats(Pout, agglom_);
328  }
329 
330  return true;
331 }
332 
333 
334 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::Random
Random number generator.
Definition: Random.H:59
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::lduSchedule
List< lduScheduleEntry > lduSchedule
Definition: lduSchedule.H:76
Foam::GAMGAgglomeration
Geometric agglomerated algebraic multigrid agglomeration class.
Definition: GAMGAgglomeration.H:64
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::pairGAMGAgglomeration::agglomerate
void agglomerate(const lduMesh &mesh, const scalarField &faceWeights)
Agglomerate all levels starting from the given face weights.
Definition: pairGAMGAgglomerate.C:34
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Definition: gatherScatter.C:150
GAMGAgglomeration.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::lduAddressing::size
label size() const
Return number of equations.
Definition: lduAddressing.H:171
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::invert
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
processorLduInterface.H
Foam::GAMGAgglomeration::calculateRegionMaster
static void calculateRegionMaster(const label comm, const labelList &procAgglomMap, labelList &masterProcs, List< label > &agglomProcIDs)
Given fine to coarse processor map determine:
Definition: GAMGAgglomerateLduAddressing.C:737
Foam::fvMesh::comm
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:316
Foam::UPstream::allocateCommunicator
static label allocateCommunicator(const label parent, const labelList &subRanks, const bool doPstream=true)
Allocate a new communicator.
Definition: UPstream.C:108
Foam::Field< label >
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:691
Foam::GAMGProcAgglomeration::agglomerate
virtual bool agglomerate()=0
Modify agglomeration. Return true if modified.
controlDict
runTime controlDict().readEntry("adjustTimeStep"
Definition: debug.C:143
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
Foam::procFacesGAMGProcAgglomeration::agglomerate
virtual bool agglomerate()
Modify agglomeration. Return true if modified.
Definition: procFacesGAMGProcAgglomeration.C:248
Foam::UPstream::freeCommunicator
static void freeCommunicator(const label communicator, const bool doPstream=true)
Free a previously allocated communicator.
Definition: UPstream.C:174
pairGAMGAgglomeration.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
processorGAMGInterface.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Random.H
procFacesGAMGProcAgglomeration.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::lduMesh::reduce
void reduce(T &Value, const BinaryOp &bop) const
Helper: reduce with current communicator.
Definition: lduMeshTemplates.C:34
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::procFacesGAMGProcAgglomeration::~procFacesGAMGProcAgglomeration
virtual ~procFacesGAMGProcAgglomeration()
Destructor.
Definition: procFacesGAMGProcAgglomeration.C:234
Foam::List< label >
Foam::lduMesh::comm
virtual label comm() const =0
Return communicator used for parallel communication.
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:309
rndGen
Random rndGen
Definition: createFields.H:23
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::GAMGProcAgglomeration
Processor agglomeration of GAMGAgglomerations.
Definition: GAMGProcAgglomeration.H:54
lduMesh.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62
Foam::lduInterfacePtrsList
UPtrList< const lduInterface > lduInterfacePtrsList
List of coupled interface fields to be used in coupling.
Definition: lduInterfacePtrsList.H:44
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58