processorColour.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) 2022 M. Janssens
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "processorColour.H"
29 #include "processorLduInterface.H"
30 #include "processorTopologyNew.H"
31 
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
44 (
45  const lduMesh& mesh,
46  labelList& procColour
47 )
48 {
49  procColour.resize_nocopy(Pstream::nProcs(mesh.comm()));
50  procColour = -1;
51 
52  // Re-use processor-topology analysis
53 
54  labelListList procNeighbours(Pstream::nProcs(mesh.comm()));
55 
56  // Fill my entry
57  {
59 
60  auto& procToProcs = procNeighbours[Pstream::myProcNo(mesh.comm())];
61  label n = 0;
62  forAll(patches, patchi)
63  {
64  if (patches.set(patchi))
65  {
66  if (isA<processorLduInterface>(patches[patchi]))
67  {
68  n++;
69  }
70  }
71  }
72 
73  procToProcs.resize_nocopy(n);
74  n = 0;
75  forAll(patches, patchi)
76  {
77  if (patches.set(patchi))
78  {
79  const auto* ppPtr = isA<processorLduInterface>(patches[patchi]);
80  if (ppPtr)
81  {
82  procToProcs[n++] = ppPtr->neighbProcNo();
83  }
84  }
85  }
86  }
87  // Send to master
88  Pstream::gatherList(procNeighbours, UPstream::msgType(), mesh.comm());
89 
90 
91  // Use greedy algorithm for now
92  // (see e.g. https://iq.opengenus.org/graph-colouring-greedy-algorithm/)
93 
94  if (Pstream::master(mesh.comm()))
95  {
96  DynamicList<label> front;
97  front.append(0); // start from processor 0
98 
99  DynamicList<label> newFront;
100  while (front.size())
101  {
102  //Pout<< "Front:" << front << endl;
103 
104  newFront.clear();
105  for (const label proci : front)
106  {
107  if (procColour[proci] == -1)
108  {
109  const labelList& nbrs = procNeighbours[proci];
110  const UIndirectList<label> nbrColour(procColour, nbrs);
111 
112  for
113  (
114  label colouri = 0;
115  colouri < Pstream::nProcs();
116  colouri++
117  )
118  {
119  if (!nbrColour.found(colouri))
120  {
121  procColour[proci] = colouri;
122  for (label nbrProci : nbrs)
123  {
124  if (procColour[nbrProci] == -1)
125  {
126  newFront.append(nbrProci);
127  }
128  }
129  break;
130  }
131  }
132  }
133  }
134 
135  front = std::move(newFront);
136  }
137  }
138 
139  Pstream::broadcast(procColour, mesh.comm());
140 
141 
142  //if (false)
143  //{
144  // volScalarField volColour
145  // (
146  // IOobject
147  // (
148  // "colour",
149  // mesh.time().timeName(),
150  // mesh,
151  // IOobject::NO_READ,
152  // IOobject::NO_WRITE,
153  // IOobject::NO_REGISTER
154  // ),
155  // mesh,
156  // dimensionedScalar(dimless, procColour[Pstream::myProcNo()]),
157  // fvPatchFieldBase::zeroGradientType()
158  // );
159  // volColour.write();
160  //}
161 
162  const label nColours = max(procColour)+1;
163 
164  if (debug)
165  {
166  Info<< typeName << " : coloured " << Pstream::nProcs(mesh.comm())
167  << " processors with in total " << nColours << " colours" << endl;
168  }
169 
170  return nColours;
171 }
172 
173 
175 (
176  const lduMesh& mesh,
177  DynamicList<label>& front,
178  labelList& cellColour
179 )
180 {
181  // Colour with the (min) coupled (global) patch
182 
183  const lduAddressing& addr = mesh.lduAddr();
184  const label* const __restrict__ uPtr = addr.upperAddr().begin();
185  const label* const __restrict__ lPtr = addr.lowerAddr().begin();
186  const label* const __restrict__ ownStartPtr = addr.ownerStartAddr().begin();
187  const label* const __restrict__ losortStartAddrPtr =
188  addr.losortStartAddr().begin();
189  const label* const __restrict__ losortAddrPtr = addr.losortAddr().begin();
190 
191  DynamicList<label> newFront;
192  while (true)
193  {
194  newFront.clear();
195  for (const label celli : front)
196  {
197  const label colouri = cellColour[celli];
198 
199  {
200  const label fStart = ownStartPtr[celli];
201  const label fEnd = ownStartPtr[celli + 1];
202 
203  for (label facei=fStart; facei<fEnd; facei++)
204  {
205  const label nbr =
206  (
207  lPtr[facei] == celli
208  ? uPtr[facei]
209  : lPtr[facei]
210  );
211  if (cellColour[nbr] == -1)
212  {
213  cellColour[nbr] = colouri;
214  newFront.append(nbr);
215  }
216  }
217  }
218  {
219  const label fStart = losortStartAddrPtr[celli];
220  const label fEnd = losortStartAddrPtr[celli + 1];
221 
222  for (label i=fStart; i<fEnd; i++)
223  {
224  label facei = losortAddrPtr[i];
225  const label nbr =
226  (
227  lPtr[facei] == celli
228  ? uPtr[facei]
229  : lPtr[facei]
230  );
231  if (cellColour[nbr] == -1)
232  {
233  cellColour[nbr] = colouri;
234  newFront.append(nbr);
235  }
236  }
237  }
238 
239  }
240 
241  if (newFront.empty())
242  {
243  break;
244  }
246  front.transfer(newFront);
247  }
248 }
249 
250 
252 (
253  const lduMesh& mesh,
254  labelList& cellColour,
255  labelList& patchToColour
256 )
257 {
258  // Colour with the (min) coupled (global) patch. Compact to have colour 0
259  // if a single region. Problematic if same patch on multiple disconnected
260  // regions.
261 
262  const lduAddressing& addr = mesh.lduAddr();
264 
265  const label nCells = addr.size();
266 
267  patchToColour.resize_nocopy(patches.size());
268  patchToColour = -1;
269  cellColour.resize_nocopy(nCells);
270  cellColour = -1;
271 
272 
273  label colouri = 0;
274 
275  // Starting front from patch faceCells
276  DynamicList<label> front;
277  forAll(patches, inti)
278  {
279  if
280  (
281  patches.set(inti)
282  && !isA<processorLduInterface>(patches[inti])
283  )
284  {
285  // 'global' interface. Seed faceCells with patch index
286 
287  if (patchToColour[inti] == -1)
288  {
289  patchToColour[inti] = colouri++;
290  }
291 
292  const auto& fc = patches[inti].faceCells();
293  forAll(fc, face)
294  {
295  const label cell = fc[face];
296  if (cellColour[cell] != patchToColour[inti])
297  {
298  cellColour[cell] = patchToColour[inti];
299  front.append(cell);
300  }
301  }
302  }
303  }
304 
305 
306  // Walk current mesh
307  walkFront(mesh, front, cellColour);
308 
309 
310  // Any unvisited cell is still labelMax. Put into colour 0.
311  for (auto& colour : cellColour)
312  {
313  if (colour == -1)
314  {
315  colour = 0;
316  }
317  }
318 
319  if (debug)
320  {
321  Info<< typeName << " : coloured " << nCells
322  << " cells with in total " << colouri << " colours" << endl;
323  }
325  return colouri;
326 }
327 
328 
329 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
330 
331 Foam::processorColour::processorColour(const lduMesh& mesh)
332 :
333  MeshObject_type(mesh)
334 {
335  nColours_ = colour(mesh, *this);
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
340 
342 {
343  auto* ptr =
344  mesh.thisDb().objectRegistry::template getObjectPtr<processorColour>
345  (
346  processorColour::typeName
347  );
348 
349  if (!ptr)
350  {
351  ptr = new processorColour(mesh);
352 
353  //regIOobject::store(static_cast<MoveableMeshObject<lduMesh>*>(ptr));
354  regIOobject::store(ptr);
355  }
356 
357  return *ptr;
358 }
359 
360 
361 // ************************************************************************* //
virtual lduInterfacePtrsList interfaces() const =0
Return a list of pointers for each patch with only those pointing to interfaces being set...
static label cellColour(const lduMesh &mesh, labelList &cellColour, labelList &patchToColour)
Calculate (locally) per cell colour according to walking from global patches. Returns number of colou...
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
T & front()
Access first element of the list, position [0].
Definition: UListI.H:230
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void gatherList(const UList< commsStruct > &comms, UList< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:168
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1252
static void walkFront(const lduMesh &mesh, DynamicList< label > &front, labelList &cellColour)
virtual label comm() const =0
Return communicator used for parallel communication.
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:1086
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...
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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:1077
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
dynamicFvMesh & mesh
const lduMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:255
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:752
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:392
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
static const processorColour & New(const lduMesh &mesh)
Should use the MeshObject provided one but that needs a.
static label colour(const lduMesh &mesh, labelList &procColour)
Calculate processor colouring from processor connectivity. Sets colour per processor such that no nei...
Define the processor-processor connection table by walking a list of patches and detecting the proces...
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:769
Colouring processors such that no neighbours have the same colour.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
const polyBoundaryMesh & patches
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
List< label > labelList
A List of labels.
Definition: List.H:62
Namespace for OpenFOAM.