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 -------------------------------------------------------------------------------
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::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 
147 Foam::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 
205 bool 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,
222  const dictionary& controlDict
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 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
procFacesGAMGProcAgglomeration(const procFacesGAMGProcAgglomeration &)=delete
No copy construct.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Random rndGen
Definition: createFields.H:23
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void reduceOr(bool &value, const label communicator=worldComm)
Logical (or) reduction (MPI_AllReduce)
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:46
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
Macros for easy insertion into run-time selection tables.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static tmp< labelField > agglomerate(label &nCoarseCells, const lduAddressing &fineMatrixAddressing, const scalarField &faceWeights)
Calculate and return agglomeration.
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
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:1065
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:415
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:741
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
runTime controlDict().readEntry("adjustTimeStep"
The central control dictionary, the contents of which are either taken directly from the FOAM_CONTROL...
Definition: debug.C:142
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
virtual bool agglomerate()
Modify agglomeration. Return true if modified.
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:29
virtual bool agglomerate()=0
Modify agglomeration.
Processor agglomeration of GAMGAgglomerations.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
constexpr label labelMax
Definition: label.H:55
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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:
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static label allocateCommunicator(const label parent, const labelRange &subRanks, const bool withComponents=true)
Allocate new communicator with contiguous sub-ranks on the parent communicator.
Definition: UPstream.C:258
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
label size() const
Return number of equations.