Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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"
31
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
37 defineTypeNameAndDebug(processorColour, 0);
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 {
58 const lduInterfacePtrsList patches = mesh.interfaces();
59
60 auto& procToProcs = procNeighbours[Pstream::myProcNo(mesh.comm())];
61 label n = 0;
62 forAll(patches, patchi)
63 {
64 if (patches.set(patchi))
65 {
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 Pout<< 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 = uPtr[facei];
206 if (cellColour[nbr] == -1)
207 {
208 cellColour[nbr] = colouri;
209 newFront.append(nbr);
210 }
211 }
212 }
213 {
214 const label fStart = losortStartAddrPtr[celli];
215 const label fEnd = losortStartAddrPtr[celli + 1];
216
217 for (label i=fStart; i<fEnd; i++)
218 {
219 label facei = losortAddrPtr[i];
220 const label nbr = lPtr[facei];
221 if (cellColour[nbr] == -1)
222 {
223 cellColour[nbr] = colouri;
224 newFront.append(nbr);
225 }
226 }
227 }
228 }
229
230 if (newFront.empty())
231 {
232 break;
234
235 front.transfer(newFront);
236 }
237}
238
239
241(
242 const lduMesh& mesh,
243 labelList& cellColour,
244 labelList& patchToColour
245)
246{
247 // Colour with the (min) coupled (global) patch. Compact to have colour 0
248 // if a single region. Problematic if same patch on multiple disconnected
249 // regions.
250
251 const lduAddressing& addr = mesh.lduAddr();
253
254 const label nCells = addr.size();
255
256 patchToColour.resize_nocopy(patches.size());
257 patchToColour = -1;
258 cellColour.resize_nocopy(nCells);
259 cellColour = -1;
260
261
262 label colouri = 0;
263
264 // Starting front from patch faceCells
265 DynamicList<label> front;
266 forAll(patches, inti)
267 {
268 if
269 (
270 patches.set(inti)
272 )
273 {
274 // 'global' interface. Seed faceCells with patch index
275
276 if (patchToColour[inti] == -1)
277 {
278 patchToColour[inti] = colouri++;
279 }
280
281 const auto& fc = patches[inti].faceCells();
282 forAll(fc, face)
283 {
284 const label cell = fc[face];
285 if (cellColour[cell] != patchToColour[inti])
286 {
287 cellColour[cell] = patchToColour[inti];
288 front.append(cell);
289 }
290 }
291 }
292 }
293
294
295 // Walk current mesh
296 walkFront(mesh, front, cellColour);
297
298
299 // Any unvisited cell is still labelMax. Put into colour 0.
300 for (auto& colour : cellColour)
301 {
302 if (colour == -1)
303 {
304 colour = 0;
305 }
306 }
307
308 if (debug)
309 {
310 Info<< typeName << " : coloured " << nCells
311 << " cells with in total " << colouri << " colours" << endl;
312 }
314 return colouri;
315}
316
317
318// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
319
320Foam::processorColour::processorColour(const lduMesh& mesh)
321:
322 MeshObject_type(mesh)
324 nColours_ = colour(mesh, *this);
325}
326
327
328// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
329
331{
332 auto* ptr =
333 mesh.thisDb().objectRegistry::template getObjectPtr<processorColour>
334 (
335 processorColour::typeName
336 );
337
338 if (!ptr)
339 {
340 ptr = new processorColour(mesh);
341
342 //regIOobject::store(static_cast<MoveableMeshObject<lduMesh>*>(ptr));
344 }
345
346 return *ptr;
347}
348
349
350// ************************************************************************* //
label n
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
bool found(const T &val, label pos=0) const
Same as contains().
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
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:171
A List with indirect addressing. Like IndirectList but does not store addressing.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
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 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
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition fvMesh.C:724
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition fvMesh.C:707
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const labelUList & ownerStartAddr() const
Return owner start addressing.
const labelUList & losortStartAddr() const
Return losort start addressing.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
label size() const noexcept
Return number of equations.
const labelUList & losortAddr() const
Return losort addressing.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
Colouring processors such that no neighbours have the same colour.
static const processorColour & New(const lduMesh &mesh)
Should use the MeshObject provided one but that needs a.
label nColours_
Max number of colours.
static void walkFront(const lduMesh &mesh, DynamicList< label > &front, labelList &cellColour)
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...
static label colour(const lduMesh &mesh, labelList &procColour)
Calculate processor colouring from processor connectivity. Sets colour per processor such that no nei...
bool store()
Register object with its registry and transfer ownership to the registry.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Define the processor-processor connection table by walking a list of patches and detecting the proces...
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299