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