FaceCellWave.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "FaceCellWave.H"
27 #include "polyMesh.H"
28 #include "processorPolyPatch.H"
29 #include "cyclicPolyPatch.H"
30 #include "cyclicAMIPolyPatch.H"
31 #include "OPstream.H"
32 #include "IPstream.H"
33 #include "PstreamReduceOps.H"
34 #include "debug.H"
35 #include "typeInfo.H"
36 #include "SubField.H"
37 #include "globalMeshData.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 template<class Type, class TrackingData>
43 
44 template<class Type, class TrackingData>
46 
47 template<class Type, class TrackingData>
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
53 template<class Type, class TrackingData>
55 (
56  const label celli,
57  const label neighbourFacei,
58  const Type& neighbourInfo,
59  const scalar tol,
60  Type& cellInfo
61 )
62 {
63  // Update info for celli, at position pt, with information from
64  // neighbouring face/cell.
65  // Updates:
66  // - changedCell_, changedCells_
67  // - statistics: nEvals_, nUnvisitedCells_
68 
69  nEvals_++;
70 
71  bool wasValid = cellInfo.valid(td_);
72 
73  bool propagate =
74  cellInfo.updateCell
75  (
76  mesh_,
77  celli,
78  neighbourFacei,
79  neighbourInfo,
80  tol,
81  td_
82  );
83 
84  if (propagate)
85  {
86  if (!changedCell_[celli])
87  {
88  changedCell_[celli] = true;
89  changedCells_.append(celli);
90  }
91  }
92 
93  if (!wasValid && cellInfo.valid(td_))
94  {
95  --nUnvisitedCells_;
96  }
97 
98  return propagate;
99 }
100 
101 
102 template<class Type, class TrackingData>
104 (
105  const label facei,
106  const label neighbourCelli,
107  const Type& neighbourInfo,
108  const scalar tol,
109  Type& faceInfo
110 )
111 {
112  // Update info for facei, at position pt, with information from
113  // neighbouring face/cell.
114  // Updates:
115  // - changedFace_, changedFaces_,
116  // - statistics: nEvals_, nUnvisitedFaces_
117 
118  nEvals_++;
119 
120  bool wasValid = faceInfo.valid(td_);
121 
122  bool propagate =
123  faceInfo.updateFace
124  (
125  mesh_,
126  facei,
127  neighbourCelli,
128  neighbourInfo,
129  tol,
130  td_
131  );
132 
133  if (propagate)
134  {
135  if (!changedFace_[facei])
136  {
137  changedFace_[facei] = true;
138  changedFaces_.append(facei);
139  }
140  }
141 
142  if (!wasValid && faceInfo.valid(td_))
143  {
144  --nUnvisitedFaces_;
145  }
146 
147  return propagate;
148 }
149 
150 
151 template<class Type, class TrackingData>
153 (
154  const label facei,
155  const Type& neighbourInfo,
156  const scalar tol,
157  Type& faceInfo
158 )
159 {
160  // Update info for facei, at position pt, with information from
161  // same face.
162  // Updates:
163  // - changedFace_, changedFaces_,
164  // - statistics: nEvals_, nUnvisitedFaces_
165 
166  nEvals_++;
167 
168  bool wasValid = faceInfo.valid(td_);
169 
170  bool propagate =
171  faceInfo.updateFace
172  (
173  mesh_,
174  facei,
175  neighbourInfo,
176  tol,
177  td_
178  );
179 
180  if (propagate)
181  {
182  if (!changedFace_[facei])
183  {
184  changedFace_[facei] = true;
185  changedFaces_.append(facei);
186  }
187  }
188 
189  if (!wasValid && faceInfo.valid(td_))
190  {
191  --nUnvisitedFaces_;
192  }
193 
194  return propagate;
195 }
196 
197 
198 template<class Type, class TrackingData>
200 (
201  const polyPatch& patch
202 ) const
203 {
204  // For debugging: check status on both sides of cyclic
205 
206  const cyclicPolyPatch& nbrPatch =
207  refCast<const cyclicPolyPatch>(patch).nbrPatch();
208 
209  forAll(patch, patchFacei)
210  {
211  label i1 = patch.start() + patchFacei;
212  label i2 = nbrPatch.start() + patchFacei;
213 
214  if
215  (
216  !allFaceInfo_[i1].sameGeometry
217  (
218  mesh_,
219  allFaceInfo_[i2],
220  geomTol_,
221  td_
222  )
223  )
224  {
226  << " faceInfo:" << allFaceInfo_[i1]
227  << " otherfaceInfo:" << allFaceInfo_[i2]
228  << abort(FatalError);
229  }
230 
231  if (changedFace_[i1] != changedFace_[i2])
232  {
234  << " faceInfo:" << allFaceInfo_[i1]
235  << " otherfaceInfo:" << allFaceInfo_[i2]
236  << " changedFace:" << changedFace_[i1]
237  << " otherchangedFace:" << changedFace_[i2]
238  << abort(FatalError);
239  }
240  }
241 }
242 
243 
244 template<class Type, class TrackingData>
245 template<class PatchType>
247 {
248  forAll(mesh_.boundaryMesh(), patchi)
249  {
250  if (isA<PatchType>(mesh_.boundaryMesh()[patchi]))
251  {
252  return true;
253  }
254  }
255  return false;
256 }
257 
258 
259 template<class Type, class TrackingData>
261 (
262  const labelList& changedFaces,
263  const List<Type>& changedFacesInfo
264 )
265 {
266  forAll(changedFaces, changedFacei)
267  {
268  label facei = changedFaces[changedFacei];
269 
270  bool wasValid = allFaceInfo_[facei].valid(td_);
271 
272  // Copy info for facei
273  allFaceInfo_[facei] = changedFacesInfo[changedFacei];
274 
275  // Maintain count of unset faces
276  if (!wasValid && allFaceInfo_[facei].valid(td_))
277  {
278  --nUnvisitedFaces_;
279  }
280 
281  // Mark facei as changed, both on list and on face itself.
282 
283  changedFace_[facei] = true;
284  changedFaces_.append(facei);
285  }
286 }
287 
288 
289 template<class Type, class TrackingData>
291 (
292  const polyPatch& patch,
293  const label nFaces,
294  const labelList& changedFaces,
295  const List<Type>& changedFacesInfo
296 )
297 {
298  // Merge face information into member data
299 
300  for (label changedFacei = 0; changedFacei < nFaces; changedFacei++)
301  {
302  const Type& neighbourWallInfo = changedFacesInfo[changedFacei];
303  label patchFacei = changedFaces[changedFacei];
304 
305  label meshFacei = patch.start() + patchFacei;
306 
307  Type& currentWallInfo = allFaceInfo_[meshFacei];
308 
309  if (!currentWallInfo.equal(neighbourWallInfo, td_))
310  {
311  updateFace
312  (
313  meshFacei,
314  neighbourWallInfo,
315  propagationTol_,
316  currentWallInfo
317  );
318  }
319  }
320 }
321 
322 
323 template<class Type, class TrackingData>
325 (
326  const polyPatch& patch,
327  const label startFacei,
328  const label nFaces,
329  labelList& changedPatchFaces,
330  List<Type>& changedPatchFacesInfo
331 ) const
332 {
333  // Construct compact patchFace change arrays for a (slice of a) single
334  // patch. changedPatchFaces in local patch numbering.
335  // Return length of arrays.
336  label nChangedPatchFaces = 0;
337 
338  for (label i = 0; i < nFaces; i++)
339  {
340  label patchFacei = i + startFacei;
341 
342  label meshFacei = patch.start() + patchFacei;
343 
344  if (changedFace_[meshFacei])
345  {
346  changedPatchFaces[nChangedPatchFaces] = patchFacei;
347  changedPatchFacesInfo[nChangedPatchFaces] = allFaceInfo_[meshFacei];
348  nChangedPatchFaces++;
349  }
350  }
351  return nChangedPatchFaces;
352 }
353 
354 
355 template<class Type, class TrackingData>
357 (
358  const polyPatch& patch,
359  const label nFaces,
360  const labelList& faceLabels,
361  const transformer& transform,
362  List<Type>& faceInfo
363 )
364 {
365  for (label i = 0; i < nFaces; i++)
366  {
367  faceInfo[i].transform(patch, faceLabels[i], transform, td_);
368  }
369 }
370 
371 
372 template<class Type, class TrackingData>
374 {
375  // Transfer all the information to/from neighbouring processors
376 
377  const globalMeshData& pData = mesh_.globalData();
378 
379  // Which patches are processor patches
380  const labelList& procPatches = pData.processorPatches();
381 
382  // Send all
383 
384  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
385 
386  forAll(procPatches, i)
387  {
388  label patchi = procPatches[i];
389 
390  const processorPolyPatch& procPatch =
391  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
392 
393  // Allocate buffers
394  label nSendFaces;
395  labelList sendFaces(procPatch.size());
396  List<Type> sendFacesInfo(procPatch.size());
397 
398  // Determine which faces changed on current patch
399  nSendFaces = getChangedPatchFaces
400  (
401  procPatch,
402  0,
403  procPatch.size(),
404  sendFaces,
405  sendFacesInfo
406  );
407 
408  if (debug & 2)
409  {
410  Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
411  << " communicating with " << procPatch.neighbProcNo()
412  << " Sending:" << nSendFaces
413  << endl;
414  }
415 
416  // Send
417  UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
418  toNeighbour
419  << SubList<label>(sendFaces, nSendFaces)
420  << SubList<Type>(sendFacesInfo, nSendFaces);
421  }
422 
423  pBufs.finishedSends();
424 
425  // Receive all
426 
427  forAll(procPatches, i)
428  {
429  label patchi = procPatches[i];
430 
431  const processorPolyPatch& procPatch =
432  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
433 
434  // Allocate buffers
435  labelList receiveFaces;
436  List<Type> receiveFacesInfo;
437 
438  // Receive
439  {
440  UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
441  fromNeighbour >> receiveFaces >> receiveFacesInfo;
442  }
443 
444  if (debug & 2)
445  {
446  Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
447  << " communicating with " << procPatch.neighbProcNo()
448  << " Receiving:" << receiveFaces.size()
449  << endl;
450  }
451 
452  // Transform info across the interface
453  transform
454  (
455  procPatch,
456  receiveFaces.size(),
457  receiveFaces,
458  procPatch.transform(),
459  receiveFacesInfo
460  );
461 
462  // Merge received info
463  mergeFaceInfo
464  (
465  procPatch,
466  receiveFaces.size(),
467  receiveFaces,
468  receiveFacesInfo
469  );
470  }
471 }
472 
473 
474 template<class Type, class TrackingData>
476 {
477  // Transfer information across cyclic halves.
478 
479  forAll(mesh_.boundaryMesh(), patchi)
480  {
481  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
482 
483  if (isA<cyclicPolyPatch>(patch))
484  {
485  const cyclicPolyPatch& cycPatch =
486  refCast<const cyclicPolyPatch>(patch);
487  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
488 
489  // Allocate buffers
490  label nReceiveFaces;
491  labelList receiveFaces(patch.size());
492  List<Type> receiveFacesInfo(patch.size());
493 
494  // Determine which faces changed
495  nReceiveFaces = getChangedPatchFaces
496  (
497  nbrPatch,
498  0,
499  nbrPatch.size(),
500  receiveFaces,
501  receiveFacesInfo
502  );
503 
504  if (debug & 2)
505  {
506  Pout<< " Cyclic patch " << patchi << ' ' << cycPatch.name()
507  << " Changed : " << nReceiveFaces
508  << endl;
509  }
510 
511  // Transform info across the interface
512  transform
513  (
514  cycPatch,
515  nReceiveFaces,
516  receiveFaces,
517  cycPatch.transform(),
518  receiveFacesInfo
519  );
520 
521  // Merge into global storage
522  mergeFaceInfo
523  (
524  cycPatch,
525  nReceiveFaces,
526  receiveFaces,
527  receiveFacesInfo
528  );
529 
530  if (debug)
531  {
532  checkCyclic(cycPatch);
533  }
534  }
535  }
536 }
537 
538 
539 template<class Type, class TrackingData>
541 {
542  // Define combine operator for AMIInterpolation
543 
544  class combine
545  {
547 
548  const cyclicAMIPolyPatch& patch_;
549 
550  public:
551 
552  combine
553  (
555  const cyclicAMIPolyPatch& patch
556  )
557  :
558  solver_(solver),
559  patch_(patch)
560  {}
561 
562  void operator()
563  (
564  Type& x,
565  const label facei,
566  const Type& y,
567  const scalar weight
568  ) const
569  {
570  if (y.valid(solver_.data()))
571  {
572  label meshFacei = -1;
573 
574  if (patch_.owner())
575  {
576  meshFacei = patch_.start() + facei;
577  }
578  else
579  {
580  meshFacei = patch_.nbrPatch().start() + facei;
581  }
582 
583  x.updateFace
584  (
585  solver_.mesh(),
586  meshFacei,
587  y,
588  solver_.propagationTol(),
589  solver_.data()
590  );
591  }
592  }
593  };
594 
595  // Transfer information across cyclicAMI halves.
596 
597  forAll(mesh_.boundaryMesh(), patchi)
598  {
599  const polyPatch& patch = mesh_.boundaryMesh()[patchi];
600 
601  if (isA<cyclicAMIPolyPatch>(patch))
602  {
603  const cyclicAMIPolyPatch& cycPatch =
604  refCast<const cyclicAMIPolyPatch>(patch);
605  const cyclicAMIPolyPatch& nbrPatch = cycPatch.nbrPatch();
606 
607  // Create a combine operator to transfer sendInfo from nbrPatch
608  // to cycPatch
609  combine cmb(*this, cycPatch);
610 
611  // Get nbrPatch data (so not just changed faces)
612  List<Type> sendInfo(nbrPatch.patchSlice(allFaceInfo_));
613 
614  // Construct default values, if necessary
615  List<Type> defVals;
616  if (cycPatch.applyLowWeightCorrection())
617  {
618  defVals = cycPatch.patchInternalList(allCellInfo_);
619  }
620 
621  // Interpolate from nbrPatch to cycPatch, applying AMI
622  // transforms as necessary
623  List<Type> receiveInfo;
624  if (cycPatch.owner())
625  {
626  forAll(cycPatch.AMIs(), i)
627  {
628  if (cycPatch.AMITransforms()[i].transformsPosition())
629  {
630  List<Type> sendInfoT(sendInfo);
631  transform
632  (
633  nbrPatch,
634  nbrPatch.size(),
635  identity(nbrPatch.size()),
636  cycPatch.AMITransforms()[i],
637  sendInfo
638  );
639  cycPatch.AMIs()[i].interpolateToSource
640  (
641  sendInfoT,
642  cmb,
643  receiveInfo,
644  defVals
645  );
646  }
647  else
648  {
649  cycPatch.AMIs()[i].interpolateToSource
650  (
651  sendInfo,
652  cmb,
653  receiveInfo,
654  defVals
655  );
656  }
657  }
658  }
659  else
660  {
661  forAll(nbrPatch.AMIs(), i)
662  {
663  if (nbrPatch.AMITransforms()[i].transformsPosition())
664  {
665  List<Type> sendInfoT(sendInfo);
666  transform
667  (
668  nbrPatch,
669  nbrPatch.size(),
670  identity(nbrPatch.size()),
671  nbrPatch.AMITransforms()[i],
672  sendInfo
673  );
674  cycPatch.nbrPatch().AMIs()[i].interpolateToTarget
675  (
676  sendInfoT,
677  cmb,
678  receiveInfo,
679  defVals
680  );
681  }
682  else
683  {
684  cycPatch.nbrPatch().AMIs()[i].interpolateToTarget
685  (
686  sendInfo,
687  cmb,
688  receiveInfo,
689  defVals
690  );
691  }
692  }
693  }
694 
695  // Shuffle up to remove invalid info
696  DynamicList<label> receiveFaces(cycPatch.size());
697  forAll(receiveInfo, patchFacei)
698  {
699  if (receiveInfo[patchFacei].valid(td_))
700  {
701  if (receiveFaces.size() != patchFacei)
702  {
703  receiveInfo[receiveFaces.size()] =
704  receiveInfo[patchFacei];
705  }
706  receiveFaces.append(patchFacei);
707  }
708  }
709 
710  // Transform info across the interface
711  if (cycPatch.transform().transformsPosition())
712  {
713  transform
714  (
715  cycPatch,
716  receiveFaces.size(),
717  receiveFaces,
718  cycPatch.transform(),
719  receiveInfo
720  );
721  }
722 
723  // Merge into global storage
724  mergeFaceInfo
725  (
726  cycPatch,
727  receiveFaces.size(),
728  receiveFaces,
729  receiveInfo
730  );
731  }
732  }
733 }
734 
735 
736 template<class Type, class TrackingData>
738 {
739  // Collect changed information
740 
741  DynamicList<label> f0Baffle(explicitConnections_.size());
742  DynamicList<Type> f0Info(explicitConnections_.size());
743 
744  DynamicList<label> f1Baffle(explicitConnections_.size());
745  DynamicList<Type> f1Info(explicitConnections_.size());
746 
747  forAll(explicitConnections_, connI)
748  {
749  const labelPair& baffle = explicitConnections_[connI];
750 
751  label f0 = baffle[0];
752  if (changedFace_[f0])
753  {
754  f0Baffle.append(connI);
755  f0Info.append(allFaceInfo_[f0]);
756  }
757 
758  label f1 = baffle[1];
759  if (changedFace_[f1])
760  {
761  f1Baffle.append(connI);
762  f1Info.append(allFaceInfo_[f1]);
763  }
764  }
765 
766 
767  // Update other side with changed information
768 
769  forAll(f1Info, index)
770  {
771  const labelPair& baffle = explicitConnections_[f1Baffle[index]];
772 
773  label f0 = baffle[0];
774  Type& currentWallInfo = allFaceInfo_[f0];
775 
776  if (!currentWallInfo.equal(f1Info[index], td_))
777  {
778  updateFace
779  (
780  f0,
781  f1Info[index],
782  propagationTol_,
783  currentWallInfo
784  );
785  }
786  }
787 
788  forAll(f0Info, index)
789  {
790  const labelPair& baffle = explicitConnections_[f0Baffle[index]];
791 
792  label f1 = baffle[1];
793  Type& currentWallInfo = allFaceInfo_[f1];
794 
795  if (!currentWallInfo.equal(f0Info[index], td_))
796  {
797  updateFace
798  (
799  f1,
800  f0Info[index],
801  propagationTol_,
802  currentWallInfo
803  );
804  }
805  }
806 }
807 
808 
809 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
810 
811 template<class Type, class TrackingData>
813 (
814  const polyMesh& mesh,
815  UList<Type>& allFaceInfo,
816  UList<Type>& allCellInfo,
817  TrackingData& td
818 )
819 :
820  mesh_(mesh),
821  explicitConnections_(0),
822  allFaceInfo_(allFaceInfo),
823  allCellInfo_(allCellInfo),
824  td_(td),
825  changedFace_(mesh_.nFaces(), false),
826  changedFaces_(mesh_.nFaces()),
827  changedCell_(mesh_.nCells(), false),
828  changedCells_(mesh_.nCells()),
829  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
830  hasCyclicAMIPatches_
831  (
832  returnReduce(hasPatch<cyclicAMIPolyPatch>(), orOp<bool>())
833  ),
834  nEvals_(0),
835  nUnvisitedCells_(mesh_.nCells()),
836  nUnvisitedFaces_(mesh_.nFaces())
837 {
838  if
839  (
840  allFaceInfo.size() != mesh_.nFaces()
841  || allCellInfo.size() != mesh_.nCells()
842  )
843  {
845  << "face and cell storage not the size of mesh faces, cells:"
846  << endl
847  << " allFaceInfo :" << allFaceInfo.size() << endl
848  << " mesh_.nFaces():" << mesh_.nFaces() << endl
849  << " allCellInfo :" << allCellInfo.size() << endl
850  << " mesh_.nCells():" << mesh_.nCells()
851  << exit(FatalError);
852  }
853 }
854 
855 
856 template<class Type, class TrackingData>
858 (
859  const polyMesh& mesh,
860  const labelList& changedFaces,
861  const List<Type>& changedFacesInfo,
862  UList<Type>& allFaceInfo,
863  UList<Type>& allCellInfo,
864  const label maxIter,
865  TrackingData& td
866 )
867 :
868  mesh_(mesh),
869  explicitConnections_(0),
870  allFaceInfo_(allFaceInfo),
871  allCellInfo_(allCellInfo),
872  td_(td),
873  changedFace_(mesh_.nFaces(), false),
874  changedFaces_(mesh_.nFaces()),
875  changedCell_(mesh_.nCells(), false),
876  changedCells_(mesh_.nCells()),
877  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
878  hasCyclicAMIPatches_
879  (
880  returnReduce(hasPatch<cyclicAMIPolyPatch>(), orOp<bool>())
881  ),
882  nEvals_(0),
883  nUnvisitedCells_(mesh_.nCells()),
884  nUnvisitedFaces_(mesh_.nFaces())
885 {
886  if
887  (
888  allFaceInfo.size() != mesh_.nFaces()
889  || allCellInfo.size() != mesh_.nCells()
890  )
891  {
893  << "face and cell storage not the size of mesh faces, cells:"
894  << endl
895  << " allFaceInfo :" << allFaceInfo.size() << endl
896  << " mesh_.nFaces():" << mesh_.nFaces() << endl
897  << " allCellInfo :" << allCellInfo.size() << endl
898  << " mesh_.nCells():" << mesh_.nCells()
899  << exit(FatalError);
900  }
901 
902  // Copy initial changed faces data
903  setFaceInfo(changedFaces, changedFacesInfo);
904 
905  // Iterate until nothing changes
906  label iter = iterate(maxIter);
907 
908  if ((maxIter > 0) && (iter >= maxIter))
909  {
911  << "Maximum number of iterations reached. Increase maxIter." << endl
912  << " maxIter:" << maxIter << endl
913  << " nChangedCells:" << changedCells_.size() << endl
914  << " nChangedFaces:" << changedFaces_.size() << endl
915  << exit(FatalError);
916  }
917 }
918 
919 
920 template<class Type, class TrackingData>
922 (
923  const polyMesh& mesh,
924  const labelPairList& explicitConnections,
925  const bool handleCyclicAMI,
926  const labelList& changedFaces,
927  const List<Type>& changedFacesInfo,
928  UList<Type>& allFaceInfo,
929  UList<Type>& allCellInfo,
930  const label maxIter,
931  TrackingData& td
932 )
933 :
934  mesh_(mesh),
935  explicitConnections_(explicitConnections),
936  allFaceInfo_(allFaceInfo),
937  allCellInfo_(allCellInfo),
938  td_(td),
939  changedFace_(mesh_.nFaces(), false),
940  changedFaces_(mesh_.nFaces()),
941  changedCell_(mesh_.nCells(), false),
942  changedCells_(mesh_.nCells()),
943  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
944  hasCyclicAMIPatches_
945  (
946  handleCyclicAMI
947  && returnReduce(hasPatch<cyclicAMIPolyPatch>(), orOp<bool>())
948  ),
949  nEvals_(0),
950  nUnvisitedCells_(mesh_.nCells()),
951  nUnvisitedFaces_(mesh_.nFaces())
952 {
953  if
954  (
955  allFaceInfo.size() != mesh_.nFaces()
956  || allCellInfo.size() != mesh_.nCells()
957  )
958  {
960  << "face and cell storage not the size of mesh faces, cells:"
961  << endl
962  << " allFaceInfo :" << allFaceInfo.size() << endl
963  << " mesh_.nFaces():" << mesh_.nFaces() << endl
964  << " allCellInfo :" << allCellInfo.size() << endl
965  << " mesh_.nCells():" << mesh_.nCells()
966  << exit(FatalError);
967  }
968 
969  // Copy initial changed faces data
970  setFaceInfo(changedFaces, changedFacesInfo);
971 
972  // Iterate until nothing changes
973  label iter = iterate(maxIter);
974 
975  if ((maxIter > 0) && (iter >= maxIter))
976  {
978  << "Maximum number of iterations reached. Increase maxIter." << endl
979  << " maxIter:" << maxIter << endl
980  << " nChangedCells:" << changedCells_.size() << endl
981  << " nChangedFaces:" << changedFaces_.size() << endl
982  << exit(FatalError);
983  }
984 }
985 
986 
987 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
988 
989 template<class Type, class TrackingData>
991 {
992  return nUnvisitedCells_;
993 }
994 
995 
996 template<class Type, class TrackingData>
998 {
999  return nUnvisitedFaces_;
1000 }
1001 
1002 
1003 template<class Type, class TrackingData>
1005 {
1006  // Propagate face to cell
1007 
1008  const labelList& owner = mesh_.faceOwner();
1009  const labelList& neighbour = mesh_.faceNeighbour();
1010  label nInternalFaces = mesh_.nInternalFaces();
1011 
1012  forAll(changedFaces_, changedFacei)
1013  {
1014  label facei = changedFaces_[changedFacei];
1015  if (!changedFace_[facei])
1016  {
1018  << "Face " << facei
1019  << " not marked as having been changed"
1020  << abort(FatalError);
1021  }
1022 
1023 
1024  const Type& neighbourWallInfo = allFaceInfo_[facei];
1025 
1026  // Evaluate all connected cells
1027 
1028  // Owner
1029  label celli = owner[facei];
1030  Type& currentWallInfo = allCellInfo_[celli];
1031 
1032  if (!currentWallInfo.equal(neighbourWallInfo, td_))
1033  {
1034  updateCell
1035  (
1036  celli,
1037  facei,
1038  neighbourWallInfo,
1039  propagationTol_,
1040  currentWallInfo
1041  );
1042  }
1043 
1044  // Neighbour.
1045  if (facei < nInternalFaces)
1046  {
1047  celli = neighbour[facei];
1048  Type& currentWallInfo2 = allCellInfo_[celli];
1049 
1050  if (!currentWallInfo2.equal(neighbourWallInfo, td_))
1051  {
1052  updateCell
1053  (
1054  celli,
1055  facei,
1056  neighbourWallInfo,
1057  propagationTol_,
1058  currentWallInfo2
1059  );
1060  }
1061  }
1062 
1063  // Reset status of face
1064  changedFace_[facei] = false;
1065  }
1066 
1067  // Handled all changed faces by now
1068  changedFaces_.clear();
1069 
1070  if (debug & 2)
1071  {
1072  Pout<< " Changed cells : " << changedCells_.size() << endl;
1073  }
1074 
1075  // Sum changedCells over all procs
1076  label totNChanged = changedCells_.size();
1077 
1078  reduce(totNChanged, sumOp<label>());
1079 
1080  return totNChanged;
1081 }
1082 
1083 
1084 template<class Type, class TrackingData>
1086 {
1087  // Propagate cell to face
1088 
1089  const cellList& cells = mesh_.cells();
1090 
1091  forAll(changedCells_, changedCelli)
1092  {
1093  label celli = changedCells_[changedCelli];
1094  if (!changedCell_[celli])
1095  {
1097  << "Cell " << celli << " not marked as having been changed"
1098  << abort(FatalError);
1099  }
1100 
1101  const Type& neighbourWallInfo = allCellInfo_[celli];
1102 
1103  // Evaluate all connected faces
1104 
1105  const labelList& faceLabels = cells[celli];
1106  forAll(faceLabels, faceLabelI)
1107  {
1108  label facei = faceLabels[faceLabelI];
1109  Type& currentWallInfo = allFaceInfo_[facei];
1110 
1111  if (!currentWallInfo.equal(neighbourWallInfo, td_))
1112  {
1113  updateFace
1114  (
1115  facei,
1116  celli,
1117  neighbourWallInfo,
1118  propagationTol_,
1119  currentWallInfo
1120  );
1121  }
1122  }
1123 
1124  // Reset status of cell
1125  changedCell_[celli] = false;
1126  }
1127 
1128  // Handled all changed cells by now
1129  changedCells_.clear();
1130 
1131 
1132  // Transfer across any explicitly provided internal connections
1133  handleExplicitConnections();
1134 
1135  if (hasCyclicPatches_)
1136  {
1137  // Transfer changed faces across cyclic halves
1138  handleCyclicPatches();
1139  }
1140 
1141  if (hasCyclicAMIPatches_)
1142  {
1143  handleAMICyclicPatches();
1144  }
1145 
1146  if (Pstream::parRun())
1147  {
1148  // Transfer changed faces from neighbouring processors.
1149  handleProcPatches();
1150  }
1151 
1152  if (debug & 2)
1153  {
1154  Pout<< " Changed faces : " << changedFaces_.size() << endl;
1155  }
1156 
1157  // Sum nChangedFaces over all procs
1158  label totNChanged = changedFaces_.size();
1159 
1160  reduce(totNChanged, sumOp<label>());
1161 
1162  return totNChanged;
1163 }
1164 
1165 
1166 // Iterate
1167 template<class Type, class TrackingData>
1169 {
1170  if (hasCyclicPatches_)
1171  {
1172  // Transfer changed faces across cyclic halves
1173  handleCyclicPatches();
1174  }
1175 
1176  if (hasCyclicAMIPatches_)
1177  {
1178  handleAMICyclicPatches();
1179  }
1180 
1181  if (Pstream::parRun())
1182  {
1183  // Transfer changed faces from neighbouring processors.
1184  handleProcPatches();
1185  }
1186 
1187  label iter = 0;
1188 
1189  while (iter < maxIter)
1190  {
1191  if (debug)
1192  {
1193  Info<< " Iteration " << iter << endl;
1194  }
1195 
1196  nEvals_ = 0;
1197 
1198  label nCells = faceToCell();
1199 
1200  if (debug)
1201  {
1202  Info<< " Total changed cells : " << nCells << endl;
1203  }
1204 
1205  if (nCells == 0)
1206  {
1207  break;
1208  }
1209 
1210  label nFaces = cellToFace();
1211 
1212  if (debug)
1213  {
1214  Info<< " Total changed faces : " << nFaces << nl
1215  << " Total evaluations : " << nEvals_ << nl
1216  << " Remaining unvisited cells: " << nUnvisitedCells_ << nl
1217  << " Remaining unvisited faces: " << nUnvisitedFaces_ << endl;
1218  }
1219 
1220  if (nFaces == 0)
1221  {
1222  break;
1223  }
1224 
1225  ++iter;
1226  }
1227 
1228  return iter;
1229 }
1230 
1231 
1232 // ************************************************************************* //
void handleAMICyclicPatches()
Merge data from across AMI cyclics.
Definition: FaceCellWave.C:540
void mergeFaceInfo(const polyPatch &patch, const label nFaces, const labelList &, const List< Type > &)
Merge received patch data into global data.
Definition: FaceCellWave.C:291
void handleExplicitConnections()
Merge data across explicitly provided local connections (usually.
Definition: FaceCellWave.C:737
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const labelList & processorPatches() const
Return list of processor patch labels.
const word & name() const
Return name.
Inter-processor communication reduction functions.
bool updateCell(const fvMesh &, const label thisCelli, const labelPair &neighbourPatchAndFacei, const FvWallInfoBase< WallInfo, Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.
Definition: FvWallInfoI.H:56
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
const UIndirectList< T > patchInternalList(const UList< T > &internalValues) const
Extract face cell data.
Definition: polyPatch.H:330
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void setFaceInfo(const labelList &changedFaces, const List< Type > &changedFacesInfo)
Set initial changed faces.
Definition: FaceCellWave.C:261
virtual label cellToFace()
Propagate from cell to face. Returns total number of faces.
error FatalError
void transform(const polyPatch &patch, const label nFaces, const labelList &faceLabels, const transformer &transform, List< Type > &faceInfo)
Transform across an interface. Implementation referred to Type.
Definition: FaceCellWave.C:357
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:297
static scalar propagationTol()
Access to tolerance.
Definition: FaceCellWave.H:234
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:76
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
virtual const transformer & transform() const
Return transformation between the coupled patches.
const cyclicPolyPatch & nbrPatch() const
const polyMesh & mesh() const
Access mesh.
Definition: FaceCellWave.H:327
void handleProcPatches()
Merge data from across processor boundaries.
Definition: FaceCellWave.C:373
scalar f1
Definition: createFields.H:15
Neighbour processor patch.
bool updateFace(const fvMesh &, const labelPair &thisPatchAndFacei, const label neighbourCelli, const FvWallInfoBase< WallInfo, Derived > &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
Definition: FvWallInfoI.H:79
const List< transformer > & AMITransforms() const
Return a reference to the AMI transforms.
virtual const cyclicAMIPolyPatch & nbrPatch() const
Return a reference to the neighbour patch.
A topoSetSource to select a faceSet from cells.
Definition: cellToFace.H:53
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:53
scalar y
label getUnsetFaces() const
Get number of unvisited faces.
Definition: FaceCellWave.C:997
A List obtained as a section of another List.
Definition: SubList.H:53
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
const cellShapeList & cells
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
bool hasPatch() const
Has cyclic patch?
Definition: FaceCellWave.C:246
const PtrList< AMIInterpolation > & AMIs() const
Return a reference to the AMI interpolators.
Cyclic plane patch.
A topoSetSource to select cells based on usage in faces.
Definition: faceToCell.H:49
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const transformer & transform() const
Return transformation between the coupled patches.
Cyclic patch for Arbitrary Mesh Interface (AMI)
int neighbProcNo() const
Return neighbour processor number.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:54
static const char nl
Definition: Ostream.H:260
bool updateFace(const label facei, const label neighbourCelli, const Type &neighbourInfo, const scalar tol, Type &faceInfo)
Updates faceInfo with information from neighbour. Updates all.
Definition: FaceCellWave.C:104
bool updateCell(const label celli, const label neighbourFacei, const Type &neighbourInfo, const scalar tol, Type &cellInfo)
Updates cellInfo with information from neighbour. Updates all.
Definition: FaceCellWave.C:55
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
FaceCellWave(const polyMesh &, UList< Type > &allFaceInfo, UList< Type > &allCellInfo, TrackingData &td=defaultTrackingData_)
Definition: FaceCellWave.C:813
virtual const transformer & transform() const
Return null transform between processor patches.
label patchi
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:325
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:315
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual bool owner() const
Does this side own the patch?
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
messageStream Info
void checkCyclic(const polyPatch &pPatch) const
Debugging: check info on both sides of cyclic.
Definition: FaceCellWave.C:200
void handleCyclicPatches()
Merge data from across cyclics.
Definition: FaceCellWave.C:475
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
treeBoundBox combine(const treeBoundBox &a, const treeBoundBox &b)
Definition: patchToPatch.C:77
virtual label faceToCell()
Propagate from face to cell. Returns total number of cells.
label getUnsetCells() const
Get number of unvisited cells, i.e. cells that were not (yet)
Definition: FaceCellWave.C:990