domainDecompositionMesh.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) 2011-2016 OpenFOAM Foundation
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 InClass
27  domainDecomposition
28 
29 Description
30  Private member of domainDecomposition.
31  Decomposes the mesh into bits
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "domainDecomposition.H"
36 #include "IOstreams.H"
37 #include "bitSet.H"
38 #include "cyclicPolyPatch.H"
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 void Foam::domainDecomposition::addInterProcFace
43 (
44  const label facei,
45  const label ownerProc,
46  const label nbrProc,
47 
48  List<Map<label>>& nbrToInterPatch,
49  List<DynamicList<DynamicList<label>>>& interPatchFaces
50 ) const
51 {
52  // Introduce turning index only for internal faces (are duplicated).
53  const label ownerIndex = facei+1;
54  const label nbrIndex = -(facei+1);
55 
56  const auto patchiter = nbrToInterPatch[ownerProc].cfind(nbrProc);
57 
58  if (patchiter.good())
59  {
60  // Existing interproc patch. Add to both sides.
61  const label toNbrProcPatchi = *patchiter;
62  interPatchFaces[ownerProc][toNbrProcPatchi].append(ownerIndex);
63 
64  if (isInternalFace(facei))
65  {
66  label toOwnerProcPatchi = nbrToInterPatch[nbrProc][ownerProc];
67  interPatchFaces[nbrProc][toOwnerProcPatchi].append(nbrIndex);
68  }
69  }
70  else
71  {
72  // Create new interproc patches.
73  const label toNbrProcPatchi = nbrToInterPatch[ownerProc].size();
74  nbrToInterPatch[ownerProc].insert(nbrProc, toNbrProcPatchi);
75 
76  DynamicList<label> oneFace;
77  oneFace.append(ownerIndex);
78  interPatchFaces[ownerProc].append(oneFace);
79 
80  if (isInternalFace(facei))
81  {
82  label toOwnerProcPatchi = nbrToInterPatch[nbrProc].size();
83  nbrToInterPatch[nbrProc].insert(ownerProc, toOwnerProcPatchi);
84  oneFace.clear();
85  oneFace.append(nbrIndex);
86  interPatchFaces[nbrProc].append(oneFace);
87  }
88  }
89 }
90 
91 
93 {
94  // Decide which cell goes to which processor
95  distributeCells();
96 
97  // Distribute the cells according to the given processor label
98 
99  // calculate the addressing information for the original mesh
100  Info<< "\nCalculating original mesh data" << endl;
101 
102  // set references to the original mesh
103  const polyBoundaryMesh& patches = boundaryMesh();
104  const faceList& fcs = faces();
105  const labelList& owner = faceOwner();
106  const labelList& neighbour = faceNeighbour();
107 
108  // loop through the list of processor labels for the cell and add the
109  // cell shape to the list of cells for the appropriate processor
110 
111  Info<< "\nDistributing cells to processors" << endl;
112 
113  // Cells per processor
114  procCellAddressing_ = invertOneToMany(nProcs_, cellToProc_);
115 
116  Info<< "\nDistributing faces to processors" << endl;
117 
118  // Loop through all internal faces and decide which processor they belong to
119  // First visit all internal faces. If cells at both sides belong to the
120  // same processor, the face is an internal face. If they are different,
121  // it belongs to both processors.
122 
123  procFaceAddressing_.setSize(nProcs_);
124 
125  // Internal faces
126  forAll(neighbour, facei)
127  {
128  if (cellToProc_[owner[facei]] == cellToProc_[neighbour[facei]])
129  {
130  // Face internal to processor. Notice no turning index.
131  procFaceAddressing_[cellToProc_[owner[facei]]].append(facei+1);
132  }
133  }
134 
135  // for all processors, set the size of start index and patch size
136  // lists to the number of patches in the mesh
137  forAll(procPatchSize_, proci)
138  {
139  procPatchSize_[proci].setSize(patches.size());
140  procPatchStartIndex_[proci].setSize(patches.size());
141  }
142 
143  forAll(patches, patchi)
144  {
145  // Reset size and start index for all processors
146  forAll(procPatchSize_, proci)
147  {
148  procPatchSize_[proci][patchi] = 0;
149  procPatchStartIndex_[proci][patchi] =
150  procFaceAddressing_[proci].size();
151  }
152 
153  const label patchStart = patches[patchi].start();
154 
155  if (!isA<cyclicPolyPatch>(patches[patchi]))
156  {
157  // Normal patch. Add faces to processor where the cell
158  // next to the face lives
159 
160  const labelUList& patchFaceCells =
161  patches[patchi].faceCells();
162 
163  forAll(patchFaceCells, facei)
164  {
165  const label curProc = cellToProc_[patchFaceCells[facei]];
166 
167  // add the face without turning index
168  procFaceAddressing_[curProc].append(patchStart+facei+1);
169 
170  // increment the number of faces for this patch
171  procPatchSize_[curProc][patchi]++;
172  }
173  }
174  else
175  {
176  const cyclicPolyPatch& pp = refCast<const cyclicPolyPatch>
177  (
178  patches[patchi]
179  );
180  // cyclic: check opposite side on this processor
181  const labelUList& patchFaceCells = pp.faceCells();
182 
183  const labelUList& nbrPatchFaceCells =
184  pp.neighbPatch().faceCells();
185 
186  forAll(patchFaceCells, facei)
187  {
188  const label curProc = cellToProc_[patchFaceCells[facei]];
189  const label nbrProc = cellToProc_[nbrPatchFaceCells[facei]];
190  if (curProc == nbrProc)
191  {
192  // add the face without turning index
193  procFaceAddressing_[curProc].append(patchStart+facei+1);
194  // increment the number of faces for this patch
195  procPatchSize_[curProc][patchi]++;
196  }
197  }
198  }
199  }
200 
201 
202  // Done internal bits of the new mesh and the ordinary patches.
203 
204 
205  // Per processor, from neighbour processor to the inter-processor patch
206  // that communicates with that neighbour
207  List<Map<label>> procNbrToInterPatch(nProcs_);
208 
209  // Per processor the faces per inter-processor patch
210  List<DynamicList<DynamicList<label>>> interPatchFaces(nProcs_);
211 
212  // Processor boundaries from internal faces
213  forAll(neighbour, facei)
214  {
215  label ownerProc = cellToProc_[owner[facei]];
216  label nbrProc = cellToProc_[neighbour[facei]];
217 
218  if (ownerProc != nbrProc)
219  {
220  // inter - processor patch face found.
221  addInterProcFace
222  (
223  facei,
224  ownerProc,
225  nbrProc,
226 
227  procNbrToInterPatch,
228  interPatchFaces
229  );
230  }
231  }
232 
233  // Add the proper processor faces to the sub information. For faces
234  // originating from internal faces this is always -1.
235  List<labelListList> subPatchIDs(nProcs_);
236  List<labelListList> subPatchStarts(nProcs_);
237  forAll(interPatchFaces, proci)
238  {
239  label nInterfaces = interPatchFaces[proci].size();
240 
241  subPatchIDs[proci].setSize(nInterfaces, labelList(1, label(-1)));
242  subPatchStarts[proci].setSize(nInterfaces, labelList(1, Zero));
243  }
244 
245 
246  // Special handling needed for the case that multiple processor cyclic
247  // patches are created on each local processor domain, e.g. if a 3x3 case
248  // is decomposed using the decomposition:
249  //
250  // | 1 | 0 | 2 |
251  // cyclic left | 2 | 0 | 1 | cyclic right
252  // | 2 | 0 | 1 |
253  //
254  // - processors 1 and 2 will both have pieces of both cyclic left- and
255  // right sub-patches present
256  // - the interface patch faces are stored in a single list, where each
257  // sub-patch is referenced into the list using a patch start index and
258  // size
259  // - if the patches are in order (in the boundary file) of left, right
260  // - processor 1 will send: left, right
261  // - processor 1 will need to receive in reverse order: right, left
262  // - similarly for processor 2
263  // - the sub-patches are therefore generated in 4 passes of the patch lists
264  // 1. add faces from owner patch where local proc i < nbr proc i
265  // 2. add faces from nbr patch where local proc i < nbr proc i
266  // 3. add faces from owner patch where local proc i > nbr proc i
267  // 4. add faces from nbr patch where local proc i > nbr proc i
268 
269  processInterCyclics
270  (
271  patches,
272  interPatchFaces,
273  procNbrToInterPatch,
274  subPatchIDs,
275  subPatchStarts,
276  true,
277  lessOp<label>()
278  );
279 
280  processInterCyclics
281  (
282  patches,
283  interPatchFaces,
284  procNbrToInterPatch,
285  subPatchIDs,
286  subPatchStarts,
287  false,
288  lessOp<label>()
289  );
290 
291  processInterCyclics
292  (
293  patches,
294  interPatchFaces,
295  procNbrToInterPatch,
296  subPatchIDs,
297  subPatchStarts,
298  false,
299  greaterOp<label>()
300  );
301 
302  processInterCyclics
303  (
304  patches,
305  interPatchFaces,
306  procNbrToInterPatch,
307  subPatchIDs,
308  subPatchStarts,
309  true,
310  greaterOp<label>()
311  );
312 
313 
314  // Sort inter-proc patch by neighbour
315  labelList order;
316  forAll(procNbrToInterPatch, proci)
317  {
318  label nInterfaces = procNbrToInterPatch[proci].size();
319 
320  procNeighbourProcessors_[proci].setSize(nInterfaces);
321  procProcessorPatchSize_[proci].setSize(nInterfaces);
322  procProcessorPatchStartIndex_[proci].setSize(nInterfaces);
323  procProcessorPatchSubPatchIDs_[proci].setSize(nInterfaces);
324  procProcessorPatchSubPatchStarts_[proci].setSize(nInterfaces);
325 
326  //Info<< "Processor " << proci << endl;
327 
328  // Get sorted neighbour processors
329  const Map<label>& curNbrToInterPatch = procNbrToInterPatch[proci];
330  labelList nbrs = curNbrToInterPatch.toc();
331 
332  sortedOrder(nbrs, order);
333 
334  DynamicList<DynamicList<label>>& curInterPatchFaces =
335  interPatchFaces[proci];
336 
337  forAll(nbrs, i)
338  {
339  const label nbrProc = nbrs[i];
340  const label interPatch = curNbrToInterPatch[nbrProc];
341 
342  procNeighbourProcessors_[proci][i] = nbrProc;
343  procProcessorPatchSize_[proci][i] =
344  curInterPatchFaces[interPatch].size();
345  procProcessorPatchStartIndex_[proci][i] =
346  procFaceAddressing_[proci].size();
347 
348  // Add size as last element to substarts and transfer
349  subPatchStarts[proci][interPatch].append
350  (
351  curInterPatchFaces[interPatch].size()
352  );
353  procProcessorPatchSubPatchIDs_[proci][i].transfer
354  (
355  subPatchIDs[proci][interPatch]
356  );
357  procProcessorPatchSubPatchStarts_[proci][i].transfer
358  (
359  subPatchStarts[proci][interPatch]
360  );
361 
362  //Info<< " nbr:" << nbrProc << endl;
363  //Info<< " interpatch:" << interPatch << endl;
364  //Info<< " size:" << procProcessorPatchSize_[proci][i] << endl;
365  //Info<< " start:" << procProcessorPatchStartIndex_[proci][i]
366  // << endl;
367  //Info<< " subPatches:"
368  // << procProcessorPatchSubPatchIDs_[proci][i]
369  // << endl;
370  //Info<< " subStarts:"
371  // << procProcessorPatchSubPatchStarts_[proci][i] << endl;
372 
373  // And add all the face labels for interPatch
374  DynamicList<label>& interPatchFaces =
375  curInterPatchFaces[interPatch];
376 
377  forAll(interPatchFaces, j)
378  {
379  procFaceAddressing_[proci].append(interPatchFaces[j]);
380  }
381  interPatchFaces.clearStorage();
382  }
383  curInterPatchFaces.clearStorage();
384  procFaceAddressing_[proci].shrink();
385  }
386 
387 
390 // forAll(procPatchStartIndex_, proci)
391 // {
392 // Info<< "Processor:" << proci << endl;
393 //
394 // Info<< " total faces:" << procFaceAddressing_[proci].size()
395 // << endl;
396 //
397 // const labelList& curProcPatchStartIndex = procPatchStartIndex_[proci];
398 //
399 // forAll(curProcPatchStartIndex, patchi)
400 // {
401 // Info<< " patch:" << patchi
402 // << "\tstart:" << curProcPatchStartIndex[patchi]
403 // << "\tsize:" << procPatchSize_[proci][patchi]
404 // << endl;
405 // }
406 // }
407 // Info<< endl;
408 //
409 // forAll(procNeighbourProcessors_, proci)
410 // {
411 // Info<< "Processor " << proci << endl;
412 //
413 // forAll(procNeighbourProcessors_[proci], i)
414 // {
415 // Info<< " nbr:" << procNeighbourProcessors_[proci][i] << endl;
416 // Info<< " size:" << procProcessorPatchSize_[proci][i] << endl;
417 // Info<< " start:" << procProcessorPatchStartIndex_[proci][i]
418 // << endl;
419 // }
420 // }
421 // Info<< endl;
422 //
423 // forAll(procFaceAddressing_, proci)
424 // {
425 // Info<< "Processor:" << proci << endl;
426 //
427 // Info<< " faces:" << procFaceAddressing_[proci] << endl;
428 // }
429 
430 
431  Info<< "\nDistributing points to processors" << endl;
432  // For every processor, loop through the list of faces for the processor.
433  // For every face, loop through the list of points and mark the point as
434  // used for the processor. Collect the list of used points for the
435  // processor.
436 
437  forAll(procPointAddressing_, proci)
438  {
439  bitSet pointsInUse(nPoints(), false);
440 
441  // For each of the faces used
442  for (const label facei : procFaceAddressing_[proci])
443  {
444  // Because of the turning index, some face labels may be -ve
445  const labelList& facePoints = fcs[mag(facei) - 1];
446 
447  // Mark the face points as being used
448  pointsInUse.set(facePoints);
449  }
450 
451  procPointAddressing_[proci] = pointsInUse.sortedToc();
452  }
453 }
454 
455 // ************************************************************************* //
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:493
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:126
label nPoints
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const polyBoundaryMesh & patches
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
void decomposeMesh()
Decompose mesh.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127