meshToMeshParallelOps.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2022 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 
29 #include "meshToMesh.H"
30 #include "OFstream.H"
31 #include "Time.H"
32 #include "globalIndex.H"
33 #include "mergePoints.H"
34 #include "processorPolyPatch.H"
35 #include "SubField.H"
36 #include "AABBTree.H"
37 #include "cellBox.H"
38 
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40 
41 Foam::label Foam::meshToMesh::calcDistribution
42 (
43  const polyMesh& src,
44  const polyMesh& tgt
45 ) const
46 {
47  label proci = 0;
48 
49  if (Pstream::parRun())
50  {
51  const bitSet hasMesh
52  (
53  UPstream::listGatherValues<bool>
54  (
55  src.nCells() > 0 || tgt.nCells() > 0
56  )
57  );
58 
59  const auto nHaveMesh = hasMesh.count();
60 
61  if (nHaveMesh == 1)
62  {
63  proci = hasMesh.find_first();
65  << "Meshes local to processor" << proci << endl;
66  }
67  else if (nHaveMesh > 1)
68  {
69  proci = -1;
71  << "Meshes split across multiple processors" << endl;
72  }
73 
74  Pstream::broadcast(proci);
75  }
76 
77  return proci;
78 }
79 
80 
81 Foam::label Foam::meshToMesh::calcOverlappingProcs
82 (
83  const List<treeBoundBoxList>& procBb,
84  const boundBox& bb,
85  boolList& overlaps
86 ) const
87 {
88  overlaps = false;
89 
90  label nOverlaps = 0;
91 
92  forAll(procBb, proci)
93  {
94  const treeBoundBoxList& bbp = procBb[proci];
95 
96  for (const treeBoundBox& b : bbp)
97  {
98  if (b.overlaps(bb))
99  {
100  overlaps[proci] = true;
101  ++nOverlaps;
102  break;
103  }
104  }
105  }
106 
107  return nOverlaps;
108 }
109 
110 
111 Foam::autoPtr<Foam::mapDistribute> Foam::meshToMesh::calcProcMap
112 (
113  const polyMesh& src,
114  const polyMesh& tgt
115 ) const
116 {
117  switch (procMapMethod_)
118  {
119  case procMapMethod::pmLOD:
120  {
121  Info<< "meshToMesh: Using processorLOD method" << endl;
122 
123  // Create processor map of overlapping faces. This map gets
124  // (possibly remote) cells from the tgt mesh such that they
125  // (together) cover all of the src mesh
126  const label nGlobalSrcCells = src.globalData().nTotalCells();
127  // Note: minimum content has to be > 1 since otherwise
128  // would keep on splitting. 16 is fairly random choice.
129  const label cellsPerBox = max(16, 0.001*nGlobalSrcCells);
130  typename processorLODs::cellBox boxLOD
131  (
132  src.cells(),
133  src.faces(),
134  src.points(),
135  tgt.cells(),
136  tgt.faces(),
137  tgt.points(),
138  cellsPerBox,
139  src.nCells()
140  );
141 
142  return boxLOD.map();
143  break;
144  }
145  default:
146  {
147  Info<< "meshToMesh: Using AABBTree method" << endl;
148 
149  // get decomposition of cells on src mesh
150  List<treeBoundBoxList> procBb(Pstream::nProcs());
151 
152  if (src.nCells() > 0)
153  {
154  procBb[Pstream::myProcNo()] = AABBTree<labelList>
155  (
156  src.cellPoints(),
157  src.points(),
158  false
159  ).boundBoxes();
160  }
161  else
162  {
163  procBb[Pstream::myProcNo()] = treeBoundBoxList();
164  }
165 
166  Pstream::allGatherList(procBb);
167 
168 
169  if (debug)
170  {
172  << "Determining extent of src mesh per processor:" << nl
173  << "\tproc\tbb" << endl;
174  forAll(procBb, proci)
175  {
176  Info<< '\t' << proci << '\t' << procBb[proci] << endl;
177  }
178  }
179 
180 
181  // Determine which cells of tgt mesh overlaps src mesh per proc
182  labelListList sendMap(Pstream::nProcs());
183 
184  {
185  // per processor indices into all segments to send
186  List<DynamicList<label>> dynSendMap(Pstream::nProcs());
187  label iniSize = floor(tgt.nCells()/Pstream::nProcs());
188 
189  forAll(dynSendMap, proci)
190  {
191  dynSendMap[proci].setCapacity(iniSize);
192  }
193 
194  // work array - whether src processor bb overlaps the tgt cell
195  // bounds
196  boolList procBbOverlaps(Pstream::nProcs());
197 
198  for (label celli = 0; celli < tgt.nCells(); ++celli)
199  {
200  // Bounding box of tgt cell
201  boundBox cellBb(tgt.cellBb(celli));
202 
203  // find the overlapping tgt cells on each src processor
204  (void)calcOverlappingProcs(procBb, cellBb, procBbOverlaps);
205 
206  forAll(procBbOverlaps, proci)
207  {
208  if (procBbOverlaps[proci])
209  {
210  dynSendMap[proci].append(celli);
211  }
212  }
213  }
214 
215  // convert dynamicList to labelList
216  forAll(sendMap, proci)
217  {
218  sendMap[proci].transfer(dynSendMap[proci]);
219  }
220 
221  // debug printing
222  if (debug)
223  {
224  Pout<< "Of my " << tgt.nCells()
225  << " target cells I need to send to:" << nl
226  << "\tproc\tcells" << endl;
227  forAll(sendMap, proci)
228  {
229  Pout<< '\t' << proci << '\t'
230  << sendMap[proci].size() << endl;
231  }
232  }
233  }
234 
235 
236  // send over how many tgt cells I need to receive from each
237  // processor
238  labelListList sendSizes(Pstream::nProcs());
239  sendSizes[Pstream::myProcNo()].resize(Pstream::nProcs());
240  forAll(sendMap, proci)
241  {
242  sendSizes[Pstream::myProcNo()][proci] = sendMap[proci].size();
243  }
244  Pstream::allGatherList(sendSizes);
245 
246 
247  // determine order of receiving
248  labelListList constructMap(Pstream::nProcs());
249 
250  label segmentI = 0;
251  forAll(constructMap, proci)
252  {
253  // what I need to receive is what other processor is sending
254  // to me
255  label nRecv = sendSizes[proci][Pstream::myProcNo()];
256  constructMap[proci].setSize(nRecv);
257 
258  for (label i = 0; i < nRecv; i++)
259  {
260  constructMap[proci][i] = segmentI++;
261  }
262  }
263 
265  (
266  segmentI, // size after construction
267  std::move(sendMap),
268  std::move(constructMap)
269  );
270 
271  break;
272  }
273  }
274 }
275 
276 
277 void Foam::meshToMesh::distributeCells
278 (
279  const mapDistribute& map,
280  const polyMesh& tgtMesh,
281  const globalIndex& globalI,
282  List<pointField>& points,
283  List<label>& nInternalFaces,
284  List<faceList>& faces,
285  List<labelList>& faceOwner,
286  List<labelList>& faceNeighbour,
287  List<labelList>& cellIDs,
288  List<labelList>& nbrProcIDs,
289  List<labelList>& procLocalFaceIDs
290 ) const
291 {
293  nInternalFaces.setSize(Pstream::nProcs(), 0);
294  faces.setSize(Pstream::nProcs());
295  faceOwner.setSize(Pstream::nProcs());
296  faceNeighbour.setSize(Pstream::nProcs());
298 
299  nbrProcIDs.setSize(Pstream::nProcs());;
300  procLocalFaceIDs.setSize(Pstream::nProcs());;
301 
302 
303  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
304 
305  for (const int domain : Pstream::allProcs())
306  {
307  const labelList& sendElems = map.subMap()[domain];
308 
309  if (sendElems.size())
310  {
311  // reverse cell map
312  labelList reverseCellMap(tgtMesh.nCells(), -1);
313  forAll(sendElems, subCelli)
314  {
315  reverseCellMap[sendElems[subCelli]] = subCelli;
316  }
317 
318  DynamicList<face> subFaces(tgtMesh.nFaces());
319  DynamicList<label> subFaceOwner(tgtMesh.nFaces());
320  DynamicList<label> subFaceNeighbour(tgtMesh.nFaces());
321 
322  DynamicList<label> subNbrProcIDs(tgtMesh.nFaces());
323  DynamicList<label> subProcLocalFaceIDs(tgtMesh.nFaces());
324 
325  label nInternal = 0;
326 
327  // internal faces
328  forAll(tgtMesh.faceNeighbour(), facei)
329  {
330  label own = tgtMesh.faceOwner()[facei];
331  label nbr = tgtMesh.faceNeighbour()[facei];
332  label subOwn = reverseCellMap[own];
333  label subNbr = reverseCellMap[nbr];
334 
335  if (subOwn != -1 && subNbr != -1)
336  {
337  nInternal++;
338 
339  if (subOwn < subNbr)
340  {
341  subFaces.append(tgtMesh.faces()[facei]);
342  subFaceOwner.append(subOwn);
343  subFaceNeighbour.append(subNbr);
344  subNbrProcIDs.append(-1);
345  subProcLocalFaceIDs.append(-1);
346  }
347  else
348  {
349  subFaces.append(tgtMesh.faces()[facei].reverseFace());
350  subFaceOwner.append(subNbr);
351  subFaceNeighbour.append(subOwn);
352  subNbrProcIDs.append(-1);
353  subProcLocalFaceIDs.append(-1);
354  }
355  }
356  }
357 
358  // boundary faces for new region
359  forAll(tgtMesh.faceNeighbour(), facei)
360  {
361  label own = tgtMesh.faceOwner()[facei];
362  label nbr = tgtMesh.faceNeighbour()[facei];
363  label subOwn = reverseCellMap[own];
364  label subNbr = reverseCellMap[nbr];
365 
366  if (subOwn != -1 && subNbr == -1)
367  {
368  subFaces.append(tgtMesh.faces()[facei]);
369  subFaceOwner.append(subOwn);
370  subFaceNeighbour.append(subNbr);
371  subNbrProcIDs.append(-1);
372  subProcLocalFaceIDs.append(-1);
373  }
374  else if (subOwn == -1 && subNbr != -1)
375  {
376  subFaces.append(tgtMesh.faces()[facei].reverseFace());
377  subFaceOwner.append(subNbr);
378  subFaceNeighbour.append(subOwn);
379  subNbrProcIDs.append(-1);
380  subProcLocalFaceIDs.append(-1);
381  }
382  }
383 
384  // boundary faces of existing region
385  forAll(tgtMesh.boundaryMesh(), patchi)
386  {
387  const polyPatch& pp = tgtMesh.boundaryMesh()[patchi];
388  const auto* procPatch = isA<processorPolyPatch>(pp);
389 
390  // Store info for faces on processor patches
391  const label nbrProci =
392  (procPatch ? procPatch->neighbProcNo() : -1);
393 
394  forAll(pp, i)
395  {
396  label facei = pp.start() + i;
397  label own = tgtMesh.faceOwner()[facei];
398 
399  if (reverseCellMap[own] != -1)
400  {
401  subFaces.append(tgtMesh.faces()[facei]);
402  subFaceOwner.append(reverseCellMap[own]);
403  subFaceNeighbour.append(-1);
404  subNbrProcIDs.append(nbrProci);
405  subProcLocalFaceIDs.append(i);
406  }
407  }
408  }
409 
410  // reverse point map
411  labelList reversePointMap(tgtMesh.nPoints(), -1);
412  DynamicList<point> subPoints(tgtMesh.nPoints());
413  forAll(subFaces, subFacei)
414  {
415  face& f = subFaces[subFacei];
416  forAll(f, fp)
417  {
418  label pointi = f[fp];
419  if (reversePointMap[pointi] == -1)
420  {
421  reversePointMap[pointi] = subPoints.size();
422  subPoints.append(tgtMesh.points()[pointi]);
423  }
424 
425  f[fp] = reversePointMap[pointi];
426  }
427  }
428 
429  // tgt cells into global numbering
430  labelList globalElems(globalI.toGlobal(sendElems));
431 
432  if (debug > 1)
433  {
434  forAll(sendElems, i)
435  {
436  Pout<< "tgtProc:" << Pstream::myProcNo()
437  << " sending tgt cell " << sendElems[i]
438  << "[" << globalElems[i] << "]"
439  << " to srcProc " << domain << endl;
440  }
441  }
442 
443  // pass data
444  if (domain == Pstream::myProcNo())
445  {
446  // allocate my own data
447  points[Pstream::myProcNo()] = subPoints;
448  nInternalFaces[Pstream::myProcNo()] = nInternal;
449  faces[Pstream::myProcNo()] = subFaces;
450  faceOwner[Pstream::myProcNo()] = subFaceOwner;
451  faceNeighbour[Pstream::myProcNo()] = subFaceNeighbour;
452  cellIDs[Pstream::myProcNo()] = globalElems;
453  nbrProcIDs[Pstream::myProcNo()] = subNbrProcIDs;
454  procLocalFaceIDs[Pstream::myProcNo()] = subProcLocalFaceIDs;
455  }
456  else
457  {
458  // send data to other processor domains
459  UOPstream toDomain(domain, pBufs);
460 
461  toDomain
462  << subPoints
463  << nInternal
464  << subFaces
465  << subFaceOwner
466  << subFaceNeighbour
467  << globalElems
468  << subNbrProcIDs
469  << subProcLocalFaceIDs;
470  }
471  }
472  }
473 
474  // Start receiving
475  pBufs.finishedSends();
476 
477  // Consume
478  for (const int domain : Pstream::allProcs())
479  {
480  const labelList& recvElems = map.constructMap()[domain];
481 
482  if (domain != Pstream::myProcNo() && recvElems.size())
483  {
484  UIPstream str(domain, pBufs);
485 
486  str >> points[domain]
487  >> nInternalFaces[domain]
488  >> faces[domain]
489  >> faceOwner[domain]
490  >> faceNeighbour[domain]
491  >> cellIDs[domain]
492  >> nbrProcIDs[domain]
493  >> procLocalFaceIDs[domain];
494  }
495 
496  if (debug)
497  {
498  Pout<< "Target mesh send sizes[" << domain << "]"
499  << ": points="<< points[domain].size()
500  << ", faces=" << faces[domain].size()
501  << ", nInternalFaces=" << nInternalFaces[domain]
502  << ", faceOwn=" << faceOwner[domain].size()
503  << ", faceNbr=" << faceNeighbour[domain].size()
504  << ", cellIDs=" << cellIDs[domain].size() << endl;
505  }
506  }
507 }
508 
509 
510 void Foam::meshToMesh::distributeAndMergeCells
511 (
512  const mapDistribute& map,
513  const polyMesh& tgt,
514  const globalIndex& globalI,
515  pointField& tgtPoints,
516  faceList& tgtFaces,
517  labelList& tgtFaceOwners,
518  labelList& tgtFaceNeighbours,
519  labelList& tgtCellIDs
520 ) const
521 {
522  // Exchange per-processor data
523  List<pointField> allPoints;
524  List<label> allNInternalFaces;
525  List<faceList> allFaces;
526  List<labelList> allFaceOwners;
527  List<labelList> allFaceNeighbours;
528  List<labelList> allTgtCellIDs;
529 
530  // Per target mesh face the neighbouring proc and index in
531  // processor patch (all -1 for normal boundary face)
532  List<labelList> allNbrProcIDs;
533  List<labelList> allProcLocalFaceIDs;
534 
535  distributeCells
536  (
537  map,
538  tgt,
539  globalI,
540  allPoints,
541  allNInternalFaces,
542  allFaces,
543  allFaceOwners,
544  allFaceNeighbours,
545  allTgtCellIDs,
546  allNbrProcIDs,
547  allProcLocalFaceIDs
548  );
549 
550  // Convert lists into format that can be used to generate a valid polyMesh
551  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
552  //
553  // Points and cells are collected into single flat lists:
554  // - i.e. proc0, proc1 ... procN
555  //
556  // Faces need to be sorted after collection to that internal faces are
557  // contiguous, followed by all boundary faces
558  //
559  // Processor patch faces between included cells on neighbouring processors
560  // are converted into internal faces
561  //
562  // Face list structure:
563  // - Per processor:
564  // - internal faces
565  // - processor faces that have been converted into internal faces
566  // - Followed by all boundary faces
567  // - from 'normal' boundary faces
568  // - from singularly-sided processor patch faces
569 
570 
571  // Number of internal+coupled faces
572  labelList allNIntCoupledFaces(allNInternalFaces);
573 
574  // Starting offset for points
575  label nPoints = 0;
576  labelList pointOffset(Pstream::nProcs(), Zero);
577  forAll(allPoints, proci)
578  {
579  pointOffset[proci] = nPoints;
580  nPoints += allPoints[proci].size();
581  }
582 
583  // Starting offset for cells
584  label nCells = 0;
585  labelList cellOffset(Pstream::nProcs(), Zero);
586  forAll(allTgtCellIDs, proci)
587  {
588  cellOffset[proci] = nCells;
589  nCells += allTgtCellIDs[proci].size();
590  }
591 
592  // Count any coupled faces
593  typedef FixedList<label, 3> label3;
594  typedef HashTable<label, label3> procCoupleInfo;
595  procCoupleInfo procFaceToGlobalCell;
596 
597  forAll(allNbrProcIDs, proci)
598  {
599  const labelList& nbrProci = allNbrProcIDs[proci];
600  const labelList& localFacei = allProcLocalFaceIDs[proci];
601 
602  forAll(nbrProci, i)
603  {
604  if (nbrProci[i] != -1 && localFacei[i] != -1)
605  {
606  label3 key;
607  key[0] = min(proci, nbrProci[i]);
608  key[1] = max(proci, nbrProci[i]);
609  key[2] = localFacei[i];
610 
611  const auto fnd = procFaceToGlobalCell.cfind(key);
612 
613  if (!fnd.found())
614  {
615  procFaceToGlobalCell.insert(key, -1);
616  }
617  else
618  {
619  if (debug > 1)
620  {
621  Pout<< "Additional internal face between procs:"
622  << key[0] << " and " << key[1]
623  << " across local face " << key[2] << endl;
624  }
625 
626  allNIntCoupledFaces[proci]++;
627  }
628  }
629  }
630  }
631 
632 
633  // Starting offset for internal faces
634  label nIntFaces = 0;
635  label nFacesTotal = 0;
636  labelList internalFaceOffset(Pstream::nProcs(), Zero);
637  forAll(allNIntCoupledFaces, proci)
638  {
639  label nCoupledFaces =
640  allNIntCoupledFaces[proci] - allNInternalFaces[proci];
641 
642  internalFaceOffset[proci] = nIntFaces;
643  nIntFaces += allNIntCoupledFaces[proci];
644  nFacesTotal += allFaceOwners[proci].size() - nCoupledFaces;
645  }
646 
647  tgtPoints.setSize(nPoints);
648  tgtFaces.setSize(nFacesTotal);
649  tgtFaceOwners.setSize(nFacesTotal);
650  tgtFaceNeighbours.setSize(nFacesTotal);
651  tgtCellIDs.setSize(nCells);
652 
653  // Insert points
654  forAll(allPoints, proci)
655  {
656  const pointField& pts = allPoints[proci];
657  SubList<point>(tgtPoints, pts.size(), pointOffset[proci]) = pts;
658  }
659 
660  // Insert cellIDs
661  forAll(allTgtCellIDs, proci)
662  {
663  const labelList& cellIDs = allTgtCellIDs[proci];
664  SubList<label>(tgtCellIDs, cellIDs.size(), cellOffset[proci]) = cellIDs;
665  }
666 
667 
668  // Insert internal faces (from internal faces)
669  forAll(allFaces, proci)
670  {
671  const faceList& fcs = allFaces[proci];
672  const labelList& faceOs = allFaceOwners[proci];
673  const labelList& faceNs = allFaceNeighbours[proci];
674 
675  SubList<face> slice
676  (
677  tgtFaces,
678  allNInternalFaces[proci],
679  internalFaceOffset[proci]
680  );
681  slice = SubList<face>(fcs, allNInternalFaces[proci]);
682  forAll(slice, i)
683  {
684  add(slice[i], pointOffset[proci]);
685  }
686 
687  SubField<label> ownSlice
688  (
689  tgtFaceOwners,
690  allNInternalFaces[proci],
691  internalFaceOffset[proci]
692  );
693  ownSlice = SubField<label>(faceOs, allNInternalFaces[proci]);
694  add(ownSlice, cellOffset[proci]);
695 
696  SubField<label> nbrSlice
697  (
698  tgtFaceNeighbours,
699  allNInternalFaces[proci],
700  internalFaceOffset[proci]
701  );
702  nbrSlice = SubField<label>(faceNs, allNInternalFaces[proci]);
703  add(nbrSlice, cellOffset[proci]);
704 
705  internalFaceOffset[proci] += allNInternalFaces[proci];
706  }
707 
708 
709  // Insert internal faces (from coupled face-pairs)
710  forAll(allNbrProcIDs, proci)
711  {
712  const labelList& nbrProci = allNbrProcIDs[proci];
713  const labelList& localFacei = allProcLocalFaceIDs[proci];
714  const labelList& faceOs = allFaceOwners[proci];
715  const faceList& fcs = allFaces[proci];
716 
717  forAll(nbrProci, i)
718  {
719  if (nbrProci[i] != -1 && localFacei[i] != -1)
720  {
721  label3 key;
722  key[0] = min(proci, nbrProci[i]);
723  key[1] = max(proci, nbrProci[i]);
724  key[2] = localFacei[i];
725 
726  auto fnd = procFaceToGlobalCell.find(key);
727 
728  if (fnd.found())
729  {
730  label tgtFacei = fnd();
731  if (tgtFacei == -1)
732  {
733  // on first visit store the new cell on this side
734  fnd() = cellOffset[proci] + faceOs[i];
735  }
736  else
737  {
738  // get owner and neighbour in new cell numbering
739  label newOwn = cellOffset[proci] + faceOs[i];
740  label newNbr = fnd();
741  label tgtFacei = internalFaceOffset[proci]++;
742 
743  if (debug > 1)
744  {
745  Pout<< " proc " << proci
746  << "\tinserting face:" << tgtFacei
747  << " connection between owner " << newOwn
748  << " and neighbour " << newNbr
749  << endl;
750  }
751 
752  if (newOwn < newNbr)
753  {
754  // we have correct orientation
755  tgtFaces[tgtFacei] = fcs[i];
756  tgtFaceOwners[tgtFacei] = newOwn;
757  tgtFaceNeighbours[tgtFacei] = newNbr;
758  }
759  else
760  {
761  // reverse orientation
762  tgtFaces[tgtFacei] = fcs[i].reverseFace();
763  tgtFaceOwners[tgtFacei] = newNbr;
764  tgtFaceNeighbours[tgtFacei] = newOwn;
765  }
766 
767  add(tgtFaces[tgtFacei], pointOffset[proci]);
768 
769  // mark with unique value
770  fnd() = -2;
771  }
772  }
773  }
774  }
775  }
776 
777 
778  forAll(allNbrProcIDs, proci)
779  {
780  const labelList& nbrProci = allNbrProcIDs[proci];
781  const labelList& localFacei = allProcLocalFaceIDs[proci];
782  const labelList& faceOs = allFaceOwners[proci];
783  const labelList& faceNs = allFaceNeighbours[proci];
784  const faceList& fcs = allFaces[proci];
785 
786  forAll(nbrProci, i)
787  {
788  // coupled boundary face
789  if (nbrProci[i] != -1 && localFacei[i] != -1)
790  {
791  label3 key;
792  key[0] = min(proci, nbrProci[i]);
793  key[1] = max(proci, nbrProci[i]);
794  key[2] = localFacei[i];
795 
796  label tgtFacei = procFaceToGlobalCell[key];
797 
798  if (tgtFacei == -1)
799  {
801  << "Unvisited " << key
802  << abort(FatalError);
803  }
804  else if (tgtFacei != -2)
805  {
806  label newOwn = cellOffset[proci] + faceOs[i];
807  label tgtFacei = nIntFaces++;
808 
809  if (debug > 1)
810  {
811  Pout<< " proc " << proci
812  << "\tinserting boundary face:" << tgtFacei
813  << " from coupled face " << key
814  << endl;
815  }
816 
817  tgtFaces[tgtFacei] = fcs[i];
818  add(tgtFaces[tgtFacei], pointOffset[proci]);
819 
820  tgtFaceOwners[tgtFacei] = newOwn;
821  tgtFaceNeighbours[tgtFacei] = -1;
822  }
823  }
824  // normal boundary face
825  else
826  {
827  label own = faceOs[i];
828  label nbr = faceNs[i];
829  if ((own != -1) && (nbr == -1))
830  {
831  label newOwn = cellOffset[proci] + faceOs[i];
832  label tgtFacei = nIntFaces++;
833 
834  tgtFaces[tgtFacei] = fcs[i];
835  add(tgtFaces[tgtFacei], pointOffset[proci]);
836 
837  tgtFaceOwners[tgtFacei] = newOwn;
838  tgtFaceNeighbours[tgtFacei] = -1;
839  }
840  }
841  }
842  }
843 
844 
845  if (debug)
846  {
847  // only merging points in debug mode
848 
849  labelList oldToNew;
850  label nChanged = Foam::inplaceMergePoints
851  (
852  tgtPoints,
853  SMALL,
854  false,
855  oldToNew
856  );
857 
858  if (nChanged)
859  {
860  Pout<< "Merged from " << oldToNew.size()
861  << " down to " << tgtPoints.size() << " points" << endl;
862 
863  for (auto& f : tgtFaces)
864  {
865  inplaceRenumber(oldToNew, f);
866  }
867  }
868  }
869 }
870 
871 
872 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:748
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:639
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
Definition: treeBoundBox.H:83
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all processes in communicator.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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) is 1 for serial run.
Definition: UPstream.H:656
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:289
const pointField & points
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
label nPoints
#define DebugInFunction
Report an information message using Foam::Info.
label inplaceMergePoints(PointList &points, const scalar mergeTol, const bool verbose, labelList &pointToUnique)
Inplace merge points, preserving the original point order. All points closer/equal mergeTol are to be...
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
int debug
Static debugging option.
Geometric merging of points. See below.
labelList f(nPoints)
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
"nonBlocking" : (MPI_Isend, MPI_Irecv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
List< label > labelList
A List of labels.
Definition: List.H:62
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
List< bool > boolList
A List of bools.
Definition: List.H:60
labelList cellIDs
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const pointField & pts
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157