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