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  // Reset buffers
317  pBufs_.clear();
318 
319  DynamicList<Type> patchInfo;
320  DynamicList<label> thisPoints;
321  DynamicList<label> nbrPoints;
322 
323 
324  // Reduce communication by only sending non-zero data,
325  // but with multiply-connected processor/processor
326  // (eg, processorCyclic) also need to send zero information
327  // to keep things synchronised
328 
329  // If data needs to be sent (index corresponding to neighbourProcs)
330  List<bool> isActiveSend(neighbourProcs.size(), false);
331 
332  for (const label patchi : procPatches)
333  {
334  const auto& procPatch =
335  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
336 
337  const label nbrProci = procPatch.neighbProcNo();
338 
339  patchInfo.clear();
340  patchInfo.reserve(procPatch.nPoints());
341 
342  thisPoints.clear();
343  thisPoints.reserve(procPatch.nPoints());
344 
345  nbrPoints.clear();
346  nbrPoints.reserve(procPatch.nPoints());
347 
348  // Get all changed points in reverse order
349  const labelList& neighbPoints = procPatch.neighbPoints();
350  forAll(neighbPoints, thisPointi)
351  {
352  label meshPointi = procPatch.meshPoints()[thisPointi];
353  if (changedPoint_.test(meshPointi))
354  {
355  patchInfo.append(allPointInfo_[meshPointi]);
356  thisPoints.append(thisPointi);
357  nbrPoints.append(neighbPoints[thisPointi]);
358  }
359  }
360 
361  // Adapt for leaving domain
362  leaveDomain(procPatch, thisPoints, patchInfo);
363 
364  // Record if send is required (non-empty data)
365  if (!patchInfo.empty())
366  {
367  const label nbrIndex = neighbourProcs.find(nbrProci);
368  if (nbrIndex >= 0) // Safety check (should be unnecessary)
369  {
370  isActiveSend[nbrIndex] = true;
371  }
372  }
373 
374  // Send to neighbour
375  {
376  UOPstream toNbr(nbrProci, pBufs_);
377  toNbr << nbrPoints << patchInfo;
378 
379  //if (debug & 2)
380  //{
381  // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
382  // << " send:" << patchInfo.size()
383  // << " to proc:" << nbrProci << endl;
384  //}
385  }
386  }
387 
388  // Eliminate unnecessary sends
389  forAll(neighbourProcs, nbrIndex)
390  {
391  if (!isActiveSend[nbrIndex])
392  {
393  pBufs_.clearSend(neighbourProcs[nbrIndex]);
394  }
395  }
396 
397  // Limit exchange to involved procs
398  pBufs_.finishedNeighbourSends(neighbourProcs);
399 
400 
401  //
402  // 2. Receive all point info on processor patches.
403  //
404 
405  for (const label patchi : procPatches)
406  {
407  const auto& procPatch =
408  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
409 
410  const label nbrProci = procPatch.neighbProcNo();
411 
412  if (!pBufs_.recvDataCount(nbrProci))
413  {
414  // Nothing to receive
415  continue;
416  }
417 
418  labelList patchPoints;
419  List<Type> patchInfo;
420  {
421  UIPstream is(nbrProci, pBufs_);
422  is >> patchPoints >> patchInfo;
423  }
424 
425  //if (debug & 2)
426  //{
427  // Pout<< "Processor patch " << patchi << ' ' << procPatch.name()
428  // << " recv:" << patchInfo.size() << " from proc:"
429  // << nbrProci << endl;
430  //}
431 
432  // Apply transform to received data for non-parallel planes
433  if (!procPatch.parallel())
434  {
435  transform(procPatch, procPatch.forwardT(), patchInfo);
436  }
437 
438  // Adapt for entering domain
439  enterDomain(procPatch, patchPoints, patchInfo);
440 
441  // Merge received info
442  const labelList& meshPoints = procPatch.meshPoints();
443  forAll(patchInfo, i)
444  {
445  label meshPointi = meshPoints[patchPoints[i]];
446 
447  if (!allPointInfo_[meshPointi].equal(patchInfo[i], td_))
448  {
449  updatePoint
450  (
451  meshPointi,
452  patchInfo[i],
453  allPointInfo_[meshPointi]
454  );
455  }
456  }
457  }
458 
459  // Collocated points should be handled by face based transfer
460  // (since that is how connectivity is worked out)
461  // They are also explicitly equalised in handleCollocatedPoints to
462  // guarantee identical values.
463 }
464 
465 
466 template<class Type, class TrackingData>
468 {
469  // 1. Send all point info on cyclic patches.
470 
471  DynamicList<Type> nbrInfo;
472  DynamicList<label> nbrPoints;
473  DynamicList<label> thisPoints;
474 
475  for (const polyPatch& patch : mesh_.boundaryMesh())
476  {
477  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(patch);
478 
479  if (cpp)
480  {
481  const auto& cycPatch = *cpp;
482  const auto& nbrPatch = cycPatch.neighbPatch();
483 
484  nbrInfo.clear();
485  nbrInfo.reserve(cycPatch.nPoints());
486  nbrPoints.clear();
487  nbrPoints.reserve(cycPatch.nPoints());
488  thisPoints.clear();
489  thisPoints.reserve(cycPatch.nPoints());
490 
491  // Collect nbrPatch points that have changed
492  {
493  const edgeList& pairs = cycPatch.coupledPoints();
494  const labelList& meshPoints = nbrPatch.meshPoints();
495 
496  forAll(pairs, pairI)
497  {
498  label thisPointi = pairs[pairI][0];
499  label nbrPointi = pairs[pairI][1];
500  label meshPointi = meshPoints[nbrPointi];
501 
502  if (changedPoint_.test(meshPointi))
503  {
504  nbrInfo.append(allPointInfo_[meshPointi]);
505  nbrPoints.append(nbrPointi);
506  thisPoints.append(thisPointi);
507  }
508  }
509 
510  // nbr : adapt for leaving domain
511  leaveDomain(nbrPatch, nbrPoints, nbrInfo);
512  }
513 
514 
515  // Apply rotation for non-parallel planes
516 
517  if (!cycPatch.parallel())
518  {
519  // received data from half1
520  transform(cycPatch, cycPatch.forwardT(), nbrInfo);
521  }
522 
523  //if (debug)
524  //{
525  // Pout<< "Cyclic patch " << patch.index() << ' ' << patch.name()
526  // << " Changed : " << nbrInfo.size()
527  // << endl;
528  //}
529 
530  // Adapt for entering domain
531  enterDomain(cycPatch, thisPoints, nbrInfo);
532 
533  // Merge received info
534  const labelList& meshPoints = cycPatch.meshPoints();
535  forAll(nbrInfo, i)
536  {
537  label meshPointi = meshPoints[thisPoints[i]];
538 
539  if (!allPointInfo_[meshPointi].equal(nbrInfo[i], td_))
540  {
541  updatePoint
542  (
543  meshPointi,
544  nbrInfo[i],
545  allPointInfo_[meshPointi]
546  );
547  }
548  }
549  }
550  }
551 }
552 
553 
554 // Guarantee collocated points have same information.
555 // Return number of points changed.
556 template<class Type, class TrackingData>
558 {
559  // Transfer onto coupled patch
560  const globalMeshData& gmd = mesh_.globalData();
561  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
562  const labelList& meshPoints = cpp.meshPoints();
563 
564  const mapDistribute& slavesMap = gmd.globalPointSlavesMap();
565  const labelListList& slaves = gmd.globalPointSlaves();
566 
567  List<Type> elems(slavesMap.constructSize());
568  forAll(meshPoints, pointi)
569  {
570  elems[pointi] = allPointInfo_[meshPoints[pointi]];
571  }
572 
573  // Pull slave data onto master (which might or might not have any
574  // initialised points). No need to update transformed slots.
575  slavesMap.distribute(elems, false);
576 
577  // Combine master data with slave data
578  combineEqOp<Type, TrackingData> cop(td_);
579 
580  forAll(slaves, pointi)
581  {
582  Type& elem = elems[pointi];
583 
584  const labelList& slavePoints = slaves[pointi];
585 
586  // Combine master with untransformed slave data
587  forAll(slavePoints, j)
588  {
589  cop(elem, elems[slavePoints[j]]);
590  }
591 
592  // Copy result back to slave slots
593  forAll(slavePoints, j)
594  {
595  elems[slavePoints[j]] = elem;
596  }
597  }
598 
599  // Push slave-slot data back to slaves
600  slavesMap.reverseDistribute(elems.size(), elems, false);
601 
602  // Extract back onto mesh
603  forAll(meshPoints, pointi)
604  {
605  if (elems[pointi].valid(td_))
606  {
607  label meshPointi = meshPoints[pointi];
608 
609  Type& elem = allPointInfo_[meshPointi];
610 
611  bool wasValid = elem.valid(td_);
612 
613  // Like updatePoint but bypass Type::updatePoint with its tolerance
614  // checking
615  //if (!elem.valid(td_) || !elem.equal(elems[pointi], td_))
616  if (!elem.equal(elems[pointi], td_))
617  {
618  nEvals_++;
619  elem = elems[pointi];
620 
621  // See if element now valid
622  if (!wasValid && elem.valid(td_))
623  {
624  --nUnvisitedPoints_;
625  }
626 
627  // Update database of changed points
628  if (changedPoint_.set(meshPointi))
629  {
630  changedPoints_.push_back(meshPointi);
631  }
632  }
633  }
634  }
635 
636  // Sum changedPoints over all procs
637  return returnReduce(nChangedPoints(), sumOp<label>());
638 }
639 
640 
641 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
643 // Iterate, propagating changedPointsInfo across mesh, until no change (or
644 // maxIter reached). Initial point values specified.
645 template<class Type, class TrackingData>
647 (
648  const polyMesh& mesh,
649  const labelList& changedPoints,
650  const List<Type>& changedPointsInfo,
651 
652  UList<Type>& allPointInfo,
653  UList<Type>& allEdgeInfo,
654 
655  const label maxIter,
656  TrackingData& td
657 )
658 :
660 
661  allPointInfo_(allPointInfo),
662  allEdgeInfo_(allEdgeInfo),
663  td_(td),
664  nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
665  nEvals_(0)
666 {
667  if (allPointInfo_.size() != mesh_.nPoints())
668  {
670  << "size of pointInfo work array is not equal to the number"
671  << " of points in the mesh" << endl
672  << " pointInfo :" << allPointInfo_.size() << endl
673  << " mesh.nPoints:" << mesh_.nPoints()
674  << exit(FatalError);
675  }
676  if (allEdgeInfo_.size() != mesh_.nEdges())
677  {
679  << "size of edgeInfo work array is not equal to the number"
680  << " of edges in the mesh" << endl
681  << " edgeInfo :" << allEdgeInfo_.size() << endl
682  << " mesh.nEdges:" << mesh_.nEdges()
683  << exit(FatalError);
684  }
685 
686 
687  // Set from initial changed points data
688  setPointInfo(changedPoints, changedPointsInfo);
689 
690  if (debug)
691  {
692  Info<< typeName << ": Seed points : "
694  }
695 
696  // Iterate until nothing changes
697  label iter = iterate(maxIter);
698 
699  if ((maxIter > 0) && (iter >= maxIter))
700  {
702  << "Maximum number of iterations reached. Increase maxIter." << endl
703  << " maxIter:" << maxIter << nl
704  << " nChangedPoints:" << nChangedPoints() << nl
705  << " nChangedEdges:" << nChangedEdges() << endl
706  << exit(FatalError);
707  }
708 }
709 
710 
711 template<class Type, class TrackingData>
713 (
714  const polyMesh& mesh,
715  UList<Type>& allPointInfo,
716  UList<Type>& allEdgeInfo,
717  TrackingData& td
718 )
719 :
721 
722  allPointInfo_(allPointInfo),
723  allEdgeInfo_(allEdgeInfo),
724  td_(td),
725  nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
726  nEvals_(0)
727 {}
728 
729 
730 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
731 
732 // Copy point information into member data
733 template<class Type, class TrackingData>
735 (
736  const labelList& changedPoints,
737  const List<Type>& changedPointsInfo
738 )
739 {
740  forAll(changedPoints, changedPointi)
741  {
742  const label pointi = changedPoints[changedPointi];
743 
744  const bool wasValid = allPointInfo_[pointi].valid(td_);
745 
746  // Copy info for pointi
747  allPointInfo_[pointi] = changedPointsInfo[changedPointi];
748 
749  // Maintain count of unset points
750  if (!wasValid && allPointInfo_[pointi].valid(td_))
751  {
752  --nUnvisitedPoints_;
753  }
754 
755  // Mark pointi as changed, both on list and on point itself.
756 
757  if (changedPoint_.set(pointi))
758  {
759  changedPoints_.push_back(pointi);
760  }
761  }
762 
763  // Sync
764  handleCollocatedPoints();
765 }
766 
767 
768 // Propagate information from edge to point. Return number of points changed.
769 template<class Type, class TrackingData>
771 {
772  for (const label edgei : changedEdges_)
773  {
774  if (!changedEdge_.test(edgei))
775  {
777  << "edge " << edgei
778  << " not marked as having been changed" << nl
779  << "This might be caused by multiple occurrences of the same"
780  << " seed point." << abort(FatalError);
781  }
782 
783 
784  const Type& neighbourWallInfo = allEdgeInfo_[edgei];
785 
786  // Evaluate all connected points (= edge endpoints)
787  const edge& e = mesh_.edges()[edgei];
788 
789  forAll(e, eI)
790  {
791  Type& currentWallInfo = allPointInfo_[e[eI]];
792 
793  if (!currentWallInfo.equal(neighbourWallInfo, td_))
794  {
795  updatePoint
796  (
797  e[eI],
798  edgei,
799  neighbourWallInfo,
800  currentWallInfo
801  );
802  }
803  }
804 
805  // Reset status of edge
806  changedEdge_.unset(edgei);
807  }
808 
809  // Handled all changed edges by now
810  changedEdges_.clear();
811 
812  if (nCyclicPatches_ > 0)
813  {
814  // Transfer changed points across cyclic halves
815  handleCyclicPatches();
816  }
817  if (Pstream::parRun())
818  {
819  // Transfer changed points from neighbouring processors.
820  handleProcPatches();
821  }
822 
823  //if (debug)
824  //{
825  // Pout<< "Changed points : " << nChangedPoints() << nl;
826  //}
827 
828  // Sum changedPoints over all procs
829  return returnReduce(nChangedPoints(), sumOp<label>());
830 }
831 
832 
833 // Propagate information from point to edge. Return number of edges changed.
834 template<class Type, class TrackingData>
836 {
837  const labelListList& pointEdges = mesh_.pointEdges();
838 
839  for (const label pointi : changedPoints_)
840  {
841  if (!changedPoint_.test(pointi))
842  {
844  << "Point " << pointi
845  << " not marked as having been changed" << nl
846  << "This might be caused by multiple occurrences of the same"
847  << " seed point." << abort(FatalError);
848  }
849 
850  const Type& neighbourWallInfo = allPointInfo_[pointi];
851 
852  // Evaluate all connected edges
853 
854  for (const label edgei : pointEdges[pointi])
855  {
856  Type& currentWallInfo = allEdgeInfo_[edgei];
857 
858  if (!currentWallInfo.equal(neighbourWallInfo, td_))
859  {
860  updateEdge
861  (
862  edgei,
863  pointi,
864  neighbourWallInfo,
865  currentWallInfo
866  );
867  }
868  }
869 
870  // Reset status of point
871  changedPoint_.unset(pointi);
872  }
873 
874  // Handled all changed points by now
875  changedPoints_.clear();
876 
877  //if (debug)
878  //{
879  // Pout<< "Changed edges : " << nChangedEdges() << endl;
880  //}
881 
882  // Sum changedEdges over all procs
883  return returnReduce(nChangedEdges(), sumOp<label>());
884 }
886 
887 // Iterate
888 template<class Type, class TrackingData>
890 (
891  const label maxIter
892 )
893 {
894  if (nCyclicPatches_ > 0)
895  {
896  // Transfer changed points across cyclic halves
897  handleCyclicPatches();
898  }
899  if (Pstream::parRun())
900  {
901  // Transfer changed points from neighbouring processors.
902  handleProcPatches();
903  }
904 
905  nEvals_ = 0;
906 
907  label iter = 0;
908 
909  while (iter < maxIter)
910  {
911  while (iter < maxIter)
912  {
913  if (debug)
914  {
915  Info<< typeName << ": Iteration " << iter << endl;
916  }
917 
918  label nEdges = pointToEdge();
919 
920  if (debug)
921  {
922  Info<< typeName << ": Total changed edges : "
923  << nEdges << endl;
924  }
925 
926  if (nEdges == 0)
927  {
928  break;
929  }
930 
931  label nPoints = edgeToPoint();
932 
933  if (debug)
934  {
935  Info<< typeName << ": Total changed points : "
936  << nPoints << nl
937  << typeName << ": Total evaluations : "
938  << returnReduce(nEvals_, sumOp<label>()) << nl
939  << typeName << ": Remaining unvisited points: "
940  << returnReduce(nUnvisitedPoints_, sumOp<label>()) << nl
941  << typeName << ": Remaining unvisited edges : "
942  << returnReduce(nUnvisitedEdges_, sumOp<label>()) << nl
943  << endl;
944  }
945 
946  if (nPoints == 0)
947  {
948  break;
949  }
950 
951  iter++;
952  }
953 
954 
955  // Enforce collocated points are exactly equal. This might still mean
956  // non-collocated points are not equal though. WIP.
957  label nPoints = handleCollocatedPoints();
958  if (debug)
959  {
960  Info<< typeName << ": Collocated point sync : "
961  << nPoints << nl << endl;
962  }
963 
964  if (nPoints == 0)
965  {
966  break;
967  }
968  }
969 
970  return iter;
971 }
972 
973 
974 // ************************************************************************* //
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:578
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:204
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:49
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:632
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1004
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:185
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:414
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:109
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:570
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:331
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:389
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:73
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.