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-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 
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  return autoPtr<mapDistribute>::New(std::move(sendMap));
236  break;
237  }
238  }
239 }
240 
241 
242 void Foam::meshToMesh::distributeCells
243 (
244  const mapDistribute& map,
245  const polyMesh& tgtMesh,
246  const globalIndex& globalI,
247  List<pointField>& points,
248  List<label>& nInternalFaces,
249  List<faceList>& faces,
250  List<labelList>& faceOwner,
251  List<labelList>& faceNeighbour,
252  List<labelList>& cellIDs,
253  List<labelList>& nbrProcIDs,
254  List<labelList>& procLocalFaceIDs
255 ) const
256 {
258  nInternalFaces.setSize(Pstream::nProcs(), 0);
259  faces.setSize(Pstream::nProcs());
260  faceOwner.setSize(Pstream::nProcs());
261  faceNeighbour.setSize(Pstream::nProcs());
263 
264  nbrProcIDs.setSize(Pstream::nProcs());;
265  procLocalFaceIDs.setSize(Pstream::nProcs());;
266 
267 
268  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
269 
270  for (const int domain : Pstream::allProcs())
271  {
272  const labelList& sendElems = map.subMap()[domain];
273 
274  if (sendElems.size())
275  {
276  // reverse cell map
277  labelList reverseCellMap(tgtMesh.nCells(), -1);
278  forAll(sendElems, subCelli)
279  {
280  reverseCellMap[sendElems[subCelli]] = subCelli;
281  }
282 
283  DynamicList<face> subFaces(tgtMesh.nFaces());
284  DynamicList<label> subFaceOwner(tgtMesh.nFaces());
285  DynamicList<label> subFaceNeighbour(tgtMesh.nFaces());
286 
287  DynamicList<label> subNbrProcIDs(tgtMesh.nFaces());
288  DynamicList<label> subProcLocalFaceIDs(tgtMesh.nFaces());
289 
290  label nInternal = 0;
291 
292  // internal faces
293  forAll(tgtMesh.faceNeighbour(), facei)
294  {
295  label own = tgtMesh.faceOwner()[facei];
296  label nbr = tgtMesh.faceNeighbour()[facei];
297  label subOwn = reverseCellMap[own];
298  label subNbr = reverseCellMap[nbr];
299 
300  if (subOwn != -1 && subNbr != -1)
301  {
302  nInternal++;
303 
304  if (subOwn < subNbr)
305  {
306  subFaces.append(tgtMesh.faces()[facei]);
307  subFaceOwner.append(subOwn);
308  subFaceNeighbour.append(subNbr);
309  subNbrProcIDs.append(-1);
310  subProcLocalFaceIDs.append(-1);
311  }
312  else
313  {
314  subFaces.append(tgtMesh.faces()[facei].reverseFace());
315  subFaceOwner.append(subNbr);
316  subFaceNeighbour.append(subOwn);
317  subNbrProcIDs.append(-1);
318  subProcLocalFaceIDs.append(-1);
319  }
320  }
321  }
322 
323  // boundary faces for new region
324  forAll(tgtMesh.faceNeighbour(), facei)
325  {
326  label own = tgtMesh.faceOwner()[facei];
327  label nbr = tgtMesh.faceNeighbour()[facei];
328  label subOwn = reverseCellMap[own];
329  label subNbr = reverseCellMap[nbr];
330 
331  if (subOwn != -1 && subNbr == -1)
332  {
333  subFaces.append(tgtMesh.faces()[facei]);
334  subFaceOwner.append(subOwn);
335  subFaceNeighbour.append(subNbr);
336  subNbrProcIDs.append(-1);
337  subProcLocalFaceIDs.append(-1);
338  }
339  else if (subOwn == -1 && subNbr != -1)
340  {
341  subFaces.append(tgtMesh.faces()[facei].reverseFace());
342  subFaceOwner.append(subNbr);
343  subFaceNeighbour.append(subOwn);
344  subNbrProcIDs.append(-1);
345  subProcLocalFaceIDs.append(-1);
346  }
347  }
348 
349  // boundary faces of existing region
350  forAll(tgtMesh.boundaryMesh(), patchi)
351  {
352  const polyPatch& pp = tgtMesh.boundaryMesh()[patchi];
353  const auto* procPatch = isA<processorPolyPatch>(pp);
354 
355  // Store info for faces on processor patches
356  const label nbrProci =
357  (procPatch ? procPatch->neighbProcNo() : -1);
358 
359  forAll(pp, i)
360  {
361  label facei = pp.start() + i;
362  label own = tgtMesh.faceOwner()[facei];
363 
364  if (reverseCellMap[own] != -1)
365  {
366  subFaces.append(tgtMesh.faces()[facei]);
367  subFaceOwner.append(reverseCellMap[own]);
368  subFaceNeighbour.append(-1);
369  subNbrProcIDs.append(nbrProci);
370  subProcLocalFaceIDs.append(i);
371  }
372  }
373  }
374 
375  // reverse point map
376  labelList reversePointMap(tgtMesh.nPoints(), -1);
377  DynamicList<point> subPoints(tgtMesh.nPoints());
378  forAll(subFaces, subFacei)
379  {
380  face& f = subFaces[subFacei];
381  forAll(f, fp)
382  {
383  label pointi = f[fp];
384  if (reversePointMap[pointi] == -1)
385  {
386  reversePointMap[pointi] = subPoints.size();
387  subPoints.append(tgtMesh.points()[pointi]);
388  }
389 
390  f[fp] = reversePointMap[pointi];
391  }
392  }
393 
394  // tgt cells into global numbering
395  labelList globalElems(globalI.toGlobal(sendElems));
396 
397  if (debug > 1)
398  {
399  forAll(sendElems, i)
400  {
401  Pout<< "tgtProc:" << Pstream::myProcNo()
402  << " sending tgt cell " << sendElems[i]
403  << "[" << globalElems[i] << "]"
404  << " to srcProc " << domain << endl;
405  }
406  }
407 
408  // pass data
409  if (domain == Pstream::myProcNo())
410  {
411  // allocate my own data
412  points[Pstream::myProcNo()] = subPoints;
413  nInternalFaces[Pstream::myProcNo()] = nInternal;
414  faces[Pstream::myProcNo()] = subFaces;
415  faceOwner[Pstream::myProcNo()] = subFaceOwner;
416  faceNeighbour[Pstream::myProcNo()] = subFaceNeighbour;
417  cellIDs[Pstream::myProcNo()] = globalElems;
418  nbrProcIDs[Pstream::myProcNo()] = subNbrProcIDs;
419  procLocalFaceIDs[Pstream::myProcNo()] = subProcLocalFaceIDs;
420  }
421  else
422  {
423  // send data to other processor domains
424  UOPstream toDomain(domain, pBufs);
425 
426  toDomain
427  << subPoints
428  << nInternal
429  << subFaces
430  << subFaceOwner
431  << subFaceNeighbour
432  << globalElems
433  << subNbrProcIDs
434  << subProcLocalFaceIDs;
435  }
436  }
437  }
438 
439  // Start receiving
440  pBufs.finishedSends();
441 
442  // Consume
443  for (const int domain : Pstream::allProcs())
444  {
445  const labelList& recvElems = map.constructMap()[domain];
446 
447  if (domain != Pstream::myProcNo() && recvElems.size())
448  {
449  UIPstream str(domain, pBufs);
450 
451  str >> points[domain]
452  >> nInternalFaces[domain]
453  >> faces[domain]
454  >> faceOwner[domain]
455  >> faceNeighbour[domain]
456  >> cellIDs[domain]
457  >> nbrProcIDs[domain]
458  >> procLocalFaceIDs[domain];
459  }
460 
461  if (debug)
462  {
463  Pout<< "Target mesh send sizes[" << domain << "]"
464  << ": points="<< points[domain].size()
465  << ", faces=" << faces[domain].size()
466  << ", nInternalFaces=" << nInternalFaces[domain]
467  << ", faceOwn=" << faceOwner[domain].size()
468  << ", faceNbr=" << faceNeighbour[domain].size()
469  << ", cellIDs=" << cellIDs[domain].size() << endl;
470  }
471  }
472 }
473 
474 
475 void Foam::meshToMesh::distributeAndMergeCells
476 (
477  const mapDistribute& map,
478  const polyMesh& tgt,
479  const globalIndex& globalI,
480  pointField& tgtPoints,
481  faceList& tgtFaces,
482  labelList& tgtFaceOwners,
483  labelList& tgtFaceNeighbours,
484  labelList& tgtCellIDs
485 ) const
486 {
487  // Exchange per-processor data
488  List<pointField> allPoints;
489  List<label> allNInternalFaces;
490  List<faceList> allFaces;
491  List<labelList> allFaceOwners;
492  List<labelList> allFaceNeighbours;
493  List<labelList> allTgtCellIDs;
494 
495  // Per target mesh face the neighbouring proc and index in
496  // processor patch (all -1 for normal boundary face)
497  List<labelList> allNbrProcIDs;
498  List<labelList> allProcLocalFaceIDs;
499 
500  distributeCells
501  (
502  map,
503  tgt,
504  globalI,
505  allPoints,
506  allNInternalFaces,
507  allFaces,
508  allFaceOwners,
509  allFaceNeighbours,
510  allTgtCellIDs,
511  allNbrProcIDs,
512  allProcLocalFaceIDs
513  );
514 
515  // Convert lists into format that can be used to generate a valid polyMesh
516  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
517  //
518  // Points and cells are collected into single flat lists:
519  // - i.e. proc0, proc1 ... procN
520  //
521  // Faces need to be sorted after collection to that internal faces are
522  // contiguous, followed by all boundary faces
523  //
524  // Processor patch faces between included cells on neighbouring processors
525  // are converted into internal faces
526  //
527  // Face list structure:
528  // - Per processor:
529  // - internal faces
530  // - processor faces that have been converted into internal faces
531  // - Followed by all boundary faces
532  // - from 'normal' boundary faces
533  // - from singularly-sided processor patch faces
534 
535 
536  // Number of internal+coupled faces
537  labelList allNIntCoupledFaces(allNInternalFaces);
538 
539  // Starting offset for points
540  label nPoints = 0;
541  labelList pointOffset(Pstream::nProcs(), Zero);
542  forAll(allPoints, proci)
543  {
544  pointOffset[proci] = nPoints;
545  nPoints += allPoints[proci].size();
546  }
547 
548  // Starting offset for cells
549  label nCells = 0;
550  labelList cellOffset(Pstream::nProcs(), Zero);
551  forAll(allTgtCellIDs, proci)
552  {
553  cellOffset[proci] = nCells;
554  nCells += allTgtCellIDs[proci].size();
555  }
556 
557  // Count any coupled faces
558  typedef FixedList<label, 3> label3;
559  typedef HashTable<label, label3> procCoupleInfo;
560  procCoupleInfo procFaceToGlobalCell;
561 
562  forAll(allNbrProcIDs, proci)
563  {
564  const labelList& nbrProci = allNbrProcIDs[proci];
565  const labelList& localFacei = allProcLocalFaceIDs[proci];
566 
567  forAll(nbrProci, i)
568  {
569  if (nbrProci[i] != -1 && localFacei[i] != -1)
570  {
571  label3 key;
572  key[0] = min(proci, nbrProci[i]);
573  key[1] = max(proci, nbrProci[i]);
574  key[2] = localFacei[i];
575 
576  const auto fnd = procFaceToGlobalCell.cfind(key);
577 
578  if (!fnd.good())
579  {
580  procFaceToGlobalCell.insert(key, -1);
581  }
582  else
583  {
584  if (debug > 1)
585  {
586  Pout<< "Additional internal face between procs:"
587  << key[0] << " and " << key[1]
588  << " across local face " << key[2] << endl;
589  }
590 
591  allNIntCoupledFaces[proci]++;
592  }
593  }
594  }
595  }
596 
597 
598  // Starting offset for internal faces
599  label nIntFaces = 0;
600  label nFacesTotal = 0;
601  labelList internalFaceOffset(Pstream::nProcs(), Zero);
602  forAll(allNIntCoupledFaces, proci)
603  {
604  label nCoupledFaces =
605  allNIntCoupledFaces[proci] - allNInternalFaces[proci];
606 
607  internalFaceOffset[proci] = nIntFaces;
608  nIntFaces += allNIntCoupledFaces[proci];
609  nFacesTotal += allFaceOwners[proci].size() - nCoupledFaces;
610  }
611 
612  tgtPoints.setSize(nPoints);
613  tgtFaces.setSize(nFacesTotal);
614  tgtFaceOwners.setSize(nFacesTotal);
615  tgtFaceNeighbours.setSize(nFacesTotal);
616  tgtCellIDs.setSize(nCells);
617 
618  // Insert points
619  forAll(allPoints, proci)
620  {
621  const pointField& pts = allPoints[proci];
622  SubList<point>(tgtPoints, pts.size(), pointOffset[proci]) = pts;
623  }
624 
625  // Insert cellIDs
626  forAll(allTgtCellIDs, proci)
627  {
628  const labelList& cellIDs = allTgtCellIDs[proci];
629  SubList<label>(tgtCellIDs, cellIDs.size(), cellOffset[proci]) = cellIDs;
630  }
631 
632 
633  // Insert internal faces (from internal faces)
634  forAll(allFaces, proci)
635  {
636  const faceList& fcs = allFaces[proci];
637  const labelList& faceOs = allFaceOwners[proci];
638  const labelList& faceNs = allFaceNeighbours[proci];
639 
640  SubList<face> slice
641  (
642  tgtFaces,
643  allNInternalFaces[proci],
644  internalFaceOffset[proci]
645  );
646  slice = SubList<face>(fcs, allNInternalFaces[proci]);
647  forAll(slice, i)
648  {
649  add(slice[i], pointOffset[proci]);
650  }
651 
652  SubField<label> ownSlice
653  (
654  tgtFaceOwners,
655  allNInternalFaces[proci],
656  internalFaceOffset[proci]
657  );
658  ownSlice = SubField<label>(faceOs, allNInternalFaces[proci]);
659  add(ownSlice, cellOffset[proci]);
660 
661  SubField<label> nbrSlice
662  (
663  tgtFaceNeighbours,
664  allNInternalFaces[proci],
665  internalFaceOffset[proci]
666  );
667  nbrSlice = SubField<label>(faceNs, allNInternalFaces[proci]);
668  add(nbrSlice, cellOffset[proci]);
669 
670  internalFaceOffset[proci] += allNInternalFaces[proci];
671  }
672 
673 
674  // Insert internal faces (from coupled face-pairs)
675  forAll(allNbrProcIDs, proci)
676  {
677  const labelList& nbrProci = allNbrProcIDs[proci];
678  const labelList& localFacei = allProcLocalFaceIDs[proci];
679  const labelList& faceOs = allFaceOwners[proci];
680  const faceList& fcs = allFaces[proci];
681 
682  forAll(nbrProci, i)
683  {
684  if (nbrProci[i] != -1 && localFacei[i] != -1)
685  {
686  label3 key;
687  key[0] = min(proci, nbrProci[i]);
688  key[1] = max(proci, nbrProci[i]);
689  key[2] = localFacei[i];
690 
691  auto fnd = procFaceToGlobalCell.find(key);
692 
693  if (fnd.good())
694  {
695  label tgtFacei = fnd();
696  if (tgtFacei == -1)
697  {
698  // on first visit store the new cell on this side
699  fnd() = cellOffset[proci] + faceOs[i];
700  }
701  else
702  {
703  // get owner and neighbour in new cell numbering
704  label newOwn = cellOffset[proci] + faceOs[i];
705  label newNbr = fnd();
706  label tgtFacei = internalFaceOffset[proci]++;
707 
708  if (debug > 1)
709  {
710  Pout<< " proc " << proci
711  << "\tinserting face:" << tgtFacei
712  << " connection between owner " << newOwn
713  << " and neighbour " << newNbr
714  << endl;
715  }
716 
717  if (newOwn < newNbr)
718  {
719  // we have correct orientation
720  tgtFaces[tgtFacei] = fcs[i];
721  tgtFaceOwners[tgtFacei] = newOwn;
722  tgtFaceNeighbours[tgtFacei] = newNbr;
723  }
724  else
725  {
726  // reverse orientation
727  tgtFaces[tgtFacei] = fcs[i].reverseFace();
728  tgtFaceOwners[tgtFacei] = newNbr;
729  tgtFaceNeighbours[tgtFacei] = newOwn;
730  }
731 
732  add(tgtFaces[tgtFacei], pointOffset[proci]);
733 
734  // mark with unique value
735  fnd() = -2;
736  }
737  }
738  }
739  }
740  }
741 
742 
743  forAll(allNbrProcIDs, proci)
744  {
745  const labelList& nbrProci = allNbrProcIDs[proci];
746  const labelList& localFacei = allProcLocalFaceIDs[proci];
747  const labelList& faceOs = allFaceOwners[proci];
748  const labelList& faceNs = allFaceNeighbours[proci];
749  const faceList& fcs = allFaces[proci];
750 
751  forAll(nbrProci, i)
752  {
753  // coupled boundary face
754  if (nbrProci[i] != -1 && localFacei[i] != -1)
755  {
756  label3 key;
757  key[0] = min(proci, nbrProci[i]);
758  key[1] = max(proci, nbrProci[i]);
759  key[2] = localFacei[i];
760 
761  label tgtFacei = procFaceToGlobalCell[key];
762 
763  if (tgtFacei == -1)
764  {
766  << "Unvisited " << key
767  << abort(FatalError);
768  }
769  else if (tgtFacei != -2)
770  {
771  label newOwn = cellOffset[proci] + faceOs[i];
772  label tgtFacei = nIntFaces++;
773 
774  if (debug > 1)
775  {
776  Pout<< " proc " << proci
777  << "\tinserting boundary face:" << tgtFacei
778  << " from coupled face " << key
779  << endl;
780  }
781 
782  tgtFaces[tgtFacei] = fcs[i];
783  add(tgtFaces[tgtFacei], pointOffset[proci]);
784 
785  tgtFaceOwners[tgtFacei] = newOwn;
786  tgtFaceNeighbours[tgtFacei] = -1;
787  }
788  }
789  // normal boundary face
790  else
791  {
792  label own = faceOs[i];
793  label nbr = faceNs[i];
794  if ((own != -1) && (nbr == -1))
795  {
796  label newOwn = cellOffset[proci] + faceOs[i];
797  label tgtFacei = nIntFaces++;
798 
799  tgtFaces[tgtFacei] = fcs[i];
800  add(tgtFaces[tgtFacei], pointOffset[proci]);
801 
802  tgtFaceOwners[tgtFacei] = newOwn;
803  tgtFaceNeighbours[tgtFacei] = -1;
804  }
805  }
806  }
807  }
808 
809 
810  if (debug)
811  {
812  // only merging points in debug mode
813 
814  labelList oldToNew;
815  label nChanged = Foam::inplaceMergePoints
816  (
817  tgtPoints,
818  SMALL,
819  false,
820  oldToNew
821  );
822 
823  if (nChanged)
824  {
825  Pout<< "Merged from " << oldToNew.size()
826  << " down to " << tgtPoints.size() << " points" << endl;
827 
828  for (auto& f : tgtFaces)
829  {
830  inplaceRenumber(oldToNew, f);
831  }
832  }
833  }
834 }
835 
836 
837 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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:1131
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
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:1004
List< treeBoundBox > treeBoundBoxList
A List of treeBoundBox.
Definition: treeBoundBox.H:83
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:1029
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
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:414
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. ...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
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:1020
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.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const pointField & pts
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133