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