FaceCellWave.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2021 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 "FaceCellWave.H"
30#include "polyMesh.H"
31#include "processorPolyPatch.H"
32#include "cyclicPolyPatch.H"
33#include "cyclicAMIPolyPatch.H"
34#include "UIPstream.H"
35#include "UOPstream.H"
36#include "PstreamReduceOps.H"
37#include "debug.H"
38#include "typeInfo.H"
39#include "SubField.H"
40#include "globalMeshData.H"
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
44template<class Type, class TrackingData>
46
47template<class Type, class TrackingData>
49
50template<class Type, class TrackingData>
52
53
54namespace Foam
55{
56 template<class Type, class TrackingData>
57 class combine
58 {
59 //- Combine operator for AMIInterpolation
60
62
63 const cyclicAMIPolyPatch& patch_;
64
65 public:
66
68 (
70 const cyclicAMIPolyPatch& patch
71 )
72 :
73 solver_(solver),
74 patch_(patch)
75 {}
76
77
78 void operator()
79 (
80 Type& x,
81 const label facei,
82 const Type& y,
83 const scalar weight
84 ) const
85 {
86 if (y.valid(solver_.data()))
87 {
88 label meshFacei = -1;
89 if (patch_.owner())
90 {
91 meshFacei = patch_.start() + facei;
92 }
93 else
94 {
95 meshFacei = patch_.neighbPatch().start() + facei;
96 }
97 x.updateFace
98 (
99 solver_.mesh(),
100 meshFacei,
101 y,
102 solver_.propagationTol(),
103 solver_.data()
104 );
105 }
106 }
107 };
108}
109
110
111// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
112
113template<class Type, class TrackingData>
115(
116 const label celli,
117 const label neighbourFacei,
118 const Type& neighbourInfo,
119 const scalar tol,
120 Type& cellInfo
121)
122{
123 // Update info for celli, at position pt, with information from
124 // neighbouring face/cell.
125 // Updates:
126 // - changedCell_, changedCells_
127 // - statistics: nEvals_, nUnvisitedCells_
128
129 ++nEvals_;
130
131 const bool wasValid = cellInfo.valid(td_);
132
133 const bool propagate =
135 (
136 mesh_,
137 celli,
138 neighbourFacei,
139 neighbourInfo,
140 tol,
141 td_
142 );
143
144 if (propagate)
145 {
146 if (changedCell_.set(celli))
147 {
148 changedCells_.append(celli);
150 }
151
152 if (!wasValid && cellInfo.valid(td_))
153 {
154 --nUnvisitedCells_;
155 }
156
157 return propagate;
158}
159
161template<class Type, class TrackingData>
163(
164 const label facei,
165 const label neighbourCelli,
166 const Type& neighbourInfo,
167 const scalar tol,
168 Type& faceInfo
169)
170{
171 // Update info for facei, at position pt, with information from
172 // neighbouring face/cell.
173 // Updates:
174 // - changedFace_, changedFaces_,
175 // - statistics: nEvals_, nUnvisitedFaces_
176
177 ++nEvals_;
178
179 const bool wasValid = faceInfo.valid(td_);
180
181 const bool propagate =
182 faceInfo.updateFace
184 mesh_,
185 facei,
186 neighbourCelli,
187 neighbourInfo,
188 tol,
189 td_
190 );
191
192 if (propagate)
193 {
194 if (changedFace_.set(facei))
195 {
196 changedFaces_.append(facei);
197 }
198 }
200 if (!wasValid && faceInfo.valid(td_))
201 {
202 --nUnvisitedFaces_;
203 }
204
205 return propagate;
206}
207
208
209template<class Type, class TrackingData>
211(
212 const label facei,
213 const Type& neighbourInfo,
214 const scalar tol,
215 Type& faceInfo
216)
217{
218 // Update info for facei, at position pt, with information from
219 // same face.
220 // Updates:
221 // - changedFace_, changedFaces_,
222 // - statistics: nEvals_, nUnvisitedFaces_
223
224 ++nEvals_;
225
226 const bool wasValid = faceInfo.valid(td_);
228 const bool propagate =
229 faceInfo.updateFace
230 (
231 mesh_,
232 facei,
233 neighbourInfo,
234 tol,
235 td_
236 );
237
238 if (propagate)
239 {
240 if (changedFace_.set(facei))
241 {
242 changedFaces_.append(facei);
243 }
244 }
246 if (!wasValid && faceInfo.valid(td_))
247 {
248 --nUnvisitedFaces_;
250
251 return propagate;
253
254
255template<class Type, class TrackingData>
257(
258 const polyPatch& patch
259) const
260{
261 // For debugging: check status on both sides of cyclic
262
263 const cyclicPolyPatch& nbrPatch =
264 refCast<const cyclicPolyPatch>(patch).neighbPatch();
265
266 forAll(patch, patchFacei)
267 {
268 const label i1 = patch.start() + patchFacei;
269 const label i2 = nbrPatch.start() + patchFacei;
270
271 if
272 (
273 !allFaceInfo_[i1].sameGeometry
274 (
275 mesh_,
276 allFaceInfo_[i2],
277 geomTol_,
278 td_
279 )
280 )
281 {
283 << " faceInfo:" << allFaceInfo_[i1]
284 << " otherfaceInfo:" << allFaceInfo_[i2]
285 << abort(FatalError);
286 }
288 if (changedFace_.test(i1) != changedFace_.test(i2))
289 {
291 << " faceInfo:" << allFaceInfo_[i1]
292 << " otherfaceInfo:" << allFaceInfo_[i2]
293 << " changedFace:" << changedFace_.test(i1)
294 << " otherchangedFace:" << changedFace_.test(i2)
295 << abort(FatalError);
296 }
297 }
299
300
301template<class Type, class TrackingData>
302template<class PatchType>
304{
305 for (const polyPatch& p : mesh_.boundaryMesh())
306 {
307 if (isA<PatchType>(p))
308 {
309 return true;
310 }
311 }
312 return false;
314
315
316template<class Type, class TrackingData>
318(
319 const label facei,
320 const Type& faceInfo
321)
322{
323 const bool wasValid = allFaceInfo_[facei].valid(td_);
324
325 // Copy info for facei
326 allFaceInfo_[facei] = faceInfo;
327
328 // Maintain count of unset faces
329 if (!wasValid && allFaceInfo_[facei].valid(td_))
330 {
331 --nUnvisitedFaces_;
332 }
333
334 // Mark facei as visited and changed (both on list and on face itself)
335 changedFace_.set(facei);
336 changedFaces_.append(facei);
337}
338
339
340template<class Type, class TrackingData>
342(
343 const labelUList& changedFaces,
344 const List<Type>& changedFacesInfo
345)
346{
347 forAll(changedFaces, changedFacei)
348 {
349 const label facei = changedFaces[changedFacei];
350
351 const bool wasValid = allFaceInfo_[facei].valid(td_);
352
353 // Copy info for facei
354 allFaceInfo_[facei] = changedFacesInfo[changedFacei];
355
356 // Maintain count of unset faces
357 if (!wasValid && allFaceInfo_[facei].valid(td_))
358 {
359 --nUnvisitedFaces_;
360 }
361
362 // Mark facei as changed, both on list and on face itself.
363 changedFace_.set(facei);
364 changedFaces_.append(facei);
366}
367
369template<class Type, class TrackingData>
371(
372 const polyPatch& patch,
373 const label nFaces,
374 const labelUList& changedFaces,
375 const List<Type>& changedFacesInfo
376)
377{
378 // Merge face information into member data
379
380 for (label changedFacei = 0; changedFacei < nFaces; ++changedFacei)
381 {
382 const Type& newInfo = changedFacesInfo[changedFacei];
383 const label patchFacei = changedFaces[changedFacei];
384
385 const label meshFacei = patch.start() + patchFacei;
387 Type& currInfo = allFaceInfo_[meshFacei];
388
389 if (!currInfo.equal(newInfo, td_))
390 {
391 updateFace
392 (
393 meshFacei,
394 newInfo,
395 propagationTol_,
396 currInfo
397 );
398 }
399 }
400}
401
402
403template<class Type, class TrackingData>
405(
406 const polyPatch& patch,
407 const label startFacei,
408 const label nFaces,
409 labelList& changedPatchFaces,
410 List<Type>& changedPatchFacesInfo
411) const
412{
413 // Construct compact patchFace change arrays for a (slice of a) single
414 // patch. changedPatchFaces in local patch numbering.
415 // Return length of arrays.
416 label nChanged = 0;
417
418 for (label i = 0; i < nFaces; ++i)
419 {
420 const label patchFacei = i + startFacei;
421 const label meshFacei = patch.start() + patchFacei;
422
423 if (changedFace_.test(meshFacei))
424 {
425 changedPatchFaces[nChanged] = patchFacei;
426 changedPatchFacesInfo[nChanged] = allFaceInfo_[meshFacei];
427 ++nChanged;
428 }
429 }
430
431 return nChanged;
432}
433
434
435template<class Type, class TrackingData>
437(
438 const polyPatch& patch,
439 const label nFaces,
440 const labelUList& faceLabels,
441 List<Type>& faceInfo
442) const
443{
444 // Handle leaving domain. Implementation referred to Type
445
446 const vectorField& fc = mesh_.faceCentres();
447
448 for (label i = 0; i < nFaces; ++i)
449 {
450 const label patchFacei = faceLabels[i];
451 const label meshFacei = patch.start() + patchFacei;
452
453 faceInfo[i].leaveDomain(mesh_, patch, patchFacei, fc[meshFacei], td_);
454 }
455}
456
457
458template<class Type, class TrackingData>
460(
461 const polyPatch& patch,
462 const label nFaces,
463 const labelUList& faceLabels,
464 List<Type>& faceInfo
465) const
466{
467 // Handle entering domain. Implementation referred to Type
468
469 const vectorField& fc = mesh_.faceCentres();
470
471 for (label i = 0; i < nFaces; ++i)
472 {
473 const label patchFacei = faceLabels[i];
474 const label meshFacei = patch.start() + patchFacei;
475
476 faceInfo[i].enterDomain(mesh_, patch, patchFacei, fc[meshFacei], td_);
477 }
478}
479
480
481template<class Type, class TrackingData>
483(
484 const tensorField& rotTensor,
485 const label nFaces,
486 List<Type>& faceInfo
487)
488{
489 // Transform. Implementation referred to Type
490
491 if (rotTensor.size() == 1)
492 {
493 const tensor& T = rotTensor[0];
494
495 for (label facei = 0; facei < nFaces; ++facei)
496 {
497 faceInfo[facei].transform(mesh_, T, td_);
498 }
499 }
500 else
501 {
502 for (label facei = 0; facei < nFaces; ++facei)
503 {
504 faceInfo[facei].transform(mesh_, rotTensor[facei], td_);
505 }
506 }
507}
508
509
510template<class Type, class TrackingData>
512(
513 const polyPatch&,
514 const label cycOffset,
515 const label nFaces,
516 labelList& faces
517)
518{
519 // Offset mesh face.
520 // Used for transferring from one cyclic half to the other.
521
522 for (label facei = 0; facei < nFaces; ++facei)
523 {
524 faces[facei] += cycOffset;
525 }
526}
527
528
529template<class Type, class TrackingData>
531{
532 // Transfer all the information to/from neighbouring processors
533
534 const globalMeshData& pData = mesh_.globalData();
535
536 // Which patches are processor patches
537 const labelList& procPatches = pData.processorPatches();
538
539 // Send all
540
541 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
542
543 for (const label patchi : procPatches)
544 {
545 const processorPolyPatch& procPatch =
546 refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
547
548 // Allocate buffers
549 label nSendFaces;
550 labelList sendFaces(procPatch.size());
551 List<Type> sendFacesInfo(procPatch.size());
552
553 // Determine which faces changed on current patch
554 nSendFaces = getChangedPatchFaces
555 (
556 procPatch,
557 0,
558 procPatch.size(),
559 sendFaces,
560 sendFacesInfo
561 );
562
563 // Adapt wallInfo for leaving domain
564 leaveDomain
565 (
566 procPatch,
567 nSendFaces,
568 sendFaces,
569 sendFacesInfo
570 );
571
572 if (debug & 2)
573 {
574 Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
575 << " communicating with " << procPatch.neighbProcNo()
576 << " Sending:" << nSendFaces
577 << endl;
578 }
579
580 UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
581 //writeFaces(nSendFaces, sendFaces, sendFacesInfo, toNeighbour);
582 toNeighbour
583 << SubList<label>(sendFaces, nSendFaces)
584 << SubList<Type>(sendFacesInfo, nSendFaces);
585 }
586
587 pBufs.finishedSends();
588
589 // Receive all
590
591 for (const label patchi : procPatches)
592 {
593 const processorPolyPatch& procPatch =
594 refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
595
596 // Allocate buffers
597 labelList receiveFaces;
598 List<Type> receiveFacesInfo;
599
600 {
601 UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
602 fromNeighbour >> receiveFaces >> receiveFacesInfo;
603 }
604
605 if (debug & 2)
606 {
607 Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
608 << " communicating with " << procPatch.neighbProcNo()
609 << " Receiving:" << receiveFaces.size()
610 << endl;
611 }
612
613 // Apply transform to received data for non-parallel planes
614 if (!procPatch.parallel())
615 {
617 (
618 procPatch.forwardT(),
619 receiveFaces.size(),
620 receiveFacesInfo
621 );
622 }
623
624 // Adapt wallInfo for entering domain
625 enterDomain
626 (
627 procPatch,
628 receiveFaces.size(),
629 receiveFaces,
630 receiveFacesInfo
631 );
632
633 // Merge received info
634 mergeFaceInfo
635 (
636 procPatch,
637 receiveFaces.size(),
638 receiveFaces,
639 receiveFacesInfo
640 );
641 }
642}
643
644
645template<class Type, class TrackingData>
647{
648 // Transfer information across cyclic halves.
649
650 for (const polyPatch& patch : mesh_.boundaryMesh())
651 {
652 const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(patch);
653
654 if (cpp)
655 {
656 const auto& cycPatch = *cpp;
657 const auto& nbrPatch = cycPatch.neighbPatch();
658
659 // Allocate buffers
660 label nReceiveFaces;
661 labelList receiveFaces(patch.size());
662 List<Type> receiveFacesInfo(patch.size());
663
664 // Determine which faces changed
665 nReceiveFaces = getChangedPatchFaces
666 (
667 nbrPatch,
668 0,
669 nbrPatch.size(),
670 receiveFaces,
671 receiveFacesInfo
672 );
673
674 // Adapt wallInfo for leaving domain
675 leaveDomain
676 (
677 nbrPatch,
678 nReceiveFaces,
679 receiveFaces,
680 receiveFacesInfo
681 );
682
683 if (!cycPatch.parallel())
684 {
685 // Received data from other half
687 (
688 cycPatch.forwardT(),
689 nReceiveFaces,
690 receiveFacesInfo
691 );
692 }
693
694 if (debug & 2)
695 {
696 Pout<< " Cyclic patch "
697 << cycPatch.index() << ' ' << cycPatch.name()
698 << " Changed : " << nReceiveFaces
699 << endl;
700 }
701
702 // Half2: Adapt wallInfo for entering domain
703 enterDomain
704 (
705 cycPatch,
706 nReceiveFaces,
707 receiveFaces,
708 receiveFacesInfo
709 );
710
711 // Merge into global storage
712 mergeFaceInfo
713 (
714 cycPatch,
715 nReceiveFaces,
716 receiveFaces,
717 receiveFacesInfo
718 );
719
720 if (debug)
721 {
722 checkCyclic(cycPatch);
723 }
724 }
725 }
726}
727
728
729template<class Type, class TrackingData>
731{
732 // Transfer information across cyclicAMI halves.
733
734 for (const polyPatch& patch : mesh_.boundaryMesh())
735 {
736 const cyclicAMIPolyPatch* cpp = isA<cyclicAMIPolyPatch>(patch);
737
738 if (cpp)
739 {
740 const auto& cycPatch = *cpp;
741 const auto& nbrPatch = cycPatch.neighbPatch();
742
743 List<Type> receiveInfo;
744
745 {
746 // Get nbrPatch data (so not just changed faces)
747 typename List<Type>::subList sendInfo
748 (
749 nbrPatch.patchSlice
750 (
751 allFaceInfo_
752 )
753 );
754
755 if (!nbrPatch.parallel() || nbrPatch.separated())
756 {
757 // Adapt sendInfo for leaving domain
758 const vectorField::subField fc = nbrPatch.faceCentres();
759 forAll(sendInfo, i)
760 {
761 sendInfo[i].leaveDomain(mesh_, nbrPatch, i, fc[i], td_);
762 }
763 }
764
765 // Transfer sendInfo to cycPatch
766 combine<Type, TrackingData> cmb(*this, cycPatch);
767
768 if (cycPatch.applyLowWeightCorrection())
769 {
770 List<Type> defVals
771 (
772 cycPatch.patchInternalList(allCellInfo_)
773 );
774
775 cycPatch.interpolate(sendInfo, cmb, receiveInfo, defVals);
776 }
777 else
778 {
779 cycPatch.interpolate(sendInfo, cmb, receiveInfo);
780 }
781 }
782
783 // Apply transform to received data for non-parallel planes
784 if (!cycPatch.parallel())
785 {
787 (
788 cycPatch.forwardT(),
789 receiveInfo.size(),
790 receiveInfo
791 );
792 }
793
794 if (!cycPatch.parallel() || cycPatch.separated())
795 {
796 // Adapt receiveInfo for entering domain
797 const vectorField::subField fc = cycPatch.faceCentres();
798 forAll(receiveInfo, i)
799 {
800 receiveInfo[i].enterDomain(mesh_, cycPatch, i, fc[i], td_);
801 }
802 }
803
804 // Merge into global storage
805 forAll(receiveInfo, i)
806 {
807 const label meshFacei = cycPatch.start()+i;
808
809 const Type& newInfo = receiveInfo[i];
810
811 Type& currInfo = allFaceInfo_[meshFacei];
812
813 if (newInfo.valid(td_) && !currInfo.equal(newInfo, td_))
814 {
815 updateFace
816 (
817 meshFacei,
818 newInfo,
819 propagationTol_,
820 currInfo
821 );
822 }
823 }
824 }
825 }
826}
827
828
829template<class Type, class TrackingData>
831{
832 changedBaffles_.clear();
833
834 // Collect all/any changed information touching a baffle
835 for (const labelPair& baffle : explicitConnections_)
836 {
837 const label f0 = baffle.first();
838 const label f1 = baffle.second();
839
840 if (changedFace_.test(f0))
841 {
842 // f0 changed. Update information on f1.
843 changedBaffles_.append(taggedInfoType(f1, allFaceInfo_[f0]));
844 }
845
846 if (changedFace_.test(f1))
847 {
848 // f1 changed. Update information on f0.
849 changedBaffles_.append(taggedInfoType(f0, allFaceInfo_[f1]));
850 }
851 }
852
853
854 // Update other side with changed information
855
856 for (const taggedInfoType& updated : changedBaffles_)
857 {
858 const label tgtFace = updated.first;
859 const Type& newInfo = updated.second;
860
861 Type& currInfo = allFaceInfo_[tgtFace];
862
863 if (!currInfo.equal(newInfo, td_))
864 {
865 updateFace
866 (
867 tgtFace,
868 newInfo,
869 propagationTol_,
870 currInfo
871 );
872 }
873 }
874
875 changedBaffles_.clear();
876}
877
878
879// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
880
881template<class Type, class TrackingData>
883(
884 const polyMesh& mesh,
885 UList<Type>& allFaceInfo,
886 UList<Type>& allCellInfo,
887 TrackingData& td
888)
889:
890 mesh_(mesh),
891 explicitConnections_(),
892 allFaceInfo_(allFaceInfo),
893 allCellInfo_(allCellInfo),
894 td_(td),
895 changedFace_(mesh_.nFaces(), false),
896 changedCell_(mesh_.nCells(), false),
897 changedFaces_(mesh_.nFaces()),
898 changedCells_(mesh_.nCells()),
899 changedBaffles_(2*explicitConnections_.size()),
900 hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
901 hasCyclicAMIPatches_
902 (
904 ),
905 nEvals_(0),
906 nUnvisitedCells_(mesh_.nCells()),
907 nUnvisitedFaces_(mesh_.nFaces())
908{
909 if
910 (
912 || allCellInfo.size() != mesh_.nCells()
913 )
914 {
916 << "face and cell storage not the size of mesh faces, cells:" << nl
917 << " allFaceInfo :" << allFaceInfo.size() << nl
918 << " mesh_.nFaces():" << mesh_.nFaces() << nl
919 << " allCellInfo :" << allCellInfo.size() << nl
920 << " mesh_.nCells():" << mesh_.nCells() << endl
921 << exit(FatalError);
922 }
923}
924
925
926template<class Type, class TrackingData>
928(
929 const polyMesh& mesh,
930 const labelUList& changedFaces,
931 const List<Type>& changedFacesInfo,
932 UList<Type>& allFaceInfo,
933 UList<Type>& allCellInfo,
934 const label maxIter,
935 TrackingData& td
936)
937:
938 mesh_(mesh),
939 explicitConnections_(),
940 allFaceInfo_(allFaceInfo),
941 allCellInfo_(allCellInfo),
942 td_(td),
943 changedFace_(mesh_.nFaces(), false),
944 changedCell_(mesh_.nCells(), false),
945 changedFaces_(mesh_.nFaces()),
946 changedCells_(mesh_.nCells()),
947 changedBaffles_(2*explicitConnections_.size()),
948 hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
949 hasCyclicAMIPatches_
950 (
952 ),
953 nEvals_(0),
954 nUnvisitedCells_(mesh_.nCells()),
955 nUnvisitedFaces_(mesh_.nFaces())
956{
957 if
958 (
960 || allCellInfo.size() != mesh_.nCells()
961 )
962 {
964 << "face and cell storage not the size of mesh faces, cells:" << nl
965 << " allFaceInfo :" << allFaceInfo.size() << nl
966 << " mesh_.nFaces():" << mesh_.nFaces() << nl
967 << " allCellInfo :" << allCellInfo.size() << nl
968 << " mesh_.nCells():" << mesh_.nCells() << endl
969 << exit(FatalError);
970 }
971
972 // Copy initial changed faces data
973 setFaceInfo(changedFaces, changedFacesInfo);
974
975 // Iterate until nothing changes
976 const label iter = iterate(maxIter);
977
978 if ((maxIter > 0) && (iter >= maxIter))
979 {
981 << "Maximum number of iterations reached. Increase maxIter." << nl
982 << " maxIter:" << maxIter << nl
983 << " nChangedCells:" << changedCells_.size() << nl
984 << " nChangedFaces:" << changedFaces_.size() << endl
985 << exit(FatalError);
986 }
987}
988
989
990template<class Type, class TrackingData>
992(
993 const polyMesh& mesh,
994 const labelPairList& explicitConnections,
995 const bool handleCyclicAMI,
996 const labelUList& changedFaces,
997 const List<Type>& changedFacesInfo,
998 UList<Type>& allFaceInfo,
999 UList<Type>& allCellInfo,
1000 const label maxIter,
1001 TrackingData& td
1002)
1003:
1004 mesh_(mesh),
1005 explicitConnections_(explicitConnections),
1006 allFaceInfo_(allFaceInfo),
1007 allCellInfo_(allCellInfo),
1008 td_(td),
1009 changedFace_(mesh_.nFaces(), false),
1010 changedCell_(mesh_.nCells(), false),
1011 changedFaces_(mesh_.nFaces()),
1012 changedCells_(mesh_.nCells()),
1013 changedBaffles_(2*explicitConnections_.size()),
1014 hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
1015 hasCyclicAMIPatches_
1016 (
1017 handleCyclicAMI
1018 && returnReduce(hasPatch<cyclicAMIPolyPatch>(), orOp<bool>())
1019 ),
1020 nEvals_(0),
1021 nUnvisitedCells_(mesh_.nCells()),
1022 nUnvisitedFaces_(mesh_.nFaces())
1023{
1024 if
1025 (
1027 || allCellInfo.size() != mesh_.nCells()
1028 )
1029 {
1031 << "face and cell storage not the size of mesh faces, cells:" << nl
1032 << " allFaceInfo :" << allFaceInfo.size() << nl
1033 << " mesh_.nFaces():" << mesh_.nFaces() << nl
1034 << " allCellInfo :" << allCellInfo.size() << nl
1035 << " mesh_.nCells():" << mesh_.nCells() << endl
1036 << exit(FatalError);
1037 }
1038
1039 // Copy initial changed faces data
1040 setFaceInfo(changedFaces, changedFacesInfo);
1041
1042 // Iterate until nothing changes
1043 const label iter = iterate(maxIter);
1044
1045 if ((maxIter > 0) && (iter >= maxIter))
1046 {
1048 << "Maximum number of iterations reached. Increase maxIter." << nl
1049 << " maxIter:" << maxIter << nl
1050 << " nChangedCells:" << changedCells_.size() << nl
1051 << " nChangedFaces:" << changedFaces_.size() << endl
1052 << exit(FatalError);
1053 }
1054}
1055
1056
1057// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1058
1059template<class Type, class TrackingData>
1061{
1062 return nUnvisitedCells_;
1063}
1064
1065
1066template<class Type, class TrackingData>
1068{
1069 return nUnvisitedFaces_;
1070}
1071
1072
1073template<class Type, class TrackingData>
1075{
1076 // Propagate face to cell
1077
1078 const labelList& owner = mesh_.faceOwner();
1079 const labelList& neighbour = mesh_.faceNeighbour();
1080 const label nInternalFaces = mesh_.nInternalFaces();
1081
1082 for (const label facei : changedFaces_)
1083 {
1084 if (!changedFace_.test(facei))
1085 {
1087 << "Face " << facei
1088 << " not marked as having been changed"
1089 << abort(FatalError);
1090 }
1091
1092 const Type& newInfo = allFaceInfo_[facei];
1093
1094 // Evaluate all connected cells
1095
1096 // Owner
1097 {
1098 const label celli = owner[facei];
1099 Type& currInfo = allCellInfo_[celli];
1100
1101 if (!currInfo.equal(newInfo, td_))
1102 {
1103 updateCell
1104 (
1105 celli,
1106 facei,
1107 newInfo,
1108 propagationTol_,
1109 currInfo
1110 );
1111 }
1112 }
1113
1114 // Neighbour.
1115 if (facei < nInternalFaces)
1116 {
1117 const label celli = neighbour[facei];
1118 Type& currInfo = allCellInfo_[celli];
1119
1120 if (!currInfo.equal(newInfo, td_))
1121 {
1122 updateCell
1123 (
1124 celli,
1125 facei,
1126 newInfo,
1127 propagationTol_,
1128 currInfo
1129 );
1130 }
1131 }
1132
1133 // Reset status of face
1134 changedFace_.unset(facei);
1135 }
1136
1137 // Handled all changed faces by now
1138 changedFaces_.clear();
1139
1140 if (debug & 2)
1141 {
1142 Pout<< " Changed cells : " << changedCells_.size() << endl;
1143 }
1144
1145 // Number of changedCells over all procs
1146 return returnReduce(changedCells_.size(), sumOp<label>());
1147}
1148
1149
1150template<class Type, class TrackingData>
1152{
1153 // Propagate cell to face
1154
1155 const cellList& cells = mesh_.cells();
1156
1157 for (const label celli : changedCells_)
1158 {
1159 if (!changedCell_.test(celli))
1160 {
1162 << "Cell " << celli << " not marked as having been changed"
1163 << abort(FatalError);
1164 }
1165
1166 const Type& newInfo = allCellInfo_[celli];
1167
1168 // Evaluate all connected faces
1169
1170 const labelList& faceLabels = cells[celli];
1171 for (const label facei : faceLabels)
1172 {
1173 Type& currInfo = allFaceInfo_[facei];
1174
1175 if (!currInfo.equal(newInfo, td_))
1176 {
1177 updateFace
1178 (
1179 facei,
1180 celli,
1181 newInfo,
1182 propagationTol_,
1183 currInfo
1184 );
1185 }
1186 }
1187
1188 // Reset status of cell
1189 changedCell_.unset(celli);
1190 }
1191
1192 // Handled all changed cells by now
1193 changedCells_.clear();
1194
1195
1196 // Transfer across any explicitly provided internal connections
1197 handleExplicitConnections();
1198
1199 if (hasCyclicPatches_)
1200 {
1201 handleCyclicPatches();
1202 }
1203
1204 if (hasCyclicAMIPatches_)
1205 {
1206 handleAMICyclicPatches();
1207 }
1208
1209 if (Pstream::parRun())
1210 {
1211 handleProcPatches();
1212 }
1213
1214 if (debug & 2)
1215 {
1216 Pout<< " Changed faces : " << changedFaces_.size() << endl;
1217 }
1218
1219
1220 // Number of changedFaces over all procs
1221 return returnReduce(changedFaces_.size(), sumOp<label>());
1222}
1223
1224
1225// Iterate
1226template<class Type, class TrackingData>
1228{
1229 if (maxIter < 0)
1230 {
1231 return 0;
1232 }
1233
1234 if (hasCyclicPatches_)
1235 {
1236 handleCyclicPatches();
1237 }
1238
1239 if (hasCyclicAMIPatches_)
1240 {
1241 handleAMICyclicPatches();
1242 }
1243
1244 if (Pstream::parRun())
1245 {
1246 handleProcPatches();
1247 }
1248
1249 label iter = 0;
1250
1251 for (/*nil*/; iter < maxIter; ++iter)
1252 {
1253 if (debug)
1254 {
1255 Info<< " Iteration " << iter << endl;
1256 }
1257
1258 nEvals_ = 0;
1259 const label nCells = faceToCell();
1260 const label nFaces = nCells ? cellToFace() : 0;
1261
1262 if (debug)
1263 {
1264 Info<< " Total evaluations : "
1265 << nEvals_ << nl
1266 << " Changed cells / faces : "
1267 << nCells << " / " << nFaces << nl
1268 << " Pending cells / faces : "
1269 << nUnvisitedCells_ << " / " << nUnvisitedFaces_ << nl;
1270 }
1271
1272 if (!nCells || !nFaces)
1273 {
1274 break;
1275 }
1276 }
1277
1278 return iter;
1279}
1280
1281
1282// ************************************************************************* //
scalar y
Inter-processor communication reduction functions.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:81
DynamicList< label > changedCells_
Definition: FaceCellWave.H:125
void handleExplicitConnections()
Merge data across explicitly provided local connections.
Definition: FaceCellWave.C:830
DynamicList< label > changedFaces_
List of changed faces.
Definition: FaceCellWave.H:122
void handleProcPatches()
Merge data from across processor boundaries.
Definition: FaceCellWave.C:530
void leaveDomain(const polyPatch &patch, const label nFaces, const labelUList &faceLabels, List< Type > &faceInfo) const
Handle leaving domain. Implementation referred to Type.
Definition: FaceCellWave.C:437
const polyMesh & mesh() const
Access mesh.
Definition: FaceCellWave.H:354
static scalar propagationTol()
Access to tolerance.
Definition: FaceCellWave.H:271
label nUnvisitedCells() const
void mergeFaceInfo(const polyPatch &patch, const label nFaces, const labelUList &changedFaces, const List< Type > &changedFacesInfo)
Merge received patch data into global data.
Definition: FaceCellWave.C:371
bool hasPatch() const
Has cyclic patch?
Definition: FaceCellWave.C:303
void handleAMICyclicPatches()
Merge data from across AMI cyclics.
Definition: FaceCellWave.C:730
void checkCyclic(const polyPatch &pPatch) const
Debugging: check info on both sides of cyclic.
Definition: FaceCellWave.C:257
void enterDomain(const polyPatch &patch, const label nFaces, const labelUList &faceLabels, List< Type > &faceInfo) const
Handle leaving domain. Implementation referred to Type.
Definition: FaceCellWave.C:460
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached.
void handleCyclicPatches()
Merge data from across cyclics.
Definition: FaceCellWave.C:646
void setFaceInfo(const label facei, const Type &faceInfo)
Set single initial changed face.
Definition: FaceCellWave.C:318
bool updateCell(const label celli, const label neighbourFacei, const Type &neighbourInfo, const scalar tol, Type &cellInfo)
Updates cellInfo with information from neighbour.
Definition: FaceCellWave.C:115
label getChangedPatchFaces(const polyPatch &patch, const label startFacei, const label nFaces, labelList &changedPatchFaces, List< Type > &changedPatchFacesInfo) const
Extract info for single patch only.
Definition: FaceCellWave.C:405
virtual label faceToCell()
Propagate from face to cell.
UList< Type > & allCellInfo()
Access allCellInfo.
Definition: FaceCellWave.H:342
const polyMesh & mesh_
Reference to mesh.
Definition: FaceCellWave.H:101
bool updateFace(const label facei, const label neighbourCelli, const Type &neighbourInfo, const scalar tol, Type &faceInfo)
Updates faceInfo with information from neighbour.
Definition: FaceCellWave.C:163
std::pair< label, Type > taggedInfoType
Information tagged with a source or destination id.
Definition: FaceCellWave.H:86
virtual label cellToFace()
Propagate from cell to face.
UList< Type > & allFaceInfo()
Access allFaceInfo.
Definition: FaceCellWave.H:336
label nUnvisitedFaces() const
Get number of unvisited faces.
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:348
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool wait=true)
Mark sends as done.
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition: cellInfo.H:65
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const cellInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: cellInfoI.H:173
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Definition: cellInfoI.H:113
A topoSetFaceSource to select all the faces from given cellSet(s).
Definition: cellToFace.H:177
combine(FaceCellWave< Type, TrackingData > &solver, const cyclicAMIPolyPatch &patch)
Definition: FaceCellWave.C:68
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrary Mesh Interface (AMI)
virtual bool owner() const
Does this side own the patch?
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Cyclic plane patch.
const cyclicPolyPatch & neighbPatch() const
A topoSetCellSource to select all cells based on usage in given faceSet(s), e.g. select cells that ar...
Definition: faceToCell.H:182
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
Default transformation behaviour.
const vector & offset() const noexcept
Offset vector (from patch faces to destination mesh objects)
const word & name() const noexcept
The patch name.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const List< T >::subList patchSlice(const UList< T > &l) const
Slice list to patch.
Definition: polyPatch.H:403
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
Neighbour processor patch.
int neighbProcNo() const
Return neighbour processor number.
Base class for solution control classes.
Definition: solver.H:54
volScalarField & p
const volScalarField & T
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const cellShapeList & cells
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333