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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "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 
42 namespace Foam
43 {
44  //- Reduction class. If x and y are not equal assign value.
45  template<class Type, class TrackingData>
46  class combineEqOp
47  {
48  TrackingData& td_;
49 
50  public:
51 
52  combineEqOp(TrackingData& td)
53  :
54  td_(td)
55  {}
56 
57  void operator()(Type& x, const Type& y) const
58  {
59  if (!x.valid(td_) && y.valid(td_))
60  {
61  x = y;
62  }
63  }
64  };
65 }
66 
67 
68 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69 
70 // Handle leaving domain. Implementation referred to Type
71 template<class Type, class TrackingData>
73 (
74  const polyPatch& patch,
75  const labelList& patchPointLabels,
76  List<Type>& pointInfo
77 ) const
78 {
79  const labelList& meshPoints = patch.meshPoints();
80 
81  forAll(patchPointLabels, i)
82  {
83  label patchPointi = patchPointLabels[i];
84 
85  const point& pt = patch.points()[meshPoints[patchPointi]];
86 
87  pointInfo[i].leaveDomain(patch, patchPointi, pt, td_);
88  }
89 }
90 
91 
92 // Handle entering domain. Implementation referred to Type
93 template<class Type, class TrackingData>
95 (
96  const polyPatch& patch,
97  const labelList& patchPointLabels,
98  List<Type>& pointInfo
99 ) const
100 {
101  const labelList& meshPoints = patch.meshPoints();
102 
103  forAll(patchPointLabels, i)
104  {
105  label patchPointi = patchPointLabels[i];
106 
107  const point& pt = patch.points()[meshPoints[patchPointi]];
108 
109  pointInfo[i].enterDomain(patch, patchPointi, pt, td_);
110  }
111 }
112 
113 
114 // Transform. Implementation referred to Type
115 template<class Type, class TrackingData>
117 (
118  const polyPatch& patch,
119  const tensorField& rotTensor,
120  List<Type>& pointInfo
121 ) const
122 {
123  if (rotTensor.size() == 1)
124  {
125  const tensor& T = rotTensor[0];
126 
127  forAll(pointInfo, i)
128  {
129  pointInfo[i].transform(T, td_);
130  }
131  }
132  else
133  {
135  << "Non-uniform transformation on patch " << patch.name()
136  << " of type " << patch.type()
137  << " not supported for point fields"
138  << abort(FatalError);
139 
140  forAll(pointInfo, i)
141  {
142  pointInfo[i].transform(rotTensor[i], td_);
143  }
144  }
145 }
146 
147 
148 // Update info for pointi, at position pt, with information from
149 // neighbouring edge.
150 // Updates:
151 // - changedPoint_, changedPoints_,
152 // - statistics: nEvals_, nUnvisitedPoints_
153 template<class Type, class TrackingData>
155 (
156  const label pointi,
157  const label neighbourEdgeI,
158  const Type& neighbourInfo,
159  Type& pointInfo
160 )
161 {
162  nEvals_++;
163 
164  bool wasValid = pointInfo.valid(td_);
165 
166  bool propagate =
167  pointInfo.updatePoint
168  (
169  mesh_,
170  pointi,
171  neighbourEdgeI,
172  neighbourInfo,
173  propagationTol_,
174  td_
175  );
176 
177  if (propagate)
178  {
179  if (changedPoint_.set(pointi))
180  {
181  changedPoints_.push_back(pointi);
182  }
183  }
184 
185  if (!wasValid && pointInfo.valid(td_))
186  {
187  --nUnvisitedPoints_;
188  }
189 
190  return propagate;
191 }
192 
193 
194 // Update info for pointi, at position pt, with information from
195 // same point.
196 // Updates:
197 // - changedPoint_, changedPoints_,
198 // - statistics: nEvals_, nUnvisitedPoints_
199 template<class Type, class TrackingData>
201 (
202  const label pointi,
203  const Type& neighbourInfo,
204  Type& pointInfo
205 )
206 {
207  nEvals_++;
208 
209  bool wasValid = pointInfo.valid(td_);
210 
211  bool propagate =
212  pointInfo.updatePoint
213  (
214  mesh_,
215  pointi,
216  neighbourInfo,
217  propagationTol_,
218  td_
219  );
220 
221  if (propagate)
222  {
223  if (changedPoint_.set(pointi))
224  {
225  changedPoints_.push_back(pointi);
226  }
227  }
228 
229  if (!wasValid && pointInfo.valid(td_))
230  {
231  --nUnvisitedPoints_;
232  }
233 
234  return propagate;
235 }
236 
237 
238 // Update info for edgei, at position pt, with information from
239 // neighbouring point.
240 // Updates:
241 // - changedEdge_, changedEdges_,
242 // - statistics: nEvals_, nUnvisitedEdge_
243 template<class Type, class TrackingData>
245 (
246  const label edgei,
247  const label neighbourPointi,
248  const Type& neighbourInfo,
249  Type& edgeInfo
250 )
251 {
252  nEvals_++;
253 
254  bool wasValid = edgeInfo.valid(td_);
255 
256  bool propagate =
257  edgeInfo.updateEdge
258  (
259  mesh_,
260  edgei,
261  neighbourPointi,
262  neighbourInfo,
263  propagationTol_,
264  td_
265  );
266 
267  if (propagate)
268  {
269  if (changedEdge_.set(edgei))
270  {
271  changedEdges_.push_back(edgei);
272  }
273  }
274 
275  if (!wasValid && edgeInfo.valid(td_))
276  {
277  --nUnvisitedEdges_;
278  }
279 
280  return propagate;
281 }
283 
284 // Check if patches of given type name are present
285 template<class Type, class TrackingData>
286 template<class PatchType>
288 {
289  label nPatches = 0;
290 
291  for (const polyPatch& p : mesh_.boundaryMesh())
292  {
293  if (isA<PatchType>(p))
294  {
295  ++nPatches;
296  }
297  }
298  return nPatches;
299 }
300 
301 
302 // Transfer all the information to/from neighbouring processors
303 template<class Type, class TrackingData>
305 {
306  // 1. Send all point info on processor patches.
307 
308  const globalMeshData& pData = mesh_.globalData();
309 
310  // Which patches are processor patches
311  const labelList& procPatches = pData.processorPatches();
312 
313  // Which processors this processor is connected to
314  const labelList& neighbourProcs = pData.topology().procNeighbours();
315 
316  // Reduce communication by only sending non-zero data,
317  // but with multiply-connected processor/processor
318  // (eg, processorCyclic) also need to send zero information
319  // to keep things synchronised
320 
321  // Reset buffers, initialize for registerSend() bookkeeping
322  pBufs_.clear();
323  pBufs_.initRegisterSend();
324 
325 
326  // Information to send
327  DynamicList<Type> patchInfo;
328  DynamicList<label> thisPoints;
329  DynamicList<label> nbrPoints;
330 
331  for (const label patchi : procPatches)
332  {
333  const auto& procPatch =
334  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
335 
336  const label nbrProci = procPatch.neighbProcNo();
337 
338  patchInfo.clear();
339  patchInfo.reserve(procPatch.nPoints());
340 
341  thisPoints.clear();
342  thisPoints.reserve(procPatch.nPoints());
343 
344  nbrPoints.clear();
345  nbrPoints.reserve(procPatch.nPoints());
346 
347  // Get all changed points in reverse order
348  const labelList& neighbPoints = procPatch.neighbPoints();
349  forAll(neighbPoints, thisPointi)
350  {
351  label meshPointi = procPatch.meshPoints()[thisPointi];
352  if (changedPoint_.test(meshPointi))
353  {
354  patchInfo.append(allPointInfo_[meshPointi]);
355  thisPoints.append(thisPointi);
356  nbrPoints.append(neighbPoints[thisPointi]);
357  }
358  }
359 
360  // Adapt for leaving domain
361  leaveDomain(procPatch, thisPoints, patchInfo);
362 
363  // Send to neighbour
364  {
365  UOPstream toNbr(nbrProci, pBufs_);
366  toNbr << nbrPoints << patchInfo;
367 
368  // Record if send is required (data are non-zero)
369  pBufs_.registerSend(nbrProci, !patchInfo.empty());
370 
371  //if (debug & 2)
372  //{
373  // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
374  // << " send:" << patchInfo.size()
375  // << " to proc:" << nbrProci << endl;
376  //}
377  }
378  }
379 
380  // Limit exchange to involved procs
381  // - automatically discards unnecessary (unregistered) sends
382  pBufs_.finishedNeighbourSends(neighbourProcs);
383 
384 
385  //
386  // 2. Receive all point info on processor patches.
387  //
388 
389  for (const label patchi : procPatches)
390  {
391  const auto& procPatch =
392  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
393 
394  const label nbrProci = procPatch.neighbProcNo();
395 
396  if (!pBufs_.recvDataCount(nbrProci))
397  {
398  // Nothing to receive
399  continue;
400  }
401 
402  labelList patchPoints;
403  List<Type> patchInfo;
404  {
405  UIPstream is(nbrProci, pBufs_);
406  is >> patchPoints >> patchInfo;
407  }
408 
409  //if (debug & 2)
410  //{
411  // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
412  // << " recv:" << patchInfo.size() << " from proc:"
413  // << nbrProci << endl;
414  //}
415 
416  // Apply transform to received data for non-parallel planes
417  if (!procPatch.parallel())
418  {
419  transform(procPatch, procPatch.forwardT(), patchInfo);
420  }
421 
422  // Adapt for entering domain
423  enterDomain(procPatch, patchPoints, patchInfo);
424 
425  // Merge received info
426  const labelList& meshPoints = procPatch.meshPoints();
427  forAll(patchInfo, i)
428  {
429  label meshPointi = meshPoints[patchPoints[i]];
430 
431  if (!allPointInfo_[meshPointi].equal(patchInfo[i], td_))
432  {
433  updatePoint
434  (
435  meshPointi,
436  patchInfo[i],
437  allPointInfo_[meshPointi]
438  );
439  }
440  }
441  }
442 
443  // Collocated points should be handled by face based transfer
444  // (since that is how connectivity is worked out)
445  // They are also explicitly equalised in handleCollocatedPoints to
446  // guarantee identical values.
447 }
448 
449 
450 template<class Type, class TrackingData>
452 {
453  // 1. Send all point info on cyclic patches.
454 
455  DynamicList<Type> nbrInfo;
456  DynamicList<label> nbrPoints;
457  DynamicList<label> thisPoints;
458 
459  for (const polyPatch& patch : mesh_.boundaryMesh())
460  {
461  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(patch);
462 
463  if (cpp)
464  {
465  const auto& cycPatch = *cpp;
466  const auto& nbrPatch = cycPatch.neighbPatch();
467 
468  nbrInfo.clear();
469  nbrInfo.reserve(cycPatch.nPoints());
470  nbrPoints.clear();
471  nbrPoints.reserve(cycPatch.nPoints());
472  thisPoints.clear();
473  thisPoints.reserve(cycPatch.nPoints());
474 
475  // Collect nbrPatch points that have changed
476  {
477  const edgeList& pairs = cycPatch.coupledPoints();
478  const labelList& meshPoints = nbrPatch.meshPoints();
479 
480  forAll(pairs, pairI)
481  {
482  label thisPointi = pairs[pairI][0];
483  label nbrPointi = pairs[pairI][1];
484  label meshPointi = meshPoints[nbrPointi];
485 
486  if (changedPoint_.test(meshPointi))
487  {
488  nbrInfo.append(allPointInfo_[meshPointi]);
489  nbrPoints.append(nbrPointi);
490  thisPoints.append(thisPointi);
491  }
492  }
493 
494  // nbr : adapt for leaving domain
495  leaveDomain(nbrPatch, nbrPoints, nbrInfo);
496  }
497 
498 
499  // Apply rotation for non-parallel planes
500 
501  if (!cycPatch.parallel())
502  {
503  // received data from half1
504  transform(cycPatch, cycPatch.forwardT(), nbrInfo);
505  }
506 
507  //if (debug)
508  //{
509  // Pout<< "Cyclic patch " << patch.index() << ' ' << patch.name()
510  // << " Changed : " << nbrInfo.size()
511  // << endl;
512  //}
513 
514  // Adapt for entering domain
515  enterDomain(cycPatch, thisPoints, nbrInfo);
516 
517  // Merge received info
518  const labelList& meshPoints = cycPatch.meshPoints();
519  forAll(nbrInfo, i)
520  {
521  label meshPointi = meshPoints[thisPoints[i]];
522 
523  if (!allPointInfo_[meshPointi].equal(nbrInfo[i], td_))
524  {
525  updatePoint
526  (
527  meshPointi,
528  nbrInfo[i],
529  allPointInfo_[meshPointi]
530  );
531  }
532  }
533  }
534  }
535 }
536 
537 
538 // Guarantee collocated points have same information.
539 // Return number of points changed.
540 template<class Type, class TrackingData>
542 {
543  // Transfer onto coupled patch
544  const globalMeshData& gmd = mesh_.globalData();
545  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
546  const labelList& meshPoints = cpp.meshPoints();
547 
548  const mapDistribute& slavesMap = gmd.globalPointSlavesMap();
549  const labelListList& slaves = gmd.globalPointSlaves();
550 
551  List<Type> elems(slavesMap.constructSize());
552  forAll(meshPoints, pointi)
553  {
554  elems[pointi] = allPointInfo_[meshPoints[pointi]];
555  }
556 
557  // Pull slave data onto master (which might or might not have any
558  // initialised points). No need to update transformed slots.
559  slavesMap.distribute(elems, false);
560 
561  // Combine master data with slave data
562  combineEqOp<Type, TrackingData> cop(td_);
563 
564  forAll(slaves, pointi)
565  {
566  Type& elem = elems[pointi];
567 
568  const labelList& slavePoints = slaves[pointi];
569 
570  // Combine master with untransformed slave data
571  forAll(slavePoints, j)
572  {
573  cop(elem, elems[slavePoints[j]]);
574  }
575 
576  // Copy result back to slave slots
577  forAll(slavePoints, j)
578  {
579  elems[slavePoints[j]] = elem;
580  }
581  }
582 
583  // Push slave-slot data back to slaves
584  slavesMap.reverseDistribute(elems.size(), elems, false);
585 
586  // Extract back onto mesh
587  forAll(meshPoints, pointi)
588  {
589  if (elems[pointi].valid(td_))
590  {
591  label meshPointi = meshPoints[pointi];
592 
593  Type& elem = allPointInfo_[meshPointi];
594 
595  bool wasValid = elem.valid(td_);
596 
597  // Like updatePoint but bypass Type::updatePoint with its tolerance
598  // checking
599  //if (!elem.valid(td_) || !elem.equal(elems[pointi], td_))
600  if (!elem.equal(elems[pointi], td_))
601  {
602  nEvals_++;
603  elem = elems[pointi];
604 
605  // See if element now valid
606  if (!wasValid && elem.valid(td_))
607  {
608  --nUnvisitedPoints_;
609  }
610 
611  // Update database of changed points
612  if (changedPoint_.set(meshPointi))
613  {
614  changedPoints_.push_back(meshPointi);
615  }
616  }
617  }
618  }
619 
620  // Sum changedPoints over all procs
621  return returnReduce(nChangedPoints(), sumOp<label>());
622 }
623 
624 
625 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
627 // Iterate, propagating changedPointsInfo across mesh, until no change (or
628 // maxIter reached). Initial point values specified.
629 template<class Type, class TrackingData>
631 (
632  const polyMesh& mesh,
633  const labelList& changedPoints,
634  const List<Type>& changedPointsInfo,
635 
636  UList<Type>& allPointInfo,
637  UList<Type>& allEdgeInfo,
638 
639  const label maxIter,
640  TrackingData& td
641 )
642 :
644 
645  allPointInfo_(allPointInfo),
646  allEdgeInfo_(allEdgeInfo),
647  td_(td),
648  nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
649  nEvals_(0)
650 {
651  if (allPointInfo_.size() != mesh_.nPoints())
652  {
654  << "size of pointInfo work array is not equal to the number"
655  << " of points in the mesh" << endl
656  << " pointInfo :" << allPointInfo_.size() << endl
657  << " mesh.nPoints:" << mesh_.nPoints()
658  << exit(FatalError);
659  }
660  if (allEdgeInfo_.size() != mesh_.nEdges())
661  {
663  << "size of edgeInfo work array is not equal to the number"
664  << " of edges in the mesh" << endl
665  << " edgeInfo :" << allEdgeInfo_.size() << endl
666  << " mesh.nEdges:" << mesh_.nEdges()
667  << exit(FatalError);
668  }
669 
670 
671  // Set from initial changed points data
672  setPointInfo(changedPoints, changedPointsInfo);
673 
674  if (debug)
675  {
676  Info<< typeName << ": Seed points : "
678  }
679 
680  // Iterate until nothing changes
681  label iter = iterate(maxIter);
682 
683  if ((maxIter > 0) && (iter >= maxIter))
684  {
686  << "Maximum number of iterations reached. Increase maxIter." << endl
687  << " maxIter:" << maxIter << nl
688  << " nChangedPoints:" << nChangedPoints() << nl
689  << " nChangedEdges:" << nChangedEdges() << endl
690  << exit(FatalError);
691  }
692 }
693 
694 
695 template<class Type, class TrackingData>
697 (
698  const polyMesh& mesh,
699  UList<Type>& allPointInfo,
700  UList<Type>& allEdgeInfo,
701  TrackingData& td
702 )
703 :
705 
706  allPointInfo_(allPointInfo),
707  allEdgeInfo_(allEdgeInfo),
708  td_(td),
709  nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
710  nEvals_(0)
711 {}
712 
713 
714 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
715 
716 // Copy point information into member data
717 template<class Type, class TrackingData>
719 (
720  const labelList& changedPoints,
721  const List<Type>& changedPointsInfo
722 )
723 {
724  forAll(changedPoints, changedPointi)
725  {
726  const label pointi = changedPoints[changedPointi];
727 
728  const bool wasValid = allPointInfo_[pointi].valid(td_);
729 
730  // Copy info for pointi
731  allPointInfo_[pointi] = changedPointsInfo[changedPointi];
732 
733  // Maintain count of unset points
734  if (!wasValid && allPointInfo_[pointi].valid(td_))
735  {
736  --nUnvisitedPoints_;
737  }
738 
739  // Mark pointi as changed, both on list and on point itself.
740 
741  if (changedPoint_.set(pointi))
742  {
743  changedPoints_.push_back(pointi);
744  }
745  }
746 
747  // Sync
748  handleCollocatedPoints();
749 }
750 
751 
752 // Propagate information from edge to point. Return number of points changed.
753 template<class Type, class TrackingData>
755 {
756  for (const label edgei : changedEdges_)
757  {
758  if (!changedEdge_.test(edgei))
759  {
761  << "edge " << edgei
762  << " not marked as having been changed" << nl
763  << "This might be caused by multiple occurrences of the same"
764  << " seed point." << abort(FatalError);
765  }
766 
767 
768  const Type& neighbourWallInfo = allEdgeInfo_[edgei];
769 
770  // Evaluate all connected points (= edge endpoints)
771  const edge& e = mesh_.edges()[edgei];
772 
773  forAll(e, eI)
774  {
775  Type& currentWallInfo = allPointInfo_[e[eI]];
776 
777  if (!currentWallInfo.equal(neighbourWallInfo, td_))
778  {
779  updatePoint
780  (
781  e[eI],
782  edgei,
783  neighbourWallInfo,
784  currentWallInfo
785  );
786  }
787  }
788 
789  // Reset status of edge
790  changedEdge_.unset(edgei);
791  }
792 
793  // Handled all changed edges by now
794  changedEdges_.clear();
795 
796  if (nCyclicPatches_ > 0)
797  {
798  // Transfer changed points across cyclic halves
799  handleCyclicPatches();
800  }
801  if (Pstream::parRun())
802  {
803  // Transfer changed points from neighbouring processors.
804  handleProcPatches();
805  }
806 
807  //if (debug)
808  //{
809  // Pout<< "Changed points : " << nChangedPoints() << nl;
810  //}
811 
812  // Sum changedPoints over all procs
813  return returnReduce(nChangedPoints(), sumOp<label>());
814 }
815 
816 
817 // Propagate information from point to edge. Return number of edges changed.
818 template<class Type, class TrackingData>
820 {
821  const labelListList& pointEdges = mesh_.pointEdges();
822 
823  for (const label pointi : changedPoints_)
824  {
825  if (!changedPoint_.test(pointi))
826  {
828  << "Point " << pointi
829  << " not marked as having been changed" << nl
830  << "This might be caused by multiple occurrences of the same"
831  << " seed point." << abort(FatalError);
832  }
833 
834  const Type& neighbourWallInfo = allPointInfo_[pointi];
835 
836  // Evaluate all connected edges
837 
838  for (const label edgei : pointEdges[pointi])
839  {
840  Type& currentWallInfo = allEdgeInfo_[edgei];
841 
842  if (!currentWallInfo.equal(neighbourWallInfo, td_))
843  {
844  updateEdge
845  (
846  edgei,
847  pointi,
848  neighbourWallInfo,
849  currentWallInfo
850  );
851  }
852  }
853 
854  // Reset status of point
855  changedPoint_.unset(pointi);
856  }
857 
858  // Handled all changed points by now
859  changedPoints_.clear();
860 
861  //if (debug)
862  //{
863  // Pout<< "Changed edges : " << nChangedEdges() << endl;
864  //}
865 
866  // Sum changedEdges over all procs
867  return returnReduce(nChangedEdges(), sumOp<label>());
868 }
870 
871 // Iterate
872 template<class Type, class TrackingData>
874 (
875  const label maxIter
876 )
877 {
878  if (nCyclicPatches_ > 0)
879  {
880  // Transfer changed points across cyclic halves
881  handleCyclicPatches();
882  }
883  if (Pstream::parRun())
884  {
885  // Transfer changed points from neighbouring processors.
886  handleProcPatches();
887  }
888 
889  nEvals_ = 0;
890 
891  label iter = 0;
892 
893  while (iter < maxIter)
894  {
895  while (iter < maxIter)
896  {
897  if (debug)
898  {
899  Info<< typeName << ": Iteration " << iter << endl;
900  }
901 
902  label nEdges = pointToEdge();
903 
904  if (debug)
905  {
906  Info<< typeName << ": Total changed edges : "
907  << nEdges << endl;
908  }
909 
910  if (nEdges == 0)
911  {
912  break;
913  }
914 
915  label nPoints = edgeToPoint();
916 
917  if (debug)
918  {
919  Info<< typeName << ": Total changed points : "
920  << nPoints << nl
921  << typeName << ": Total evaluations : "
922  << returnReduce(nEvals_, sumOp<label>()) << nl
923  << typeName << ": Remaining unvisited points: "
924  << returnReduce(nUnvisitedPoints_, sumOp<label>()) << nl
925  << typeName << ": Remaining unvisited edges : "
926  << returnReduce(nUnvisitedEdges_, sumOp<label>()) << nl
927  << endl;
928  }
929 
930  if (nPoints == 0)
931  {
932  break;
933  }
934 
935  iter++;
936  }
937 
938 
939  // Enforce collocated points are exactly equal. This might still mean
940  // non-collocated points are not equal though. WIP.
941  label nPoints = handleCollocatedPoints();
942  if (debug)
943  {
944  Info<< typeName << ": Collocated point sync : "
945  << nPoints << nl << endl;
946  }
947 
948  if (nPoints == 0)
949  {
950  break;
951  }
952  }
953 
954  return iter;
955 }
956 
957 
958 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Reduction class. If x and y are not equal assign value.
Definition: PointEdgeWave.C:41
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
label pointToEdge()
Propagate from point to edge. Returns total number of edges.
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const polyMesh & mesh_
Reference to mesh.
Definition: PointEdgeWave.H:97
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
label edgeToPoint()
Propagate from edge to point. Returns total number of points.
void setPointInfo(const labelList &changedPoints, const List< Type > &changedPointsInfo)
Copy initial data into allPointInfo_.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
combineEqOp(TrackingData &td)
Definition: PointEdgeWave.C:47
void push_back(const T &val)
Append an element at the end of the list.
Definition: ListI.H:227
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
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:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
label nPoints
Cyclic plane patch.
label nChangedEdges() const noexcept
Current number of changed edges.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
errorManip< error > abort(error &err)
Definition: errorManip.H:139
label nEdges() const
Number of mesh edges.
void reserve(const label len)
Reserve allocation space for at least this size, allocating new space if required and retaining old c...
Definition: DynamicListI.H:333
int debug
Static debugging option.
void operator()(Type &x, const Type &y) const
Definition: PointEdgeWave.C:52
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
const volScalarField & T
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
const labelList & procNeighbours() const
The neighbour processor connections (ascending order) associated with the local rank.
vector point
Point is a vector.
Definition: point.H:37
const processorTopology & topology() const noexcept
The processor to processor topology.
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Wave propagation of information through grid. Every iteration information goes through one layer of e...
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
label nChangedPoints() const noexcept
Current number of changed points.
Namespace for OpenFOAM.