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-------------------------------------------------------------------------------
11License
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
41Foam::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
81Foam::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
111Foam::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 {
164 }
165
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 const cellList& cells = tgt.cells();
183 const faceList& faces = tgt.faces();
184 const pointField& points = tgt.points();
185
186 labelListList sendMap;
187
188 {
189 // per processor indices into all segments to send
190 List<DynamicList<label>> dynSendMap(Pstream::nProcs());
191 label iniSize = floor(tgt.nCells()/Pstream::nProcs());
192
193 forAll(dynSendMap, proci)
194 {
195 dynSendMap[proci].setCapacity(iniSize);
196 }
197
198 // work array - whether src processor bb overlaps the tgt cell
199 // bounds
200 boolList procBbOverlaps(Pstream::nProcs());
201 forAll(cells, celli)
202 {
203 const cell& c = cells[celli];
204
205 // determine bounding box of tgt cell
206 boundBox cellBb(boundBox::invertedBox);
207 forAll(c, facei)
208 {
209 const face& f = faces[c[facei]];
210 forAll(f, fp)
211 {
212 cellBb.add(points, f);
213 }
214 }
215
216 // find the overlapping tgt cells on each src processor
217 (void)calcOverlappingProcs(procBb, cellBb, procBbOverlaps);
218
219 forAll(procBbOverlaps, proci)
220 {
221 if (procBbOverlaps[proci])
222 {
223 dynSendMap[proci].append(celli);
224 }
225 }
226 }
227
228 // convert dynamicList to labelList
229 sendMap.setSize(Pstream::nProcs());
230 forAll(sendMap, proci)
231 {
232 sendMap[proci].transfer(dynSendMap[proci]);
233 }
234 }
235
236 // debug printing
237 if (debug)
238 {
239 Pout<< "Of my " << cells.size()
240 << " target cells I need to send to:" << nl
241 << "\tproc\tcells" << endl;
242 forAll(sendMap, proci)
243 {
244 Pout<< '\t' << proci << '\t' << sendMap[proci].size()
245 << endl;
246 }
247 }
248
249
250 // send over how many tgt cells I need to receive from each
251 // processor
252 labelListList sendSizes(Pstream::nProcs());
253 sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
254 forAll(sendMap, proci)
255 {
256 sendSizes[Pstream::myProcNo()][proci] = sendMap[proci].size();
257 }
258 Pstream::allGatherList(sendSizes);
259
260
261 // determine order of receiving
262 labelListList constructMap(Pstream::nProcs());
263
264 label segmentI = 0;
265 forAll(constructMap, proci)
266 {
267 // what I need to receive is what other processor is sending
268 // to me
269 label nRecv = sendSizes[proci][Pstream::myProcNo()];
270 constructMap[proci].setSize(nRecv);
271
272 for (label i = 0; i < nRecv; i++)
273 {
274 constructMap[proci][i] = segmentI++;
275 }
276 }
277
279 (
280 segmentI, // size after construction
281 std::move(sendMap),
282 std::move(constructMap)
283 );
284
285 break;
286 }
287 }
288}
289
290
291void Foam::meshToMesh::distributeCells
292(
293 const mapDistribute& map,
294 const polyMesh& tgtMesh,
295 const globalIndex& globalI,
296 List<pointField>& points,
297 List<label>& nInternalFaces,
298 List<faceList>& faces,
299 List<labelList>& faceOwner,
300 List<labelList>& faceNeighbour,
301 List<labelList>& cellIDs,
302 List<labelList>& nbrProcIDs,
303 List<labelList>& procLocalFaceIDs
304) const
305{
307 nInternalFaces.setSize(Pstream::nProcs(), 0);
308 faces.setSize(Pstream::nProcs());
309 faceOwner.setSize(Pstream::nProcs());
310 faceNeighbour.setSize(Pstream::nProcs());
312
313 nbrProcIDs.setSize(Pstream::nProcs());;
314 procLocalFaceIDs.setSize(Pstream::nProcs());;
315
316
317 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
318
319 for (const int domain : Pstream::allProcs())
320 {
321 const labelList& sendElems = map.subMap()[domain];
322
323 if (sendElems.size())
324 {
325 // reverse cell map
326 labelList reverseCellMap(tgtMesh.nCells(), -1);
327 forAll(sendElems, subCelli)
328 {
329 reverseCellMap[sendElems[subCelli]] = subCelli;
330 }
331
332 DynamicList<face> subFaces(tgtMesh.nFaces());
333 DynamicList<label> subFaceOwner(tgtMesh.nFaces());
334 DynamicList<label> subFaceNeighbour(tgtMesh.nFaces());
335
336 DynamicList<label> subNbrProcIDs(tgtMesh.nFaces());
337 DynamicList<label> subProcLocalFaceIDs(tgtMesh.nFaces());
338
339 label nInternal = 0;
340
341 // internal faces
342 forAll(tgtMesh.faceNeighbour(), facei)
343 {
344 label own = tgtMesh.faceOwner()[facei];
345 label nbr = tgtMesh.faceNeighbour()[facei];
346 label subOwn = reverseCellMap[own];
347 label subNbr = reverseCellMap[nbr];
348
349 if (subOwn != -1 && subNbr != -1)
350 {
351 nInternal++;
352
353 if (subOwn < subNbr)
354 {
355 subFaces.append(tgtMesh.faces()[facei]);
356 subFaceOwner.append(subOwn);
357 subFaceNeighbour.append(subNbr);
358 subNbrProcIDs.append(-1);
359 subProcLocalFaceIDs.append(-1);
360 }
361 else
362 {
363 subFaces.append(tgtMesh.faces()[facei].reverseFace());
364 subFaceOwner.append(subNbr);
365 subFaceNeighbour.append(subOwn);
366 subNbrProcIDs.append(-1);
367 subProcLocalFaceIDs.append(-1);
368 }
369 }
370 }
371
372 // boundary faces for new region
373 forAll(tgtMesh.faceNeighbour(), facei)
374 {
375 label own = tgtMesh.faceOwner()[facei];
376 label nbr = tgtMesh.faceNeighbour()[facei];
377 label subOwn = reverseCellMap[own];
378 label subNbr = reverseCellMap[nbr];
379
380 if (subOwn != -1 && subNbr == -1)
381 {
382 subFaces.append(tgtMesh.faces()[facei]);
383 subFaceOwner.append(subOwn);
384 subFaceNeighbour.append(subNbr);
385 subNbrProcIDs.append(-1);
386 subProcLocalFaceIDs.append(-1);
387 }
388 else if (subOwn == -1 && subNbr != -1)
389 {
390 subFaces.append(tgtMesh.faces()[facei].reverseFace());
391 subFaceOwner.append(subNbr);
392 subFaceNeighbour.append(subOwn);
393 subNbrProcIDs.append(-1);
394 subProcLocalFaceIDs.append(-1);
395 }
396 }
397
398 // boundary faces of existing region
399 forAll(tgtMesh.boundaryMesh(), patchi)
400 {
401 const polyPatch& pp = tgtMesh.boundaryMesh()[patchi];
402 const auto* procPatch = isA<processorPolyPatch>(pp);
403
404 // Store info for faces on processor patches
405 const label nbrProci =
406 (procPatch ? procPatch->neighbProcNo() : -1);
407
408 forAll(pp, i)
409 {
410 label facei = pp.start() + i;
411 label own = tgtMesh.faceOwner()[facei];
412
413 if (reverseCellMap[own] != -1)
414 {
415 subFaces.append(tgtMesh.faces()[facei]);
416 subFaceOwner.append(reverseCellMap[own]);
417 subFaceNeighbour.append(-1);
418 subNbrProcIDs.append(nbrProci);
419 subProcLocalFaceIDs.append(i);
420 }
421 }
422 }
423
424 // reverse point map
425 labelList reversePointMap(tgtMesh.nPoints(), -1);
426 DynamicList<point> subPoints(tgtMesh.nPoints());
427 forAll(subFaces, subFacei)
428 {
429 face& f = subFaces[subFacei];
430 forAll(f, fp)
431 {
432 label pointi = f[fp];
433 if (reversePointMap[pointi] == -1)
434 {
435 reversePointMap[pointi] = subPoints.size();
436 subPoints.append(tgtMesh.points()[pointi]);
437 }
438
439 f[fp] = reversePointMap[pointi];
440 }
441 }
442
443 // tgt cells into global numbering
444 labelList globalElems(globalI.toGlobal(sendElems));
445
446 if (debug > 1)
447 {
448 forAll(sendElems, i)
449 {
450 Pout<< "tgtProc:" << Pstream::myProcNo()
451 << " sending tgt cell " << sendElems[i]
452 << "[" << globalElems[i] << "]"
453 << " to srcProc " << domain << endl;
454 }
455 }
456
457 // pass data
458 if (domain == Pstream::myProcNo())
459 {
460 // allocate my own data
461 points[Pstream::myProcNo()] = subPoints;
462 nInternalFaces[Pstream::myProcNo()] = nInternal;
463 faces[Pstream::myProcNo()] = subFaces;
464 faceOwner[Pstream::myProcNo()] = subFaceOwner;
465 faceNeighbour[Pstream::myProcNo()] = subFaceNeighbour;
466 cellIDs[Pstream::myProcNo()] = globalElems;
467 nbrProcIDs[Pstream::myProcNo()] = subNbrProcIDs;
468 procLocalFaceIDs[Pstream::myProcNo()] = subProcLocalFaceIDs;
469 }
470 else
471 {
472 // send data to other processor domains
473 UOPstream toDomain(domain, pBufs);
474
475 toDomain
476 << subPoints
477 << nInternal
478 << subFaces
479 << subFaceOwner
480 << subFaceNeighbour
481 << globalElems
482 << subNbrProcIDs
483 << subProcLocalFaceIDs;
484 }
485 }
486 }
487
488 // Start receiving
489 pBufs.finishedSends();
490
491 // Consume
492 for (const int domain : Pstream::allProcs())
493 {
494 const labelList& recvElems = map.constructMap()[domain];
495
496 if (domain != Pstream::myProcNo() && recvElems.size())
497 {
498 UIPstream str(domain, pBufs);
499
500 str >> points[domain]
501 >> nInternalFaces[domain]
502 >> faces[domain]
503 >> faceOwner[domain]
504 >> faceNeighbour[domain]
505 >> cellIDs[domain]
506 >> nbrProcIDs[domain]
507 >> procLocalFaceIDs[domain];
508 }
509
510 if (debug)
511 {
512 Pout<< "Target mesh send sizes[" << domain << "]"
513 << ": points="<< points[domain].size()
514 << ", faces=" << faces[domain].size()
515 << ", nInternalFaces=" << nInternalFaces[domain]
516 << ", faceOwn=" << faceOwner[domain].size()
517 << ", faceNbr=" << faceNeighbour[domain].size()
518 << ", cellIDs=" << cellIDs[domain].size() << endl;
519 }
520 }
521}
522
523
524void Foam::meshToMesh::distributeAndMergeCells
525(
526 const mapDistribute& map,
527 const polyMesh& tgt,
528 const globalIndex& globalI,
529 pointField& tgtPoints,
530 faceList& tgtFaces,
531 labelList& tgtFaceOwners,
532 labelList& tgtFaceNeighbours,
533 labelList& tgtCellIDs
534) const
535{
536 // Exchange per-processor data
537 List<pointField> allPoints;
538 List<label> allNInternalFaces;
539 List<faceList> allFaces;
540 List<labelList> allFaceOwners;
541 List<labelList> allFaceNeighbours;
542 List<labelList> allTgtCellIDs;
543
544 // Per target mesh face the neighbouring proc and index in
545 // processor patch (all -1 for normal boundary face)
546 List<labelList> allNbrProcIDs;
547 List<labelList> allProcLocalFaceIDs;
548
549 distributeCells
550 (
551 map,
552 tgt,
553 globalI,
554 allPoints,
555 allNInternalFaces,
556 allFaces,
557 allFaceOwners,
558 allFaceNeighbours,
559 allTgtCellIDs,
560 allNbrProcIDs,
561 allProcLocalFaceIDs
562 );
563
564 // Convert lists into format that can be used to generate a valid polyMesh
565 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
566 //
567 // Points and cells are collected into single flat lists:
568 // - i.e. proc0, proc1 ... procN
569 //
570 // Faces need to be sorted after collection to that internal faces are
571 // contiguous, followed by all boundary faces
572 //
573 // Processor patch faces between included cells on neighbouring processors
574 // are converted into internal faces
575 //
576 // Face list structure:
577 // - Per processor:
578 // - internal faces
579 // - processor faces that have been converted into internal faces
580 // - Followed by all boundary faces
581 // - from 'normal' boundary faces
582 // - from singularly-sided processor patch faces
583
584
585 // Number of internal+coupled faces
586 labelList allNIntCoupledFaces(allNInternalFaces);
587
588 // Starting offset for points
589 label nPoints = 0;
590 labelList pointOffset(Pstream::nProcs(), Zero);
591 forAll(allPoints, proci)
592 {
593 pointOffset[proci] = nPoints;
594 nPoints += allPoints[proci].size();
595 }
596
597 // Starting offset for cells
598 label nCells = 0;
599 labelList cellOffset(Pstream::nProcs(), Zero);
600 forAll(allTgtCellIDs, proci)
601 {
602 cellOffset[proci] = nCells;
603 nCells += allTgtCellIDs[proci].size();
604 }
605
606 // Count any coupled faces
607 typedef FixedList<label, 3> label3;
608 typedef HashTable<label, label3> procCoupleInfo;
609 procCoupleInfo procFaceToGlobalCell;
610
611 forAll(allNbrProcIDs, proci)
612 {
613 const labelList& nbrProci = allNbrProcIDs[proci];
614 const labelList& localFacei = allProcLocalFaceIDs[proci];
615
616 forAll(nbrProci, i)
617 {
618 if (nbrProci[i] != -1 && localFacei[i] != -1)
619 {
620 label3 key;
621 key[0] = min(proci, nbrProci[i]);
622 key[1] = max(proci, nbrProci[i]);
623 key[2] = localFacei[i];
624
625 const auto fnd = procFaceToGlobalCell.cfind(key);
626
627 if (!fnd.found())
628 {
629 procFaceToGlobalCell.insert(key, -1);
630 }
631 else
632 {
633 if (debug > 1)
634 {
635 Pout<< "Additional internal face between procs:"
636 << key[0] << " and " << key[1]
637 << " across local face " << key[2] << endl;
638 }
639
640 allNIntCoupledFaces[proci]++;
641 }
642 }
643 }
644 }
645
646
647 // Starting offset for internal faces
648 label nIntFaces = 0;
649 label nFacesTotal = 0;
650 labelList internalFaceOffset(Pstream::nProcs(), Zero);
651 forAll(allNIntCoupledFaces, proci)
652 {
653 label nCoupledFaces =
654 allNIntCoupledFaces[proci] - allNInternalFaces[proci];
655
656 internalFaceOffset[proci] = nIntFaces;
657 nIntFaces += allNIntCoupledFaces[proci];
658 nFacesTotal += allFaceOwners[proci].size() - nCoupledFaces;
659 }
660
661 tgtPoints.setSize(nPoints);
662 tgtFaces.setSize(nFacesTotal);
663 tgtFaceOwners.setSize(nFacesTotal);
664 tgtFaceNeighbours.setSize(nFacesTotal);
665 tgtCellIDs.setSize(nCells);
666
667 // Insert points
668 forAll(allPoints, proci)
669 {
670 const pointField& pts = allPoints[proci];
671 SubList<point>(tgtPoints, pts.size(), pointOffset[proci]) = pts;
672 }
673
674 // Insert cellIDs
675 forAll(allTgtCellIDs, proci)
676 {
677 const labelList& cellIDs = allTgtCellIDs[proci];
678 SubList<label>(tgtCellIDs, cellIDs.size(), cellOffset[proci]) = cellIDs;
679 }
680
681
682 // Insert internal faces (from internal faces)
683 forAll(allFaces, proci)
684 {
685 const faceList& fcs = allFaces[proci];
686 const labelList& faceOs = allFaceOwners[proci];
687 const labelList& faceNs = allFaceNeighbours[proci];
688
689 SubList<face> slice
690 (
691 tgtFaces,
692 allNInternalFaces[proci],
693 internalFaceOffset[proci]
694 );
695 slice = SubList<face>(fcs, allNInternalFaces[proci]);
696 forAll(slice, i)
697 {
698 add(slice[i], pointOffset[proci]);
699 }
700
701 SubField<label> ownSlice
702 (
703 tgtFaceOwners,
704 allNInternalFaces[proci],
705 internalFaceOffset[proci]
706 );
707 ownSlice = SubField<label>(faceOs, allNInternalFaces[proci]);
708 add(ownSlice, cellOffset[proci]);
709
710 SubField<label> nbrSlice
711 (
712 tgtFaceNeighbours,
713 allNInternalFaces[proci],
714 internalFaceOffset[proci]
715 );
716 nbrSlice = SubField<label>(faceNs, allNInternalFaces[proci]);
717 add(nbrSlice, cellOffset[proci]);
718
719 internalFaceOffset[proci] += allNInternalFaces[proci];
720 }
721
722
723 // Insert internal faces (from coupled face-pairs)
724 forAll(allNbrProcIDs, proci)
725 {
726 const labelList& nbrProci = allNbrProcIDs[proci];
727 const labelList& localFacei = allProcLocalFaceIDs[proci];
728 const labelList& faceOs = allFaceOwners[proci];
729 const faceList& fcs = allFaces[proci];
730
731 forAll(nbrProci, i)
732 {
733 if (nbrProci[i] != -1 && localFacei[i] != -1)
734 {
735 label3 key;
736 key[0] = min(proci, nbrProci[i]);
737 key[1] = max(proci, nbrProci[i]);
738 key[2] = localFacei[i];
739
740 auto fnd = procFaceToGlobalCell.find(key);
741
742 if (fnd.found())
743 {
744 label tgtFacei = fnd();
745 if (tgtFacei == -1)
746 {
747 // on first visit store the new cell on this side
748 fnd() = cellOffset[proci] + faceOs[i];
749 }
750 else
751 {
752 // get owner and neighbour in new cell numbering
753 label newOwn = cellOffset[proci] + faceOs[i];
754 label newNbr = fnd();
755 label tgtFacei = internalFaceOffset[proci]++;
756
757 if (debug > 1)
758 {
759 Pout<< " proc " << proci
760 << "\tinserting face:" << tgtFacei
761 << " connection between owner " << newOwn
762 << " and neighbour " << newNbr
763 << endl;
764 }
765
766 if (newOwn < newNbr)
767 {
768 // we have correct orientation
769 tgtFaces[tgtFacei] = fcs[i];
770 tgtFaceOwners[tgtFacei] = newOwn;
771 tgtFaceNeighbours[tgtFacei] = newNbr;
772 }
773 else
774 {
775 // reverse orientation
776 tgtFaces[tgtFacei] = fcs[i].reverseFace();
777 tgtFaceOwners[tgtFacei] = newNbr;
778 tgtFaceNeighbours[tgtFacei] = newOwn;
779 }
780
781 add(tgtFaces[tgtFacei], pointOffset[proci]);
782
783 // mark with unique value
784 fnd() = -2;
785 }
786 }
787 }
788 }
789 }
790
791
792 forAll(allNbrProcIDs, proci)
793 {
794 const labelList& nbrProci = allNbrProcIDs[proci];
795 const labelList& localFacei = allProcLocalFaceIDs[proci];
796 const labelList& faceOs = allFaceOwners[proci];
797 const labelList& faceNs = allFaceNeighbours[proci];
798 const faceList& fcs = allFaces[proci];
799
800 forAll(nbrProci, i)
801 {
802 // coupled boundary face
803 if (nbrProci[i] != -1 && localFacei[i] != -1)
804 {
805 label3 key;
806 key[0] = min(proci, nbrProci[i]);
807 key[1] = max(proci, nbrProci[i]);
808 key[2] = localFacei[i];
809
810 label tgtFacei = procFaceToGlobalCell[key];
811
812 if (tgtFacei == -1)
813 {
815 << "Unvisited " << key
816 << abort(FatalError);
817 }
818 else if (tgtFacei != -2)
819 {
820 label newOwn = cellOffset[proci] + faceOs[i];
821 label tgtFacei = nIntFaces++;
822
823 if (debug > 1)
824 {
825 Pout<< " proc " << proci
826 << "\tinserting boundary face:" << tgtFacei
827 << " from coupled face " << key
828 << endl;
829 }
830
831 tgtFaces[tgtFacei] = fcs[i];
832 add(tgtFaces[tgtFacei], pointOffset[proci]);
833
834 tgtFaceOwners[tgtFacei] = newOwn;
835 tgtFaceNeighbours[tgtFacei] = -1;
836 }
837 }
838 // normal boundary face
839 else
840 {
841 label own = faceOs[i];
842 label nbr = faceNs[i];
843 if ((own != -1) && (nbr == -1))
844 {
845 label newOwn = cellOffset[proci] + faceOs[i];
846 label tgtFacei = nIntFaces++;
847
848 tgtFaces[tgtFacei] = fcs[i];
849 add(tgtFaces[tgtFacei], pointOffset[proci]);
850
851 tgtFaceOwners[tgtFacei] = newOwn;
852 tgtFaceNeighbours[tgtFacei] = -1;
853 }
854 }
855 }
856 }
857
858
859 if (debug)
860 {
861 // only merging points in debug mode
862
863 labelList oldToNew;
864 label nChanged = Foam::inplaceMergePoints
865 (
866 tgtPoints,
867 SMALL,
868 false,
869 oldToNew
870 );
871
872 if (nChanged)
873 {
874 Pout<< "Merged from " << oldToNew.size()
875 << " down to " << tgtPoints.size() << " points" << endl;
876
877 for (auto& f : tgtFaces)
878 {
879 inplaceRenumber(oldToNew, f);
880 }
881 }
882 }
883}
884
885
886// ************************************************************************* //
labelList cellIDs
void setSize(const label n)
Alias for resize()
Definition: List.H:218
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void allGatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void broadcast(Type &value, const label comm=UPstream::worldComm)
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
int myProcNo() const noexcept
Return processor number.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
label nPoints
const cellShapeList & cells
Geometric merging of points. See below.
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
const dimensionedScalar c
Speed of light in a vacuum.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:108
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
List< label > labelList
A List of labels.
Definition: List.H:66
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
List< treeBoundBox > treeBoundBoxList
List of bounding boxes.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
label inplaceMergePoints(PointList &points, const scalar mergeTol, const bool verbose, labelList &pointToUnique)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333