Loading...
Searching...
No Matches
cyclicACMIGAMGInterface.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2025 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 "AMIInterpolation.H"
32#include "Map.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
44 );
46 (
50 );
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
56Foam::cyclicACMIGAMGInterface::cyclicACMIGAMGInterface
57(
58 const label index,
59 const lduInterfacePtrsList& coarseInterfaces,
60 const lduInterface& fineInterface,
61 const labelField& localRestrictAddressing,
62 const labelField& neighbourRestrictAddressing,
63 const label fineLevelIndex,
64 const label coarseComm
65)
66:
68 (
69 index,
70 coarseInterfaces
71 ),
72 neighbPatchID_
73 (
74 refCast<const cyclicACMILduInterface>(fineInterface).neighbPatchID()
75 ),
76 owner_
77 (
78 refCast<const cyclicACMILduInterface>(fineInterface).owner()
79 ),
80 forwardT_
81 (
82 refCast<const cyclicACMILduInterface>(fineInterface).forwardT()
83 ),
84 reverseT_
85 (
86 refCast<const cyclicACMILduInterface>(fineInterface).reverseT()
87 ),
88 myProcNo_(-1)
89{
90 const auto& fineCyclicACMIInterface =
92
93 // Construct face agglomeration from cell agglomeration
94 {
95 // From coarse face to cell
96 DynamicList<label> dynFaceCells(localRestrictAddressing.size());
97
98 // From face to coarse face
99 DynamicList<label> dynFaceRestrictAddressing
100 (
101 localRestrictAddressing.size()
102 );
103
104 Map<label> masterToCoarseFace(localRestrictAddressing.size());
105
106 for (const label curMaster : localRestrictAddressing)
107 {
108 const auto iter = masterToCoarseFace.cfind(curMaster);
109
110 if (iter.good())
111 {
112 // Already have coarse face
113 dynFaceRestrictAddressing.append(iter.val());
114 }
115 else
116 {
117 // New coarse face
118 const label coarseI = dynFaceCells.size();
119 dynFaceRestrictAddressing.append(coarseI);
120 dynFaceCells.append(curMaster);
121 masterToCoarseFace.insert(curMaster, coarseI);
122 }
123 }
124
125 faceCells_.transfer(dynFaceCells);
126 faceRestrictAddressing_.transfer(dynFaceRestrictAddressing);
127 }
128
129
130 // On the owner side construct the AMI
131
132 if (fineCyclicACMIInterface.owner())
133 {
134 // Construct the neighbour side agglomeration (as the neighbour would
135 // do it so it the exact loop above using neighbourRestrictAddressing
136 // instead of localRestrictAddressing)
137
138 labelList nbrFaceRestrictAddressing;
139 {
140 // From face to coarse face
141 DynamicList<label> dynNbrFaceRestrictAddressing
142 (
143 neighbourRestrictAddressing.size()
144 );
145
146 Map<label> masterToCoarseFace(neighbourRestrictAddressing.size());
147
148 for (const label curMaster : neighbourRestrictAddressing)
149 {
150 const auto iter = masterToCoarseFace.cfind(curMaster);
151
152 if (iter.good())
153 {
154 // Already have coarse face
155 dynNbrFaceRestrictAddressing.append(iter.val());
156 }
157 else
158 {
159 // New coarse face
160 const label coarseI = masterToCoarseFace.size();
161 dynNbrFaceRestrictAddressing.append(coarseI);
162 masterToCoarseFace.insert(curMaster, coarseI);
163 }
164 }
165
166 nbrFaceRestrictAddressing.transfer(dynNbrFaceRestrictAddressing);
167 }
168
169
170 amiPtr_.reset
171 (
173 (
174 fineCyclicACMIInterface.AMI(),
176 nbrFaceRestrictAddressing
177 )
178 );
179 }
180}
181
182
183Foam::cyclicACMIGAMGInterface::cyclicACMIGAMGInterface
184(
185 const label index,
186 const lduInterfacePtrsList& coarseInterfaces,
187 Istream& is
188)
189:
190 GAMGInterface(index, coarseInterfaces, is),
191 neighbPatchID_(readLabel(is)),
192 owner_(readBool(is)),
193 forwardT_(is),
194 reverseT_(is),
195 myProcNo_(-1)
196{
197 const bool hasAMI(readBool(is));
198
199 if (hasAMI)
200 {
201 amiPtr_.reset(new AMIPatchToPatchInterpolation(is));
202
203 // Store originating ranks locally - used when processor agglomerating
204 // onto a processor that wasn't in the communicator originally (since
205 // it had no faces)
206 const label comm = AMI().comm();
207
208 if (comm != -1)
210 is >> myProcNo_;
211 }
212 }
213}
214
215
216Foam::cyclicACMIGAMGInterface::cyclicACMIGAMGInterface
217(
218 const label index,
219 const lduInterfacePtrsList& coarseInterfaces,
220 const lduInterface& fineInterface,
221 const labelList& interfaceMap,
222 const labelUList& faceCells,
223 const labelUList& faceRestrictAddresssing,
224 const labelUList& faceOffsets,
225 const lduInterfacePtrsList& allInterfaces,
226 const label coarseComm,
227 const label myProcNo,
228 const labelList& procAgglomMap
229)
230:
232 (
233 index,
234 coarseInterfaces,
235 faceCells,
236 faceRestrictAddresssing
237 ),
238 neighbPatchID_
239 (
240 interfaceMap.find
241 (
242 refCast
243 <
245 >(fineInterface).neighbPatchID()
246 )
247 ),
248 owner_
249 (
250 refCast<const cyclicACMILduInterface>(fineInterface).owner()
251 ),
252 forwardT_
253 (
254 refCast<const cyclicACMILduInterface>(fineInterface).forwardT()
255 ),
256 reverseT_
257 (
258 refCast<const cyclicACMILduInterface>(fineInterface).reverseT()
259 ),
260 myProcNo_(-1)
261{
262 if (!owner_)
263 {
264 return;
265 }
266
267
268 // Get stats, sizes from the input interfaces. For the global settings
269 // the problem is that the
270 // local processor might not have any valid interfaces so here just
271 // collect and do a global reduction afterwards.
272
273 // Structure to pack all. First element is used to decide who has the
274 // valid AMI.
275 typedef
276 Tuple2
277 <
278 label,
279 Tuple2
280 <
281 Tuple2
282 <
283 FixedList<bool, 4>,
284 scalar
285 >,
286 label
287 >
288 > AMIType;
289
290 AMIType globalInfo;
291 FixedList<bool, 4>& bools = globalInfo.second().first().first();
292
293 // Define aliases to make our life easier
294 label& firstValidAMI = globalInfo.first();
295 bool& requireMatch = bools[0];
296 bool& reverseTarget = bools[1];
297 bool& srcHasFlip = bools[2];
298 bool& tgtHasFlip = bools[3];
299 scalar& lowWeightCorrection = globalInfo.second().first().second();
300 label& singlePatchProc = globalInfo.second().second();
301
302 // Initialise all global variables
303 firstValidAMI = labelMax;
304 requireMatch = false;
305 reverseTarget = false;
306 srcHasFlip = false;
307 tgtHasFlip = false;
308 lowWeightCorrection = -1;
309 singlePatchProc = -1;
310
311 // Initialise all local variables
312 bool hasSrcMagSf = false;
313 bool hasSrcCentroids = false;
314 bool hasTgtMagSf = false;
315 label nSrc = 0;
316 label nTgt = 0;
317
318 forAll(allInterfaces, inti)
319 {
320 if (allInterfaces.set(inti))
321 {
322 const auto& intf =
323 refCast<const cyclicACMIGAMGInterface>(allInterfaces[inti]);
324
325 if (!intf.amiPtr_)
326 {
327 continue;
328 }
329
330 if (firstValidAMI == labelMax)
331 {
332 firstValidAMI = inti;
333 }
334
335 const auto& AMI = intf.AMI();
336
337 // Note:
338 // - AMI.singlePatchProc or AMI.distributed : global properties.
339 // Should be already parallel synced.
340 // - AMI.comm() : -1 on processors that are not involved
341 // - AMI.srcMap(),tgtMap() : not valid if AMI.comm() is -1
342
343 if (!AMI.distributed())
344 {
345 singlePatchProc = AMI.singlePatchProc();
346 }
347 else if (AMI.comm() != -1)
348 {
349 singlePatchProc = -1;
350 srcHasFlip =
351 srcHasFlip || AMI.srcMap().constructHasFlip();
352 tgtHasFlip =
353 tgtHasFlip || AMI.tgtMap().constructHasFlip();
354 }
355 requireMatch = AMI.requireMatch();
356 reverseTarget = AMI.reverseTarget();
357 lowWeightCorrection = AMI.lowWeightCorrection();
358
359 nSrc += AMI.srcAddress().size();
360 nTgt += AMI.tgtAddress().size();
361
362 if (AMI.srcMagSf().size())
363 {
364 hasSrcMagSf = true;
365 if (AMI.srcMagSf().size() != AMI.srcAddress().size())
366 {
368 << "srcMagSf size:" << AMI.srcMagSf().size()
369 << "srcAddress size:" << AMI.srcAddress().size()
370 << exit(FatalError);
371 }
372 }
373 if (AMI.srcCentroids().size())
374 {
375 hasSrcCentroids = true;
376 if (AMI.srcCentroids().size() != AMI.srcAddress().size())
377 {
379 << "srcCentroids size:" << AMI.srcCentroids().size()
380 << "srcAddress size:" << AMI.srcAddress().size()
381 << exit(FatalError);
382 }
383 }
384 if (AMI.tgtMagSf().size())
385 {
386 hasTgtMagSf = true;
387 if (AMI.tgtMagSf().size() != AMI.tgtAddress().size())
388 {
390 << "tgtMagSf size:" << AMI.tgtMagSf().size()
391 << "tgtAddress size:" << AMI.tgtAddress().size()
392 << exit(FatalError);
393 }
394 }
395 }
396 }
397
398
399 // Reduce global information in case one of the coarse ranks does not
400 // have an input AMI to get data from. Could use minFirstEqOp from Tuple2
401 // instead ...
403 (
404 globalInfo,
405 [](AMIType& x, const AMIType& y)
406 {
407 if (y.first() < x.first())
408 {
409 x = y;
410 }
411 },
413 coarseComm
414 );
415
417 << "Input amis :"
418 << " singlePatchProc:" << singlePatchProc
419 << " srcHasFlip:" << srcHasFlip
420 << " tgtHasFlip:" << tgtHasFlip
421 << " requireMatch:" << requireMatch
422 << " reverseTarget:" << reverseTarget
423 << " lowWeightCorrection:" << lowWeightCorrection
424 << " hasSrcMagSf:" << hasSrcMagSf
425 << " hasSrcCentroids:" << hasSrcCentroids
426 << " hasTgtMagSf:" << hasTgtMagSf
427 << " nSrc:" << nSrc
428 << " nTgt:" << nTgt
429 << endl;
430
431
432 labelListList srcAddress;
433 scalarListList srcWeights;
434 scalarList srcMagSf;
435 // Needed?
436 pointListList srcCentroids;
437
438 labelListList tgtAddress;
439 scalarListList tgtWeights;
440 scalarList tgtMagSf;
441
442
443 // Map to send src side data to tgt side
444 autoPtr<mapDistribute> srcToTgtMap;
445
446 // Map to send tgt side data to src side
447 autoPtr<mapDistribute> tgtToSrcMap;
448
449 if (singlePatchProc == -1)
450 {
451 // Find ranks that agglomerate together
452 const label myAgglom = UPstream::myProcNo(coarseComm);
453
454 // Per input map either -1 or the index in the maps that is local
455 // data.
456 labelList localRanks(allInterfaces.size(), -1);
457 // From rank in coarse communicator back to rank in original (fine)
458 // communicator.
459 labelListList newToOldRanks;
460 {
461 // Pass 1: count number of valid maps
462 label nOldRanks = 0;
463 forAll(allInterfaces, inti)
464 {
465 if (allInterfaces.set(inti))
466 {
468 (
469 allInterfaces[inti]
470 );
471
472 if (!intf.amiPtr_ || intf.AMI().comm() == -1)
473 {
474 continue;
475 }
476 nOldRanks++;
477 }
478 }
479
480 // Pass 2: collect
481 DynamicList<label> oldRanks(nOldRanks);
482 forAll(allInterfaces, inti)
483 {
484 if (allInterfaces.set(inti))
485 {
487 (
488 allInterfaces[inti]
489 );
490
491 if (!intf.amiPtr_ || intf.AMI().comm() == -1)
492 {
493 continue;
494 }
495
496 label fineRank = -1;
497 if (intf.myProcNo() == -1)
498 {
499 // The interface was already local so got never
500 // sent across so myProcNo_ is never set ...
501 fineRank = UPstream::myProcNo(intf.AMI().comm());
502 }
503 else
504 {
505 fineRank = intf.myProcNo();
506 }
507
508 oldRanks.append(fineRank);
509 localRanks[inti] = fineRank;
510 }
511 }
512
513 // Pull individual parts together - this is the only communication
514 // needed.
515 newToOldRanks = Pstream::listGatherValues
516 (
517 labelList(std::move(oldRanks)),
518 coarseComm
519 );
520 Pstream::broadcast(newToOldRanks, coarseComm);
521 }
522
523
524 // Create combined maps
525 UPtrList<const mapDistribute> srcMaps(allInterfaces.size());
526 UPtrList<const mapDistribute> tgtMaps(allInterfaces.size());
527 forAll(allInterfaces, inti)
528 {
529 if (allInterfaces.set(inti))
530 {
532 (
533 allInterfaces[inti]
534 );
535
536 if (!intf.amiPtr_)
537 {
538 // Should not be in allInterfaces?
539 continue;
540 }
541
542 const auto& AMI = intf.AMI();
543
544 if (AMI.comm() != -1)
545 {
546 srcMaps.set(inti, &AMI.srcMap());
547 tgtMaps.set(inti, &AMI.tgtMap());
548 }
549 }
550 }
551
552
553 // Offsets for slots into results of srcToTgtMap
554 labelList srcStartOfLocal;
555 List<Map<label>> srcCompactMaps;
556
557 srcToTgtMap.reset
558 (
559 new mapDistribute
560 (
561 srcMaps,
562 localRanks, // per src map which rank represents local data
563 coarseComm,
564 newToOldRanks, // destination rank to source ranks
565 srcStartOfLocal,
566 srcCompactMaps
567 )
568 );
569
570
571 // Assemble tgtAddress
572 tgtAddress.setSize(nTgt);
573 if (tgtAddress.size())
574 {
575 label alli = 0;
576 forAll(allInterfaces, inti)
577 {
578 if (allInterfaces.set(inti))
579 {
581 (
582 allInterfaces[inti]
583 );
584
585 if (!intf.amiPtr_)
586 {
587 continue;
588 }
589
590 const auto& AMI = intf.AMI();
591 const auto& tgtSlots = AMI.tgtAddress();
592 const label localSize =
593 srcStartOfLocal[inti+1]
594 - srcStartOfLocal[inti];
595
596 forAll(tgtSlots, tgti)
597 {
598 // Append old slots: copy old values and adapt
599 auto& newSlots = tgtAddress[alli++];
600 newSlots = tgtSlots[tgti];
601
602 // Renumber to new indices
604 (
605 newSlots,
606 localSize,
607 srcStartOfLocal[inti],
608 srcCompactMaps[inti],
609 srcHasFlip //hasFlip
610 );
611
612 for (const label slot : newSlots)
613 {
614 if
615 (
616 slot < 0
617 || slot >= srcToTgtMap().constructSize()
618 )
619 {
620 FatalErrorInFunction << " newSlots:"
621 << newSlots << exit(FatalError);
622 }
623 }
624 }
625 }
626 }
627
628 if (nTgt != alli)
629 {
630 FatalErrorInFunction << "nTgt:" << nTgt
631 << " alli:" << alli << exit(FatalError);
632 }
633 }
634
635 // Offsets for slots into results of tgtToSrcMap
636 labelList tgtStartOfLocal;
637 List<Map<label>> tgtCompactMaps;
638
639 tgtToSrcMap.reset
640 (
641 new mapDistribute
642 (
643 tgtMaps,
644 localRanks,
645 coarseComm,
646 newToOldRanks,
647 tgtStartOfLocal,
648 tgtCompactMaps
649 )
650 );
651
652
653 // Assemble srcAddress
654 srcAddress.setSize(nSrc);
655 if (srcAddress.size())
656 {
657 label alli = 0;
658 forAll(allInterfaces, inti)
659 {
660 if (allInterfaces.set(inti))
661 {
663 (
664 allInterfaces[inti]
665 );
666
667 if (!intf.amiPtr_)
668 {
669 continue;
670 }
671
672 const auto& AMI = intf.AMI();
673 const auto& srcSlots = AMI.srcAddress();
674 const label localSize =
675 tgtStartOfLocal[inti+1]
676 - tgtStartOfLocal[inti];
677
678 forAll(srcSlots, srci)
679 {
680 // Append old slots: copy old values and adapt
681 auto& newSlots = srcAddress[alli++];
682 newSlots = srcSlots[srci];
683 // Renumber to new indices
685 (
686 newSlots,
687 localSize,
688 tgtStartOfLocal[inti],
689 tgtCompactMaps[inti],
690 tgtHasFlip
691 );
692
693 for (const label slot : newSlots)
694 {
695 if
696 (
697 slot < 0
698 || slot >= tgtToSrcMap().constructSize()
699 )
700 {
701 FatalErrorInFunction << " newSlots:"
702 << newSlots << exit(FatalError);
703 }
704 }
705 }
706 }
707 }
708
709 if (nSrc != alli)
710 {
711 FatalErrorInFunction << "nSrc:" << nSrc
712 << " alli:" << alli << exit(FatalError);
713 }
714 }
715
716
717 // Clean up: if no remote elements sent/received mark as
718 // non-distributed. We could do this at the start but this
719 // needs to take all the internal transport into account. Easier
720 // (but less efficient) to do afterwards now all is compacted.
721 {
722 const auto& map = srcToTgtMap().subMap();
723
724 bool usesRemote = false;
725 forAll(map, proci)
726 {
727 if (proci != myAgglom)
728 {
729 const auto& ss = srcToTgtMap().subMap()[proci];
730 const auto& sc = srcToTgtMap().constructMap()[proci];
731 const auto& ts = tgtToSrcMap().subMap()[proci];
732 const auto& tc = tgtToSrcMap().constructMap()[proci];
733
734 if (ss.size() || sc.size() || ts.size() || tc.size())
735 {
736 usesRemote = true;
737 break;
738 }
739 }
740 }
741
742 // We can't have a single rank become fully-local since we
743 // expect singlePatchProc to be synchronised. So make sure all
744 // have become local
745
746 if (!returnReduceOr(usesRemote, coarseComm))
747 {
748 DebugPout<< "** making fully local on new rank "
749 << myAgglom << " in comm:" << coarseComm << endl;
750 singlePatchProc = myAgglom;
751 srcToTgtMap.clear();
752 tgtToSrcMap.clear();
753 }
754 }
755 }
756 else
757 {
758 // src/tgt address are straight indices
759
760 srcAddress.setSize(nSrc);
761 tgtAddress.setSize(nTgt);
762
763 nSrc = 0;
764 nTgt = 0;
765 forAll(allInterfaces, inti)
766 {
767 if (allInterfaces.set(inti))
768 {
770 (
771 allInterfaces[inti]
772 );
773
774 if (!intf.amiPtr_)
775 {
776 continue;
777 }
778
779 const auto& AMI = intf.AMI();
780
781 const auto& srcA = AMI.srcAddress();
782 if (srcAddress.size())
783 {
784 label srci = nSrc;
785 forAll(srcA, i)
786 {
787 srcAddress[srci++] = srcA[i]+nTgt;
788 }
789 }
790
791 const auto& tgtA = AMI.tgtAddress();
792 if (tgtAddress.size())
793 {
794 label tgti = nTgt;
795 forAll(tgtA, i)
796 {
797 tgtAddress[tgti++] = tgtA[i]+nSrc;
798 }
799 }
800
801 nSrc += srcA.size();
802 nTgt += tgtA.size();
803 }
804 }
805 }
806
807 srcWeights.setSize(nSrc);
808 if (hasSrcMagSf)
809 {
810 srcMagSf.setSize(nSrc);
811 }
812 if (hasSrcCentroids)
813 {
814 srcCentroids.setSize(nSrc);
815 }
816 tgtWeights.setSize(nTgt);
817 if (hasTgtMagSf)
818 {
819 tgtMagSf.setSize(nTgt);
820 }
821
822
823 // Append individual data
824 nSrc = 0;
825 nTgt = 0;
826 forAll(allInterfaces, inti)
827 {
828 if (allInterfaces.set(inti))
829 {
831 (
832 allInterfaces[inti]
833 );
834
835 if (!intf.amiPtr_)
836 {
837 continue;
838 }
839
840 const auto& AMI = intf.AMI();
841
842 const auto& srcA = AMI.srcAddress();
843 {
844 // weights
845 SubList<scalarList>(srcWeights, srcA.size(), nSrc) =
846 AMI.srcWeights();
847
848 // magSf
849 if (hasSrcMagSf)
850 {
851 SubList<scalar>(srcMagSf, srcA.size(), nSrc) =
852 AMI.srcMagSf();
853 }
854
855 // centroids
856 if (hasSrcCentroids)
857 {
858 SubList<pointList>(srcCentroids, srcA.size(), nSrc) =
859 AMI.srcCentroids();
860 }
861 }
862
863 const auto& tgtA = AMI.tgtAddress();
864 {
865 // weights
866 SubList<scalarList>(tgtWeights, tgtA.size(), nTgt) =
867 AMI.tgtWeights();
868
869 if (hasTgtMagSf)
870 {
871 SubList<scalar>(tgtMagSf, tgtA.size(), nTgt) =
872 AMI.tgtMagSf();
873 }
874 }
875
876 nSrc += srcA.size();
877 nTgt += tgtA.size();
878 }
879 }
880
881
882 // Construct with same arguments as original
883 amiPtr_.reset
884 (
886 (
887 requireMatch,
888 reverseTarget,
889 lowWeightCorrection
890 )
891 );
892 amiPtr_().comm(coarseComm),
893 amiPtr_().reset
894 (
895 std::move(srcToTgtMap),
896 std::move(tgtToSrcMap),
897 std::move(srcAddress),
898 std::move(srcWeights),
899 std::move(tgtAddress),
900 std::move(tgtWeights),
901 singlePatchProc
902 );
903 amiPtr_().srcMagSf() = std::move(srcMagSf);
904 amiPtr_().srcCentroids() = std::move(srcCentroids);
905 amiPtr_().tgtMagSf() = std::move(tgtMagSf);
906
907 //Pout<< "** constructed new ami:"
908 // << " comm:" << amiPtr_().comm()
909 // << " srcMap.comm:" << amiPtr_().srcMap().comm()
910 // << " tgtMap.comm:" << amiPtr_().tgtMap().comm()
911 // << endl;
912}
913
914
915// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
916
918(
919 const Pstream::commsTypes commsType,
920 const labelUList& iF
921) const
922{
923 const cyclicACMIGAMGInterface& nbr =
924 dynamic_cast<const cyclicACMIGAMGInterface&>(neighbPatch());
925 const labelUList& nbrFaceCells = nbr.faceCells();
926
927 auto tpnf = tmp<labelField>::New(nbrFaceCells.size());
928 labelField& pnf = tpnf.ref();
929
930 forAll(pnf, facei)
931 {
932 pnf[facei] = iF[nbrFaceCells[facei]];
933 }
934
935 return tpnf;
936}
937
938
939void Foam::cyclicACMIGAMGInterface::write(Ostream& os) const
940{
942
943 const bool hasAMI = bool(amiPtr_);
944
945 os << token::SPACE << neighbPatchID_
946 << token::SPACE << owner_
947 << token::SPACE << forwardT_
948 << token::SPACE << reverseT_
949 << token::SPACE << hasAMI;
950
951 if (hasAMI)
952 {
953 os << token::SPACE;
954 AMI().writeData(os);
955
956 // Write processors in communicator
957 const label comm = AMI().comm();
958
959 if (comm != -1)
960 {
962 << UPstream::myProcNo(comm);
963 }
964 }
965}
966
967
968// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label comm() const noexcept
Communicator (local or otherwise) for parallel operations.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
Abstract base class for GAMG agglomerated interfaces.
labelList faceRestrictAddressing_
Face restrict addressing.
virtual label index() const
virtual const lduInterfacePtrsList & coarseInterfaces() const
virtual const labelUList & faceCells() const
Return faceCell addressing.
labelList faceCells_
Face-cell addressing.
GAMGInterface(const GAMGInterface &)=delete
No copy construct.
virtual void write(Ostream &) const =0
Write to stream.
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static void combineReduce(T &value, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors.
static List< T > listGatherValues(const T &localValue, const int communicator=UPstream::worldComm, const int tag=UPstream::msgType())
Gather individual values into list locations.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
commsTypes
Communications types.
Definition UPstream.H:81
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void clear() noexcept
Same as reset(nullptr).
Definition autoPtr.H:255
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
GAMG agglomerated cyclic ACMI interface.
virtual const cyclicACMIGAMGInterface & neighbPatch() const
Return processor number.
virtual void write(Ostream &) const
Write to stream.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual label neighbPatchID() const
Cyclic interface functions.
virtual const AMIPatchToPatchInterpolation & AMI() const
virtual label myProcNo() const
-1 or old local rank
virtual tmp< labelField > internalFieldTransfer(const Pstream::commsTypes commsType, const labelUList &iF) const
Transfer and return internal field adjacent to the interface.
virtual const tensorField & forwardT() const
Return face transformation tensor.
An abstract base class for cyclic ACMI coupled interfaces.
cyclicACMILduInterface() noexcept=default
Default construct.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
An abstract base class for implicitly-coupled interfaces e.g. processor and cyclic patches.
static label renumberMap(labelListList &mapElements, const labelUList &oldToNew, const bool hasFlip)
Helper for renumbering the (compacted) map elements using the supplied old-to-new mapping.
Class containing processor-to-processor mapping information.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
@ SPACE
Space [isspace].
Definition token.H:144
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define DebugPout
Report an information message using Foam::Pout.
List< bool > bools(const labelHashSet &locations)
Transform the on locations to a boolList, with true for each non-negative location and false for all ...
Definition HashOps.C:72
label find(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as ListOps::find_if.
Definition ListOps.H:855
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition label.H:63
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
List< pointList > pointListList
List of pointList.
Definition pointList.H:35
bool readBool(Istream &is)
Read bool from stream using Foam::Switch(Istream&).
Definition bool.C:72
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299