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-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 "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  // Which processors this processor is connected to
530  const labelList& neighbourProcs = pData.topology().procNeighbours();
531 
532  // Reduce communication by only sending non-zero data,
533  // but with multiply-connected processor/processor
534  // (eg, processorCyclic) also need to send zero information
535  // to keep things synchronised
536 
537  // Reset buffers, initialize for registerSend() bookkeeping
538  pBufs_.clear();
539  pBufs_.initRegisterSend();
540 
541 
542  // Information to send
543  DynamicList<Type> sendFacesInfo;
544  DynamicList<label> sendFaces;
545 
546  for (const label patchi : procPatches)
547  {
548  const auto& procPatch =
549  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
550 
551  const label nbrProci = procPatch.neighbProcNo();
552 
553  // Resize buffers
554  sendFaces.resize_nocopy(procPatch.size());
555  sendFacesInfo.resize_nocopy(procPatch.size());
556 
557  // Determine which faces changed on current patch
558  const label nSendFaces = getChangedPatchFaces
559  (
560  procPatch,
561  0,
562  procPatch.size(),
563  sendFaces,
564  sendFacesInfo
565  );
566 
567  // Shrink
568  sendFaces.resize(nSendFaces);
569  sendFacesInfo.resize(nSendFaces);
570 
571  // Adapt wallInfo for leaving domain
572  leaveDomain
573  (
574  procPatch,
575  nSendFaces,
576  sendFaces,
577  sendFacesInfo
578  );
579 
580  // Send to neighbour
581  {
582  UOPstream toNbr(nbrProci, pBufs_);
583  toNbr << sendFaces << sendFacesInfo;
584 
585  // Record if send is required (data are non-zero)
586  pBufs_.registerSend(nbrProci, !sendFaces.empty());
587 
588  if (debug & 2)
589  {
590  Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
591  << " send:" << sendFaces.size() << " to proc:" << nbrProci
592  << endl;
593  }
594  }
595  }
596 
597  // Limit exchange to involved procs
598  // - automatically discards unnecessary (unregistered) sends
599  pBufs_.finishedNeighbourSends(neighbourProcs);
600 
601 
602  for (const label patchi : procPatches)
603  {
604  const auto& procPatch =
605  refCast<const processorPolyPatch>(mesh_.boundaryMesh()[patchi]);
606 
607  const label nbrProci = procPatch.neighbProcNo();
608 
609  if (!pBufs_.recvDataCount(nbrProci))
610  {
611  // Nothing to receive
612  continue;
613  }
614 
615 
616  labelList receiveFaces;
617  List<Type> receiveFacesInfo;
618  {
619  UIPstream is(nbrProci, pBufs_);
620  is >> receiveFaces >> receiveFacesInfo;
621  }
622 
623  const label nReceiveFaces = receiveFaces.size();
624 
625  if (debug & 2)
626  {
627  Pout<< " Processor patch " << patchi << ' ' << procPatch.name()
628  << " recv:" << receiveFaces.size() << " from proci:"
629  << nbrProci << endl;
630  }
631 
632  // Apply transform to received data for non-parallel planes
633  if (!procPatch.parallel())
634  {
635  transform
636  (
637  procPatch.forwardT(),
638  nReceiveFaces,
639  receiveFacesInfo
640  );
641  }
642 
643  // Adapt wallInfo for entering domain
644  enterDomain
645  (
646  procPatch,
647  nReceiveFaces,
648  receiveFaces,
649  receiveFacesInfo
650  );
651 
652  // Merge received info
653  mergeFaceInfo
654  (
655  procPatch,
656  nReceiveFaces,
657  receiveFaces,
658  receiveFacesInfo
659  );
660  }
661 }
662 
663 
664 template<class Type, class TrackingData>
666 {
667  // Transfer information across cyclic halves.
668 
669  for (const polyPatch& patch : mesh_.boundaryMesh())
670  {
671  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(patch);
672 
673  if (cpp)
674  {
675  const auto& cycPatch = *cpp;
676  const auto& nbrPatch = cycPatch.neighbPatch();
677 
678  // Allocate buffers
679  labelList receiveFaces(patch.size());
680  List<Type> receiveFacesInfo(patch.size());
681 
682  // Determine which faces changed
683  const label nReceiveFaces = getChangedPatchFaces
684  (
685  nbrPatch,
686  0,
687  nbrPatch.size(),
688  receiveFaces,
689  receiveFacesInfo
690  );
691 
692  // Adapt wallInfo for leaving domain
693  leaveDomain
694  (
695  nbrPatch,
696  nReceiveFaces,
697  receiveFaces,
698  receiveFacesInfo
699  );
700 
701  if (!cycPatch.parallel())
702  {
703  // Received data from other half
704  transform
705  (
706  cycPatch.forwardT(),
707  nReceiveFaces,
708  receiveFacesInfo
709  );
710  }
711 
712  if (debug & 2)
713  {
714  Pout<< " Cyclic patch "
715  << cycPatch.index() << ' ' << cycPatch.name()
716  << " Changed : " << nReceiveFaces
717  << endl;
718  }
719 
720  // Half2: Adapt wallInfo for entering domain
721  enterDomain
722  (
723  cycPatch,
724  nReceiveFaces,
725  receiveFaces,
726  receiveFacesInfo
727  );
728 
729  // Merge into global storage
730  mergeFaceInfo
731  (
732  cycPatch,
733  nReceiveFaces,
734  receiveFaces,
735  receiveFacesInfo
736  );
737 
738  if (debug)
739  {
740  checkCyclic(cycPatch);
741  }
742  }
743  }
744 }
745 
746 
747 template<class Type, class TrackingData>
749 {
750  // Transfer information across cyclicAMI halves.
751 
752  for (const polyPatch& patch : mesh_.boundaryMesh())
753  {
754  const cyclicAMIPolyPatch* cpp = isA<cyclicAMIPolyPatch>(patch);
755 
756  // Note:
757  // - can either do owner and neighbour separately or have owner
758  // do both
759  // - separately means that neighbour will receive updated owner
760  // properties which might be beneficial or involve extra work?
761 
762  if (cpp)
763  {
764  const auto& cycPatch = *cpp;
765  const auto& nbrPatch = cycPatch.neighbPatch();
766 
767  List<Type> receiveInfo;
768 
769  {
770  // Get nbrPatch data (so not just changed faces). Do not use
771  // a slice here since the leaveDomain might change the values
772  List<Type> sendInfo(nbrPatch.patchSlice(allFaceInfo_));
773 
774  if (!nbrPatch.parallel() || nbrPatch.separated())
775  {
776  // Adapt sendInfo for leaving domain
777  const vectorField::subField fc = nbrPatch.faceCentres();
778  forAll(sendInfo, i)
779  {
780  sendInfo[i].leaveDomain(mesh_, nbrPatch, i, fc[i], td_);
781  }
782  }
783 
784  // Transfer sendInfo to cycPatch
785  combine<Type, TrackingData> cmb(*this, cycPatch);
786 
787  if (cycPatch.applyLowWeightCorrection())
788  {
789  List<Type> defVals
790  (
791  cycPatch.patchInternalList(allCellInfo_)
792  );
793 
794  cycPatch.interpolate(sendInfo, cmb, receiveInfo, defVals);
795  }
796  else
797  {
798  cycPatch.interpolate(sendInfo, cmb, receiveInfo);
799  }
800  }
801 
802  // Apply transform to received data for non-parallel planes
803  if (!cycPatch.parallel())
804  {
805  transform
806  (
807  cycPatch.forwardT(),
808  receiveInfo.size(),
809  receiveInfo
810  );
811  }
812 
813  if (!cycPatch.parallel() || cycPatch.separated())
814  {
815  // Adapt receiveInfo for entering domain
816  const vectorField::subField fc = cycPatch.faceCentres();
817  forAll(receiveInfo, i)
818  {
819  receiveInfo[i].enterDomain(mesh_, cycPatch, i, fc[i], td_);
820  }
821  }
822 
823  // Merge into global storage
824  forAll(receiveInfo, i)
825  {
826  const label meshFacei = cycPatch.start()+i;
827 
828  const Type& newInfo = receiveInfo[i];
829 
830  Type& currInfo = allFaceInfo_[meshFacei];
831 
832  if (newInfo.valid(td_) && !currInfo.equal(newInfo, td_))
833  {
834  updateFace
835  (
836  meshFacei,
837  newInfo,
838  propagationTol_,
839  currInfo
840  );
841  }
842  }
843  }
844  }
845 }
846 
847 
848 template<class Type, class TrackingData>
850 {
851  changedBaffles_.clear();
852 
853  // Collect all/any changed information touching a baffle
854  for (const labelPair& baffle : explicitConnections_)
855  {
856  const label f0 = baffle.first();
857  const label f1 = baffle.second();
858 
859  if (changedFace_.test(f0))
860  {
861  // f0 changed. Update information on f1.
862  changedBaffles_.push_back(taggedInfoType(f1, allFaceInfo_[f0]));
863  }
864 
865  if (changedFace_.test(f1))
866  {
867  // f1 changed. Update information on f0.
868  changedBaffles_.push_back(taggedInfoType(f0, allFaceInfo_[f1]));
869  }
870  }
871 
872 
873  // Update other side with changed information
874 
875  for (const taggedInfoType& updated : changedBaffles_)
876  {
877  const label tgtFace = updated.first;
878  const Type& newInfo = updated.second;
879 
880  Type& currInfo = allFaceInfo_[tgtFace];
881 
882  if (!currInfo.equal(newInfo, td_))
883  {
884  updateFace
885  (
886  tgtFace,
887  newInfo,
888  propagationTol_,
889  currInfo
890  );
891  }
892  }
893 
894  changedBaffles_.clear();
895 }
896 
898 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
899 
900 template<class Type, class TrackingData>
902 (
903  const polyMesh& mesh,
904  UList<Type>& allFaceInfo,
905  UList<Type>& allCellInfo,
906  TrackingData& td
907 )
908 :
910 
911  explicitConnections_(),
912  allFaceInfo_(allFaceInfo),
913  allCellInfo_(allCellInfo),
914  td_(td),
915  changedBaffles_(2*explicitConnections_.size()),
916  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
917  hasCyclicAMIPatches_
918  (
920  ),
921  nEvals_(0)
922 {
923  if
924  (
925  allFaceInfo.size() != mesh_.nFaces()
926  || allCellInfo.size() != mesh_.nCells()
927  )
928  {
930  << "face and cell storage not the size of mesh faces, cells:" << nl
931  << " allFaceInfo :" << allFaceInfo.size() << nl
932  << " mesh_.nFaces():" << mesh_.nFaces() << nl
933  << " allCellInfo :" << allCellInfo.size() << nl
934  << " mesh_.nCells():" << mesh_.nCells() << endl
935  << exit(FatalError);
936  }
937 }
938 
939 
940 template<class Type, class TrackingData>
942 (
943  const polyMesh& mesh,
944  const labelUList& changedFaces,
945  const List<Type>& changedFacesInfo,
946  UList<Type>& allFaceInfo,
947  UList<Type>& allCellInfo,
948  const label maxIter,
949  TrackingData& td
950 )
951 :
953 
954  explicitConnections_(),
955  allFaceInfo_(allFaceInfo),
956  allCellInfo_(allCellInfo),
957  td_(td),
958  changedBaffles_(2*explicitConnections_.size()),
959  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
960  hasCyclicAMIPatches_
961  (
963  ),
964  nEvals_(0)
965 {
966  if
967  (
968  allFaceInfo.size() != mesh_.nFaces()
969  || allCellInfo.size() != mesh_.nCells()
970  )
971  {
973  << "face and cell storage not the size of mesh faces, cells:" << nl
974  << " allFaceInfo :" << allFaceInfo.size() << nl
975  << " mesh_.nFaces():" << mesh_.nFaces() << nl
976  << " allCellInfo :" << allCellInfo.size() << nl
977  << " mesh_.nCells():" << mesh_.nCells() << endl
978  << exit(FatalError);
979  }
980 
981  // Copy initial changed faces data
982  setFaceInfo(changedFaces, changedFacesInfo);
983 
984  // Iterate until nothing changes
985  const label iter = iterate(maxIter);
986 
987  if ((maxIter > 0) && (iter >= maxIter))
988  {
990  << "Maximum number of iterations reached. Increase maxIter." << nl
991  << " maxIter:" << maxIter << nl
992  << " nChangedCells:" << nChangedCells() << nl
993  << " nChangedFaces:" << nChangedFaces() << endl
994  << exit(FatalError);
995  }
996 }
997 
998 
999 template<class Type, class TrackingData>
1001 (
1002  const polyMesh& mesh,
1003  const labelPairList& explicitConnections,
1004  const bool handleCyclicAMI,
1005  const labelUList& changedFaces,
1006  const List<Type>& changedFacesInfo,
1007  UList<Type>& allFaceInfo,
1008  UList<Type>& allCellInfo,
1009  const label maxIter,
1010  TrackingData& td
1011 )
1012 :
1014 
1015  explicitConnections_(explicitConnections),
1016  allFaceInfo_(allFaceInfo),
1017  allCellInfo_(allCellInfo),
1018  td_(td),
1019  changedBaffles_(2*explicitConnections_.size()),
1020  hasCyclicPatches_(hasPatch<cyclicPolyPatch>()),
1021  hasCyclicAMIPatches_
1022  (
1023  handleCyclicAMI
1024  && returnReduceOr(hasPatch<cyclicAMIPolyPatch>())
1025  ),
1026  nEvals_(0)
1027 {
1028  if
1029  (
1030  allFaceInfo.size() != mesh_.nFaces()
1031  || allCellInfo.size() != mesh_.nCells()
1032  )
1033  {
1035  << "face and cell storage not the size of mesh faces, cells:" << nl
1036  << " allFaceInfo :" << allFaceInfo.size() << nl
1037  << " mesh_.nFaces():" << mesh_.nFaces() << nl
1038  << " allCellInfo :" << allCellInfo.size() << nl
1039  << " mesh_.nCells():" << mesh_.nCells() << endl
1040  << exit(FatalError);
1041  }
1042 
1043  // Copy initial changed faces data
1044  setFaceInfo(changedFaces, changedFacesInfo);
1045 
1046  // Iterate until nothing changes
1047  const label iter = iterate(maxIter);
1048 
1049  if ((maxIter > 0) && (iter >= maxIter))
1050  {
1052  << "Maximum number of iterations reached. Increase maxIter." << nl
1053  << " maxIter:" << maxIter << nl
1054  << " nChangedCells:" << nChangedCells() << nl
1055  << " nChangedFaces:" << nChangedFaces() << endl
1056  << exit(FatalError);
1057  }
1058 }
1060 
1061 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1062 
1063 template<class Type, class TrackingData>
1065 {
1066  // Propagate face to cell
1067 
1068  const labelList& owner = mesh_.faceOwner();
1069  const labelList& neighbour = mesh_.faceNeighbour();
1070  const label nInternalFaces = mesh_.nInternalFaces();
1071 
1072  for (const label facei : changedFaces_)
1073  {
1074  if (!changedFace_.test(facei))
1075  {
1077  << "Face " << facei
1078  << " not marked as having been changed"
1079  << abort(FatalError);
1080  }
1081 
1082  const Type& newInfo = allFaceInfo_[facei];
1083 
1084  // Evaluate all connected cells
1085 
1086  // Owner
1087  {
1088  const label celli = owner[facei];
1089  Type& currInfo = allCellInfo_[celli];
1090 
1091  if (!currInfo.equal(newInfo, td_))
1092  {
1093  updateCell
1094  (
1095  celli,
1096  facei,
1097  newInfo,
1098  propagationTol_,
1099  currInfo
1100  );
1101  }
1102  }
1103 
1104  // Neighbour.
1105  if (facei < nInternalFaces)
1106  {
1107  const label celli = neighbour[facei];
1108  Type& currInfo = allCellInfo_[celli];
1109 
1110  if (!currInfo.equal(newInfo, td_))
1111  {
1112  updateCell
1113  (
1114  celli,
1115  facei,
1116  newInfo,
1117  propagationTol_,
1118  currInfo
1119  );
1120  }
1121  }
1122 
1123  // Reset status of face
1124  changedFace_.unset(facei);
1125  }
1126 
1127  // Handled all changed faces by now
1128  changedFaces_.clear();
1129 
1130  if (debug & 2)
1131  {
1132  Pout<< " Changed cells : " << nChangedCells() << endl;
1133  }
1134 
1135  // Number of changedCells over all procs
1136  return returnReduce(nChangedCells(), sumOp<label>());
1137 }
1138 
1139 
1140 template<class Type, class TrackingData>
1142 {
1143  // Propagate cell to face
1144 
1145  const cellList& cells = mesh_.cells();
1146 
1147  for (const label celli : changedCells_)
1148  {
1149  if (!changedCell_.test(celli))
1150  {
1152  << "Cell " << celli << " not marked as having been changed"
1153  << abort(FatalError);
1154  }
1155 
1156  const Type& newInfo = allCellInfo_[celli];
1157 
1158  // Evaluate all connected faces
1159 
1160  const labelList& faceLabels = cells[celli];
1161  for (const label facei : faceLabels)
1162  {
1163  Type& currInfo = allFaceInfo_[facei];
1164 
1165  if (!currInfo.equal(newInfo, td_))
1166  {
1167  updateFace
1168  (
1169  facei,
1170  celli,
1171  newInfo,
1172  propagationTol_,
1173  currInfo
1174  );
1175  }
1176  }
1177 
1178  // Reset status of cell
1179  changedCell_.unset(celli);
1180  }
1181 
1182  // Handled all changed cells by now
1183  changedCells_.clear();
1184 
1185 
1186  // Transfer across any explicitly provided internal connections
1187  handleExplicitConnections();
1188 
1189  if (hasCyclicPatches_)
1190  {
1191  handleCyclicPatches();
1192  }
1193 
1194  if (hasCyclicAMIPatches_)
1195  {
1196  handleAMICyclicPatches();
1197  }
1198 
1199  if (Pstream::parRun())
1200  {
1201  handleProcPatches();
1202  }
1203 
1204  if (debug & 2)
1205  {
1206  Pout<< " Changed faces : " << nChangedFaces() << endl;
1207  }
1208 
1209 
1210  // Number of changedFaces over all procs
1211  return returnReduce(nChangedFaces(), sumOp<label>());
1213 
1214 
1215 // Iterate
1216 template<class Type, class TrackingData>
1217 Foam::label Foam::FaceCellWave<Type, TrackingData>::iterate(const label maxIter)
1218 {
1219  if (maxIter < 0)
1220  {
1221  return 0;
1222  }
1223 
1224  if (hasCyclicPatches_)
1225  {
1226  handleCyclicPatches();
1227  }
1228 
1229  if (hasCyclicAMIPatches_)
1230  {
1231  handleAMICyclicPatches();
1232  }
1233 
1234  if (Pstream::parRun())
1235  {
1236  handleProcPatches();
1237  }
1238 
1239  label iter = 0;
1240 
1241  for (/*nil*/; iter < maxIter; ++iter)
1242  {
1243  if (debug)
1244  {
1245  Info<< " Iteration " << iter << endl;
1246  }
1247 
1248  nEvals_ = 0;
1249  const label nCells = faceToCell();
1250  const label nFaces = nCells ? cellToFace() : 0;
1251 
1252  if (debug)
1253  {
1254  Info<< " Total evaluations : "
1255  << nEvals_ << nl
1256  << " Changed cells / faces : "
1257  << nCells << " / " << nFaces << nl
1258  << " Pending cells / faces : "
1259  << nUnvisitedCells_ << " / " << nUnvisitedFaces_ << nl;
1260  }
1261 
1262  if (!nCells || !nFaces)
1263  {
1264  break;
1265  }
1266  }
1267 
1268  return iter;
1269 }
1270 
1271 
1272 // ************************************************************************* //
void handleAMICyclicPatches()
Merge data from across AMI cyclics.
Definition: FaceCellWave.C:743
void handleExplicitConnections()
Merge data across explicitly provided local connections.
Definition: FaceCellWave.C:844
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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 solver class.
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...
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
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:489
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:206
const TrackingData & data() const noexcept
Additional data to be passed into container.
Definition: FaceCellWave.H:506
Holds information regarding type of cell. Used in inside/outside determination in cellClassification...
Definition: cellInfo.H:60
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
labelList faceLabels(nFaceLabels)
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:147
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
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:168
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:51
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
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual bool owner() const
Does this side own the patch?
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:660
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
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
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.