PointEdgeWave.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) 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 "PointEdgeWave.H"
30#include "polyMesh.H"
31#include "processorPolyPatch.H"
32#include "cyclicPolyPatch.H"
33#include "UIPstream.H"
34#include "UOPstream.H"
35#include "debug.H"
36#include "typeInfo.H"
37#include "globalMeshData.H"
38#include "pointFields.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42template<class Type, class TrackingData>
44
45template<class Type, class TrackingData>
47
48namespace Foam
49{
50 //- Reduction class. If x and y are not equal assign value.
51 template<class Type, class TrackingData>
53 {
54 TrackingData& td_;
55
56 public:
57
58 combineEqOp(TrackingData& td)
59 :
60 td_(td)
61 {}
62
63 void operator()(Type& x, const Type& y) const
64 {
65 if (!x.valid(td_) && y.valid(td_))
66 {
67 x = y;
68 }
69 }
70 };
71}
72
73
74// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
75
76// Handle leaving domain. Implementation referred to Type
77template<class Type, class TrackingData>
79(
80 const polyPatch& patch,
81 const labelList& patchPointLabels,
82 List<Type>& pointInfo
83) const
84{
85 const labelList& meshPoints = patch.meshPoints();
86
87 forAll(patchPointLabels, i)
88 {
89 label patchPointi = patchPointLabels[i];
90
91 const point& pt = patch.points()[meshPoints[patchPointi]];
92
93 pointInfo[i].leaveDomain(patch, patchPointi, pt, td_);
94 }
95}
96
97
98// Handle entering domain. Implementation referred to Type
99template<class Type, class TrackingData>
101(
102 const polyPatch& patch,
103 const labelList& patchPointLabels,
104 List<Type>& pointInfo
105) const
106{
107 const labelList& meshPoints = patch.meshPoints();
108
109 forAll(patchPointLabels, i)
110 {
111 label patchPointi = patchPointLabels[i];
112
113 const point& pt = patch.points()[meshPoints[patchPointi]];
114
115 pointInfo[i].enterDomain(patch, patchPointi, pt, td_);
116 }
117}
118
119
120// Transform. Implementation referred to Type
121template<class Type, class TrackingData>
123(
124 const polyPatch& patch,
125 const tensorField& rotTensor,
126 List<Type>& pointInfo
127) const
128{
129 if (rotTensor.size() == 1)
130 {
131 const tensor& T = rotTensor[0];
132
133 forAll(pointInfo, i)
134 {
135 pointInfo[i].transform(T, td_);
136 }
137 }
138 else
139 {
141 << "Non-uniform transformation on patch " << patch.name()
142 << " of type " << patch.type()
143 << " not supported for point fields"
144 << abort(FatalError);
145
146 forAll(pointInfo, i)
147 {
148 pointInfo[i].transform(rotTensor[i], td_);
149 }
150 }
151}
152
153
154// Update info for pointi, at position pt, with information from
155// neighbouring edge.
156// Updates:
157// - changedPoint_, changedPoints_, nChangedPoints_,
158// - statistics: nEvals_, nUnvisitedPoints_
159template<class Type, class TrackingData>
161(
162 const label pointi,
163 const label neighbourEdgeI,
164 const Type& neighbourInfo,
165 Type& pointInfo
166)
167{
168 nEvals_++;
169
170 bool wasValid = pointInfo.valid(td_);
171
172 bool propagate =
173 pointInfo.updatePoint
174 (
175 mesh_,
176 pointi,
177 neighbourEdgeI,
178 neighbourInfo,
179 propagationTol_,
180 td_
181 );
182
183 if (propagate)
184 {
185 if (!changedPoint_[pointi])
186 {
187 changedPoint_[pointi] = true;
188 changedPoints_[nChangedPoints_++] = pointi;
189 }
190 }
191
192 if (!wasValid && pointInfo.valid(td_))
193 {
194 --nUnvisitedPoints_;
195 }
196
197 return propagate;
198}
199
200
201// Update info for pointi, at position pt, with information from
202// same point.
203// Updates:
204// - changedPoint_, changedPoints_, nChangedPoints_,
205// - statistics: nEvals_, nUnvisitedPoints_
206template<class Type, class TrackingData>
208(
209 const label pointi,
210 const Type& neighbourInfo,
211 Type& pointInfo
212)
213{
214 nEvals_++;
215
216 bool wasValid = pointInfo.valid(td_);
217
218 bool propagate =
219 pointInfo.updatePoint
220 (
221 mesh_,
222 pointi,
223 neighbourInfo,
224 propagationTol_,
225 td_
226 );
227
228 if (propagate)
229 {
230 if (!changedPoint_[pointi])
231 {
232 changedPoint_[pointi] = true;
233 changedPoints_[nChangedPoints_++] = pointi;
234 }
235 }
236
237 if (!wasValid && pointInfo.valid(td_))
238 {
239 --nUnvisitedPoints_;
240 }
241
242 return propagate;
243}
244
245
246// Update info for edgeI, at position pt, with information from
247// neighbouring point.
248// Updates:
249// - changedEdge_, changedEdges_, nChangedEdges_,
250// - statistics: nEvals_, nUnvisitedEdge_
251template<class Type, class TrackingData>
253(
254 const label edgeI,
255 const label neighbourPointi,
256 const Type& neighbourInfo,
257 Type& edgeInfo
258)
259{
260 nEvals_++;
261
262 bool wasValid = edgeInfo.valid(td_);
263
264 bool propagate =
265 edgeInfo.updateEdge
266 (
267 mesh_,
268 edgeI,
269 neighbourPointi,
270 neighbourInfo,
271 propagationTol_,
272 td_
273 );
274
275 if (propagate)
276 {
277 if (!changedEdge_[edgeI])
278 {
279 changedEdge_[edgeI] = true;
280 changedEdges_[nChangedEdges_++] = edgeI;
281 }
282 }
283
284 if (!wasValid && edgeInfo.valid(td_))
285 {
286 --nUnvisitedEdges_;
287 }
288
289 return propagate;
290}
291
292
293// Check if patches of given type name are present
294template<class Type, class TrackingData>
295template<class PatchType>
297{
298 label nPatches = 0;
299
300 forAll(mesh_.boundaryMesh(), patchi)
301 {
302 if (isA<PatchType>(mesh_.boundaryMesh()[patchi]))
303 {
304 nPatches++;
305 }
306 }
307 return nPatches;
308}
309
310
311// Transfer all the information to/from neighbouring processors
312template<class Type, class TrackingData>
314{
315 // 1. Send all point info on processor patches.
316
317 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
318
319 DynamicList<Type> patchInfo;
320 DynamicList<label> thisPoints;
321 DynamicList<label> nbrPoints;
322
323 forAll(mesh_.globalData().processorPatches(), i)
324 {
325 label patchi = mesh_.globalData().processorPatches()[i];
326 const processorPolyPatch& procPatch =
327 refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
328
329 patchInfo.clear();
330 patchInfo.reserve(procPatch.nPoints());
331 thisPoints.clear();
332 thisPoints.reserve(procPatch.nPoints());
333 nbrPoints.clear();
334 nbrPoints.reserve(procPatch.nPoints());
335
336 // Get all changed points in reverse order
337 const labelList& neighbPoints = procPatch.neighbPoints();
338 forAll(neighbPoints, thisPointi)
339 {
340 label meshPointi = procPatch.meshPoints()[thisPointi];
341 if (changedPoint_[meshPointi])
342 {
343 patchInfo.append(allPointInfo_[meshPointi]);
344 thisPoints.append(thisPointi);
345 nbrPoints.append(neighbPoints[thisPointi]);
346 }
347 }
348
349 // Adapt for leaving domain
350 leaveDomain(procPatch, thisPoints, patchInfo);
351
352 //if (debug)
353 //{
354 // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
355 // << " communicating with " << procPatch.neighbProcNo()
356 // << " Sending:" << patchInfo.size() << endl;
357 //}
358
359 UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
360 toNeighbour << nbrPoints << patchInfo;
361 }
362
363
364 pBufs.finishedSends();
365
366 //
367 // 2. Receive all point info on processor patches.
368 //
369
370 forAll(mesh_.globalData().processorPatches(), i)
371 {
372 label patchi = mesh_.globalData().processorPatches()[i];
373 const processorPolyPatch& procPatch =
374 refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
375
376 List<Type> patchInfo;
377 labelList patchPoints;
378
379 {
380 UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
381 fromNeighbour >> patchPoints >> patchInfo;
382 }
383
384 //if (debug)
385 //{
386 // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
387 // << " communicating with " << procPatch.neighbProcNo()
388 // << " Received:" << patchInfo.size() << endl;
389 //}
390
391 // Apply transform to received data for non-parallel planes
392 if (!procPatch.parallel())
393 {
394 transform(procPatch, procPatch.forwardT(), patchInfo);
395 }
396
397 // Adapt for entering domain
398 enterDomain(procPatch, patchPoints, patchInfo);
399
400 // Merge received info
401 const labelList& meshPoints = procPatch.meshPoints();
402 forAll(patchInfo, i)
403 {
404 label meshPointi = meshPoints[patchPoints[i]];
405
406 if (!allPointInfo_[meshPointi].equal(patchInfo[i], td_))
407 {
408 updatePoint
409 (
410 meshPointi,
411 patchInfo[i],
412 allPointInfo_[meshPointi]
413 );
414 }
415 }
416 }
417
418 // Collocated points should be handled by face based transfer
419 // (since that is how connectivity is worked out)
420 // They are also explicitly equalised in handleCollocatedPoints to
421 // guarantee identical values.
422}
423
424
425template<class Type, class TrackingData>
427{
428 // 1. Send all point info on cyclic patches.
429
430 DynamicList<Type> nbrInfo;
431 DynamicList<label> nbrPoints;
432 DynamicList<label> thisPoints;
433
434 forAll(mesh_.boundaryMesh(), patchi)
435 {
436 const polyPatch& patch = mesh_.boundaryMesh()[patchi];
437
438 const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(patch);
439
440 if (cpp)
441 {
442 const auto& cycPatch = *cpp;
443 const auto& nbrPatch = cycPatch.neighbPatch();
444
445 nbrInfo.clear();
446 nbrInfo.reserve(cycPatch.nPoints());
447 nbrPoints.clear();
448 nbrPoints.reserve(cycPatch.nPoints());
449 thisPoints.clear();
450 thisPoints.reserve(cycPatch.nPoints());
451
452 // Collect nbrPatch points that have changed
453 {
454 const edgeList& pairs = cycPatch.coupledPoints();
455 const labelList& meshPoints = nbrPatch.meshPoints();
456
457 forAll(pairs, pairI)
458 {
459 label thisPointi = pairs[pairI][0];
460 label nbrPointi = pairs[pairI][1];
461 label meshPointi = meshPoints[nbrPointi];
462
463 if (changedPoint_[meshPointi])
464 {
465 nbrInfo.append(allPointInfo_[meshPointi]);
466 nbrPoints.append(nbrPointi);
467 thisPoints.append(thisPointi);
468 }
469 }
470
471 // nbr : adapt for leaving domain
472 leaveDomain(nbrPatch, nbrPoints, nbrInfo);
473 }
474
475
476 // Apply rotation for non-parallel planes
477
478 if (!cycPatch.parallel())
479 {
480 // received data from half1
481 transform(cycPatch, cycPatch.forwardT(), nbrInfo);
482 }
483
484 //if (debug)
485 //{
486 // Pout<< "Cyclic patch " << patchi << ' ' << patch.name()
487 // << " Changed : " << nbrInfo.size()
488 // << endl;
489 //}
490
491 // Adapt for entering domain
492 enterDomain(cycPatch, thisPoints, nbrInfo);
493
494 // Merge received info
495 const labelList& meshPoints = cycPatch.meshPoints();
496 forAll(nbrInfo, i)
497 {
498 label meshPointi = meshPoints[thisPoints[i]];
499
500 if (!allPointInfo_[meshPointi].equal(nbrInfo[i], td_))
501 {
502 updatePoint
503 (
504 meshPointi,
505 nbrInfo[i],
506 allPointInfo_[meshPointi]
507 );
508 }
509 }
510 }
511 }
512}
513
514
515// Guarantee collocated points have same information.
516// Return number of points changed.
517template<class Type, class TrackingData>
519{
520 // Transfer onto coupled patch
521 const globalMeshData& gmd = mesh_.globalData();
522 const indirectPrimitivePatch& cpp = gmd.coupledPatch();
523 const labelList& meshPoints = cpp.meshPoints();
524
525 const mapDistribute& slavesMap = gmd.globalPointSlavesMap();
526 const labelListList& slaves = gmd.globalPointSlaves();
527
528 List<Type> elems(slavesMap.constructSize());
529 forAll(meshPoints, pointi)
530 {
531 elems[pointi] = allPointInfo_[meshPoints[pointi]];
532 }
533
534 // Pull slave data onto master (which might or might not have any
535 // initialised points). No need to update transformed slots.
536 slavesMap.distribute(elems, false);
537
538 // Combine master data with slave data
539 combineEqOp<Type, TrackingData> cop(td_);
540
541 forAll(slaves, pointi)
542 {
543 Type& elem = elems[pointi];
544
545 const labelList& slavePoints = slaves[pointi];
546
547 // Combine master with untransformed slave data
548 forAll(slavePoints, j)
549 {
550 cop(elem, elems[slavePoints[j]]);
551 }
552
553 // Copy result back to slave slots
554 forAll(slavePoints, j)
555 {
556 elems[slavePoints[j]] = elem;
557 }
558 }
559
560 // Push slave-slot data back to slaves
561 slavesMap.reverseDistribute(elems.size(), elems, false);
562
563 // Extract back onto mesh
564 forAll(meshPoints, pointi)
565 {
566 if (elems[pointi].valid(td_))
567 {
568 label meshPointi = meshPoints[pointi];
569
570 Type& elem = allPointInfo_[meshPointi];
571
572 bool wasValid = elem.valid(td_);
573
574 // Like updatePoint but bypass Type::updatePoint with its tolerance
575 // checking
576 //if (!elem.valid(td_) || !elem.equal(elems[pointi], td_))
577 if (!elem.equal(elems[pointi], td_))
578 {
579 nEvals_++;
580 elem = elems[pointi];
581
582 // See if element now valid
583 if (!wasValid && elem.valid(td_))
584 {
585 --nUnvisitedPoints_;
586 }
587
588 // Update database of changed points
589 if (!changedPoint_[meshPointi])
590 {
591 changedPoint_[meshPointi] = true;
592 changedPoints_[nChangedPoints_++] = meshPointi;
593 }
594 }
595 }
596 }
597
598 // Sum nChangedPoints over all procs
599 label totNChanged = nChangedPoints_;
600
601 reduce(totNChanged, sumOp<label>());
602
603 return totNChanged;
604}
605
606
607// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
608
609// Iterate, propagating changedPointsInfo across mesh, until no change (or
610// maxIter reached). Initial point values specified.
611template<class Type, class TrackingData>
613(
614 const polyMesh& mesh,
615 const labelList& changedPoints,
616 const List<Type>& changedPointsInfo,
617
618 UList<Type>& allPointInfo,
619 UList<Type>& allEdgeInfo,
620
621 const label maxIter,
622 TrackingData& td
623)
624:
625 mesh_(mesh),
626 allPointInfo_(allPointInfo),
627 allEdgeInfo_(allEdgeInfo),
628 td_(td),
629 changedPoint_(mesh_.nPoints(), false),
630 changedPoints_(mesh_.nPoints()),
631 nChangedPoints_(0),
632 changedEdge_(mesh_.nEdges(), false),
633 changedEdges_(mesh_.nEdges()),
634 nChangedEdges_(0),
635 nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
636 nEvals_(0),
637 nUnvisitedPoints_(mesh_.nPoints()),
638 nUnvisitedEdges_(mesh_.nEdges())
639{
640 if (allPointInfo_.size() != mesh_.nPoints())
641 {
643 << "size of pointInfo work array is not equal to the number"
644 << " of points in the mesh" << endl
645 << " pointInfo :" << allPointInfo_.size() << endl
646 << " mesh.nPoints:" << mesh_.nPoints()
647 << exit(FatalError);
648 }
649 if (allEdgeInfo_.size() != mesh_.nEdges())
650 {
652 << "size of edgeInfo work array is not equal to the number"
653 << " of edges in the mesh" << endl
654 << " edgeInfo :" << allEdgeInfo_.size() << endl
655 << " mesh.nEdges:" << mesh_.nEdges()
656 << exit(FatalError);
657 }
658
659
660 // Set from initial changed points data
661 setPointInfo(changedPoints, changedPointsInfo);
662
663 if (debug)
664 {
665 Info<< typeName << ": Seed points : "
666 << returnReduce(nChangedPoints_, sumOp<label>()) << endl;
667 }
668
669 // Iterate until nothing changes
670 label iter = iterate(maxIter);
671
672 if ((maxIter > 0) && (iter >= maxIter))
673 {
675 << "Maximum number of iterations reached. Increase maxIter." << endl
676 << " maxIter:" << maxIter << endl
677 << " nChangedPoints:" << nChangedPoints_ << endl
678 << " nChangedEdges:" << nChangedEdges_ << endl
679 << exit(FatalError);
680 }
681}
682
683
684template<class Type, class TrackingData>
686(
687 const polyMesh& mesh,
688 UList<Type>& allPointInfo,
689 UList<Type>& allEdgeInfo,
690 TrackingData& td
691)
692:
693 mesh_(mesh),
694 allPointInfo_(allPointInfo),
695 allEdgeInfo_(allEdgeInfo),
696 td_(td),
697 changedPoint_(mesh_.nPoints(), false),
698 changedPoints_(mesh_.nPoints()),
699 nChangedPoints_(0),
700 changedEdge_(mesh_.nEdges(), false),
701 changedEdges_(mesh_.nEdges()),
702 nChangedEdges_(0),
703 nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
704 nEvals_(0),
705 nUnvisitedPoints_(mesh_.nPoints()),
706 nUnvisitedEdges_(mesh_.nEdges())
707{}
708
709
710// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
711
712template<class Type, class TrackingData>
714{
715 return nUnvisitedPoints_;
716}
717
718
719template<class Type, class TrackingData>
721{
722 return nUnvisitedEdges_;
723}
724
725
726// Copy point information into member data
727template<class Type, class TrackingData>
729(
730 const labelList& changedPoints,
731 const List<Type>& changedPointsInfo
732)
733{
734 forAll(changedPoints, changedPointi)
735 {
736 const label pointi = changedPoints[changedPointi];
737
738 const bool wasValid = allPointInfo_[pointi].valid(td_);
739
740 // Copy info for pointi
741 allPointInfo_[pointi] = changedPointsInfo[changedPointi];
742
743 // Maintain count of unset points
744 if (!wasValid && allPointInfo_[pointi].valid(td_))
745 {
746 --nUnvisitedPoints_;
747 }
748
749 // Mark pointi as changed, both on list and on point itself.
750
751 if (!changedPoint_[pointi])
752 {
753 changedPoint_[pointi] = true;
754 changedPoints_[nChangedPoints_++] = pointi;
755 }
756 }
757
758 // Sync
759 handleCollocatedPoints();
760}
761
762
763// Propagate information from edge to point. Return number of points changed.
764template<class Type, class TrackingData>
766{
767 for
768 (
769 label changedEdgeI = 0;
770 changedEdgeI < nChangedEdges_;
771 changedEdgeI++
772 )
773 {
774 label edgeI = changedEdges_[changedEdgeI];
775
776 if (!changedEdge_[edgeI])
777 {
779 << "edge " << edgeI
780 << " not marked as having been changed" << nl
781 << "This might be caused by multiple occurrences of the same"
782 << " seed point." << abort(FatalError);
783 }
784
785
786 const Type& neighbourWallInfo = allEdgeInfo_[edgeI];
787
788 // Evaluate all connected points (= edge endpoints)
789 const edge& e = mesh_.edges()[edgeI];
790
791 forAll(e, eI)
792 {
793 Type& currentWallInfo = allPointInfo_[e[eI]];
794
795 if (!currentWallInfo.equal(neighbourWallInfo, td_))
796 {
797 updatePoint
798 (
799 e[eI],
800 edgeI,
801 neighbourWallInfo,
802 currentWallInfo
803 );
804 }
805 }
806
807 // Reset status of edge
808 changedEdge_[edgeI] = false;
809 }
810
811 // Handled all changed edges by now
812 nChangedEdges_ = 0;
813
814 if (nCyclicPatches_ > 0)
815 {
816 // Transfer changed points across cyclic halves
817 handleCyclicPatches();
818 }
819 if (Pstream::parRun())
820 {
821 // Transfer changed points from neighbouring processors.
822 handleProcPatches();
823 }
824
825 //if (debug)
826 //{
827 // Pout<< "Changed points : " << nChangedPoints_ << endl;
828 //}
829
830 // Sum nChangedPoints over all procs
831 label totNChanged = nChangedPoints_;
832
833 reduce(totNChanged, sumOp<label>());
834
835 return totNChanged;
836}
837
838
839// Propagate information from point to edge. Return number of edges changed.
840template<class Type, class TrackingData>
842{
843 const labelListList& pointEdges = mesh_.pointEdges();
844
845 for
846 (
847 label changedPointi = 0;
848 changedPointi < nChangedPoints_;
849 changedPointi++
850 )
851 {
852 label pointi = changedPoints_[changedPointi];
853
854 if (!changedPoint_[pointi])
855 {
857 << "Point " << pointi
858 << " not marked as having been changed" << nl
859 << "This might be caused by multiple occurrences of the same"
860 << " seed point." << abort(FatalError);
861 }
862
863 const Type& neighbourWallInfo = allPointInfo_[pointi];
864
865 // Evaluate all connected edges
866
867 const labelList& edgeLabels = pointEdges[pointi];
868 forAll(edgeLabels, edgeLabelI)
869 {
870 label edgeI = edgeLabels[edgeLabelI];
871
872 Type& currentWallInfo = allEdgeInfo_[edgeI];
873
874 if (!currentWallInfo.equal(neighbourWallInfo, td_))
875 {
876 updateEdge
877 (
878 edgeI,
879 pointi,
880 neighbourWallInfo,
881 currentWallInfo
882 );
883 }
884 }
885
886 // Reset status of point
887 changedPoint_[pointi] = false;
888 }
889
890 // Handled all changed points by now
891 nChangedPoints_ = 0;
892
893 //if (debug)
894 //{
895 // Pout<< "Changed edges : " << nChangedEdges_ << endl;
896 //}
897
898 // Sum nChangedPoints over all procs
899 label totNChanged = nChangedEdges_;
900
901 reduce(totNChanged, sumOp<label>());
902
903 return totNChanged;
904}
905
906
907// Iterate
908template<class Type, class TrackingData>
910(
911 const label maxIter
912)
913{
914 if (nCyclicPatches_ > 0)
915 {
916 // Transfer changed points across cyclic halves
917 handleCyclicPatches();
918 }
919 if (Pstream::parRun())
920 {
921 // Transfer changed points from neighbouring processors.
922 handleProcPatches();
923 }
924
925 nEvals_ = 0;
926
927 label iter = 0;
928
929 while (iter < maxIter)
930 {
931 while (iter < maxIter)
932 {
933 if (debug)
934 {
935 Info<< typeName << ": Iteration " << iter << endl;
936 }
937
938 label nEdges = pointToEdge();
939
940 if (debug)
941 {
942 Info<< typeName << ": Total changed edges : "
943 << nEdges << endl;
944 }
945
946 if (nEdges == 0)
947 {
948 break;
949 }
950
951 label nPoints = edgeToPoint();
952
953 if (debug)
954 {
955 Info<< typeName << ": Total changed points : "
956 << nPoints << nl
957 << typeName << ": Total evaluations : "
958 << returnReduce(nEvals_, sumOp<label>()) << nl
959 << typeName << ": Remaining unvisited points: "
960 << returnReduce(nUnvisitedPoints_, sumOp<label>()) << nl
961 << typeName << ": Remaining unvisited edges : "
962 << returnReduce(nUnvisitedEdges_, sumOp<label>()) << nl
963 << endl;
964 }
965
966 if (nPoints == 0)
967 {
968 break;
969 }
970
971 iter++;
972 }
973
974
975 // Enforce collocated points are exactly equal. This might still mean
976 // non-collocated points are not equal though. WIP.
977 label nPoints = handleCollocatedPoints();
978 if (debug)
979 {
980 Info<< typeName << ": Collocated point sync : "
981 << nPoints << nl << endl;
982 }
983
984 if (nPoints == 0)
985 {
986 break;
987 }
988 }
989
990 return iter;
991}
992
993
994// ************************************************************************* //
scalar y
reduce(hasMovingMesh, orOp< bool >())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
void reserve(const label len)
Definition: DynamicListI.H:333
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:91
label nUnvisitedEdges() const
Number of unvisited edges, i.e. edges that were not (yet)
void setPointInfo(const labelList &changedPoints, const List< Type > &changedPointsInfo)
Copy initial data into allPointInfo_.
label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
label pointToEdge()
Propagate from point to edge. Returns total number of edges.
label edgeToPoint()
Propagate from edge to point. Returns total number of points.
label nUnvisitedPoints() const
label nPoints() const
Number of points supporting patch faces.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
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
Reduction class. If x and y are not equal assign value.
Definition: PointEdgeWave.C:53
combineEqOp(TrackingData &td)
Definition: PointEdgeWave.C:58
void operator()(Type &x, const Type &y) const
Definition: PointEdgeWave.C:63
Cyclic plane patch.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Default transformation behaviour.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
label nPoints() const noexcept
Number of mesh points.
label nEdges() const
Number of mesh edges.
Neighbour processor patch.
int neighbProcNo() const
Return neighbour processor number.
const labelList & neighbPoints() const
Return neighbour point labels. WIP.
Tensor of scalars, i.e. Tensor<scalar>.
const volScalarField & T
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const label nPatches
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
label nPoints
bool equal(const Matrix< Form1, Type > &A, const Matrix< Form2, Type > &B, const bool verbose=false, const label maxDiffs=10, const scalar relTol=1e-5, const scalar absTol=1e-8)
Compare matrix elements for absolute or relative equality.
Definition: MatrixTools.C:34
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
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 reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
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