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