cyclicACMIPolyPatch.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 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 "cyclicACMIPolyPatch.H"
30 #include "SubField.H"
31 #include "Time.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(cyclicACMIPolyPatch, 0);
39 
40  addToRunTimeSelectionTable(polyPatch, cyclicACMIPolyPatch, word);
41  addToRunTimeSelectionTable(polyPatch, cyclicACMIPolyPatch, dictionary);
42 }
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
47 {
48  const polyMesh& mesh = boundaryMesh().mesh();
49 
50  bool updated = false;
51 
52  if (!owner())
53  {
54  return updated;
55  }
56 
57  // Check if underlying AMI up to date
58  if (!mesh.upToDatePoints(AMITime_))
59  {
60  // This should not happen normally since resetAMI is triggered
61  // by any point motion.
62  FatalErrorInFunction << "Problem : AMI is up to event:"
63  << AMITime_.eventNo()
64  << " mesh points are up to time " << mesh.pointsInstance()
65  << " patch:" << this->name()
66  << exit(FatalError);
67  }
68 
69  // Check if scaling enabled (and necessary)
70  if
71  (
72  srcScalePtr_.valid()
73  && (updated || prevTimeIndex_ != mesh.time().timeIndex())
74  )
75  {
76  if (debug)
77  {
78  Pout<< "cyclicACMIPolyPatch::updateAreas() :"
79  << " patch:" << this->name()
80  << " neighbPatch:" << this->neighbPatch().name()
81  << " AMITime_:" << AMITime_.eventNo()
82  << " uptodate:" << mesh.upToDatePoints(AMITime_)
83  << " mesh.time().timeIndex():" << mesh.time().timeIndex()
84  << " prevTimeIndex_:" << prevTimeIndex_
85  << endl;
86  }
87 
88  if (createAMIFaces_)
89  {
91  << "Topology changes and scaling currently not supported."
92  << " Patch " << this->name() << endl;
93  }
94 
95  const scalar t = mesh.time().timeOutputValue();
96 
97  // Note: ideally preserve src/tgtMask before clipping to tolerance ...
98  srcScaledMask_ =
99  min
100  (
101  scalar(1) - tolerance_,
102  max(tolerance_, srcScalePtr_->value(t)*srcMask_)
103  );
104 
105 
106  if (!tgtScalePtr_.valid())
107  {
108  tgtScalePtr_= srcScalePtr_.clone(neighbPatch());
109  }
110 
111  tgtScaledMask_ =
112  min
113  (
114  scalar(1) - tolerance_,
115  max(tolerance_, tgtScalePtr_->value(t)*tgtMask_)
116  );
117 
118  if (debug)
119  {
120  Pout<< "cyclicACMIPolyPatch::updateAreas : scaling masks"
121  << " for " << name() << " mask " << gAverage(srcScaledMask_)
122  << " and " << nonOverlapPatch().name()
123  << " mask " << gAverage(srcScaledMask_) << endl;
124  }
125 
126  // Calculate areas from the masks
127  cyclicACMIPolyPatch& cpp = const_cast<cyclicACMIPolyPatch&>(*this);
128  const cyclicACMIPolyPatch& nbrCpp = neighbPatch();
129 
130  cpp.scalePatchFaceAreas(*this, srcScaledMask_, thisSf_, thisNoSf_);
131  cpp.scalePatchFaceAreas(nbrCpp, tgtScaledMask_, nbrSf_, nbrNoSf_);
132 
133  prevTimeIndex_ = mesh.time().timeIndex();
134  AMITime_.setUpToDate();
135  updated = true;
136  }
137 
138  return updated;
139 }
140 
141 
143 {
144  // Is io up to date with
145  // - underlying AMI
146  // - scaling
147  return io.upToDate(AMITime_);
148 }
149 
150 
152 {
153  io.setUpToDate();
154 }
155 
156 
158 (
159  const word& name,
160  const scalarField& weightSum
161 ) const
162 {
163  label nUncovered = 0;
164  label nCovered = 0;
165  for (const scalar sum : weightSum)
166  {
167  if (sum < tolerance_)
168  {
169  ++nUncovered;
170  }
171  else if (sum > scalar(1) - tolerance_)
172  {
173  ++nCovered;
174  }
175  }
176  reduce(nUncovered, sumOp<label>());
177  reduce(nCovered, sumOp<label>());
178  label nTotal = returnReduce(weightSum.size(), sumOp<label>());
179 
180  Info<< "ACMI: Patch " << name << " uncovered/blended/covered = "
181  << nUncovered << ", " << nTotal-nUncovered-nCovered
182  << ", " << nCovered << endl;
183 }
184 
185 
187 (
188  const cyclicACMIPolyPatch& acmipp,
189  const scalarField& mask, // srcMask_
190  const vectorList& faceArea, // this->faceAreas();
191  const vectorList& noFaceArea // nonOverlapPatch.faceAreas()
192 )
193 {
194  // Primitive patch face areas have been cleared/reset based on the raw
195  // points - need to reset to avoid double-accounting of face areas
196 
197  const scalar maxTol = scalar(1) - tolerance_;
198 
199  const polyPatch& nonOverlapPatch = acmipp.nonOverlapPatch();
200  vectorField::subField noSf = nonOverlapPatch.faceAreas();
201 
202  DebugPout
203  << "rescaling non-overlap patch areas for: "
204  << nonOverlapPatch.name() << endl;
205 
206  if (mask.size() != noSf.size())
207  {
209  << "Inconsistent sizes for patch: " << acmipp.name()
210  << " - not manipulating patches" << nl
211  << " - size: " << size() << nl
212  << " - non-overlap patch size: " << noSf.size() << nl
213  << " - mask size: " << mask.size() << nl
214  << "This is OK for decomposition but"
215  << " should be considered fatal at run-time" << endl;
216 
217  return;
218  }
219 
220  forAll(noSf, facei)
221  {
222  const scalar w = min(maxTol, max(tolerance_, mask[facei]));
223  noSf[facei] = noFaceArea[facei]*(scalar(1) - w);
224  }
225 
226  if (!createAMIFaces_)
227  {
228  // Note: for topological update (createAMIFaces_ = true)
229  // AMI coupled patch face areas are updated as part of the topological
230  // updates, e.g. by the calls to cyclicAMIPolyPatch's setTopology and
231  // initMovePoints
232  DebugPout
233  << "scaling coupled patch areas for: " << acmipp.name() << endl;
234 
235  // Scale the coupled patch face areas
236  vectorField::subField Sf = acmipp.faceAreas();
237 
238  forAll(Sf, facei)
239  {
240  Sf[facei] = faceArea[facei]*max(tolerance_, mask[facei]);
241  }
242 
243  // Re-normalise the weights since the effect of overlap is already
244  // accounted for in the area
245  auto& weights = const_cast<scalarListList&>(acmipp.weights());
246  auto& weightsSum = const_cast<scalarField&>(acmipp.weightsSum());
247  forAll(weights, i)
248  {
249  scalarList& wghts = weights[i];
250  if (wghts.size())
251  {
252  scalar& sum = weightsSum[i];
253 
254  forAll(wghts, j)
255  {
256  wghts[j] /= sum;
257  }
258  sum = 1.0;
259  }
260  }
261  }
262 }
263 
264 
266 {
267  resetAMI(boundaryMesh().mesh().points());
268 }
269 
270 
272 {
273  if (!owner())
274  {
275  return;
276  }
277 
278  const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
279 
280  DebugPout
281  << "cyclicACMIPolyPatch::resetAMI : recalculating weights"
282  << " for " << name() << " and " << nonOverlapPatch.name()
283  << endl;
284 
285  const polyMesh& mesh = boundaryMesh().mesh();
286 
287  if (!createAMIFaces_ && mesh.hasCellCentres())
288  {
289  DebugPout
290  << "cyclicACMIPolyPatch::resetAMI : clearing cellCentres"
291  << " for " << name() << " and " << nonOverlapPatch.name() << nl
292  << "The mesh already has cellCentres calculated when"
293  << " resetting ACMI " << name() << "." << nl
294  << "This is a problem since ACMI adapts the face areas"
295  << " (to close cells) so this has" << nl
296  << "to be done before cell centre calculation." << nl
297  << "This can happen if e.g. the cyclicACMI is after"
298  << " any processor patches in the boundary." << endl;
299 
301  }
302 
303  // At this point we want face geometry but not cell geometry since we want
304  // correct the face area on duplicate baffles before calculating the cell
305  // centres and volumes.
306  if (!mesh.hasFaceAreas())
307  {
309  << "primitiveMesh must already have face geometry"
310  << abort(FatalError);
311  }
312 
313  // Calculate the AMI using partial face-area-weighted. This leaves
314  // the weights as fractions of local areas (sum(weights) = 1 means
315  // face is fully covered)
317 
318  const AMIPatchToPatchInterpolation& AMI = this->AMI();
319 
320  // Output some statistics
321  reportCoverage("source", AMI.srcWeightsSum());
322  reportCoverage("target", AMI.tgtWeightsSum());
323 
324  // Set the mask fields
325  // Note:
326  // - assumes that the non-overlap patches are decomposed using the same
327  // decomposition as the coupled patches (per side)
328  srcMask_ = min(scalar(1), max(scalar(0), AMI.srcWeightsSum()));
329  tgtMask_ = min(scalar(1), max(scalar(0), AMI.tgtWeightsSum()));
330 
331  if (debug)
332  {
333  Pout<< "resetAMI" << endl;
334  {
335  const cyclicACMIPolyPatch& patch = *this;
336  Pout<< "patch:" << patch.name() << " size:" << patch.size()
337  << " non-overlap patch: " << patch.nonOverlapPatch().name()
338  << " size:" << patch.nonOverlapPatch().size()
339  << " mask size:" << patch.srcMask().size() << endl;
340  }
341  {
342  const cyclicACMIPolyPatch& patch = this->neighbPatch();
343  Pout<< "patch:" << patch.name() << " size:" << patch.size()
344  << " non-overlap patch: " << patch.nonOverlapPatch().name()
345  << " size:" << patch.nonOverlapPatch().size()
346  << " mask size:" << patch.neighbPatch().tgtMask().size()
347  << endl;
348  }
349  }
350 }
351 
352 
354 {
355  if (!owner() || !canResetAMI())
356  {
357  return;
358  }
359 
360  const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
361  const cyclicACMIPolyPatch& nbrPatch = this->neighbPatch();
362  const polyPatch& nbrNonOverlapPatch = nbrPatch.nonOverlapPatch();
363 
364  if (srcScalePtr_.valid())
365  {
366  // Save overlap geometry for later scaling
367  thisSf_ = this->faceAreas();
368  thisNoSf_ = nonOverlapPatch.faceAreas();
369  nbrSf_ = nbrPatch.faceAreas();
370  nbrNoSf_ = nbrNonOverlapPatch.faceAreas();
371  }
372 
373  // In-place scale the patch areas
374  scalePatchFaceAreas
375  (
376  *this,
377  srcMask_, // unscaled mask
378  this->faceAreas(),
379  nonOverlapPatch.faceAreas()
380  );
381  scalePatchFaceAreas
382  (
383  nbrPatch,
384  tgtMask_, // unscaled mask
385  nbrPatch.faceAreas(),
386  nbrNonOverlapPatch.faceAreas()
387  );
388 
389  // Mark current AMI as up to date with points
390  boundaryMesh().mesh().setUpToDatePoints(AMITime_);
391 }
392 
393 
395 {
396  DebugPout << "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
397 
398  // Note: calculates transformation and triggers face centre calculation
400 
401  // Initialise the AMI early to make sure we adapt the face areas before the
402  // cell centre calculation gets triggered.
403  if (!createAMIFaces_ && canResetAMI())
404  {
405  resetAMI();
406  }
407 
408  scalePatchFaceAreas();
409 }
410 
411 
413 {
414  DebugPout << "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
415 
417 }
418 
419 
421 (
422  PstreamBuffers& pBufs,
423  const pointField& p
424 )
425 {
426  DebugPout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
427 
428  // Note: calculates transformation and triggers face centre calculation
430 
431  if (!createAMIFaces_ && canResetAMI())
432  {
433  resetAMI();
434  }
435 
436  scalePatchFaceAreas();
437 }
438 
439 
441 (
442  PstreamBuffers& pBufs,
443  const pointField& p
444 )
445 {
446  DebugPout << "cyclicACMIPolyPatch::movePoints : " << name() << endl;
447 
448  // When topology is changing, this will scale the duplicate AMI faces
450 }
451 
452 
454 {
455  DebugPout << "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
456 
458 }
459 
460 
462 {
463  DebugPout << "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
464 
466 }
467 
468 
470 {
471  DebugPout << "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
472 
474 }
475 
476 
478 {
479  if (srcScalePtr_.valid())
480  {
481  // Make sure areas are up-to-date
482  updateAreas();
483 
484  return srcScaledMask_;
485  }
486  else
487  {
488  return srcMask_;
489  }
490 }
491 
492 
494 {
495  if (tgtScalePtr_.valid())
496  {
497  // Make sure areas are up-to-date
498  updateAreas();
499 
500  return tgtScaledMask_;
501  }
502  else
503  {
504  return tgtMask_;
505  }
506 }
507 
508 
509 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
510 
512 (
513  const word& name,
514  const label size,
515  const label start,
516  const label index,
517  const polyBoundaryMesh& bm,
518  const word& patchType,
519  const transformType transform,
520  const word& defaultAMIMethod
521 )
522 :
524  (
525  name,
526  size,
527  start,
528  index,
529  bm,
530  patchType,
531  transform,
532  defaultAMIMethod
533  ),
534  nonOverlapPatchName_(word::null),
535  nonOverlapPatchID_(-1),
536  srcMask_(),
537  tgtMask_(),
538  AMITime_
539  (
540  IOobject
541  (
542  "AMITime",
543  boundaryMesh().mesh().pointsInstance(),
544  boundaryMesh().mesh(),
547  false
548  ),
549  dimensionedScalar("time", dimTime, -GREAT)
550  ),
551  prevTimeIndex_(-1)
552 {
553  AMIPtr_->setRequireMatch(false);
554 
555  // Non-overlapping patch might not be valid yet so cannot determine
556  // associated patchID
557 }
558 
559 
561 (
562  const word& name,
563  const dictionary& dict,
564  const label index,
565  const polyBoundaryMesh& bm,
566  const word& patchType,
567  const word& defaultAMIMethod
568 )
569 :
570  cyclicAMIPolyPatch(name, dict, index, bm, patchType, defaultAMIMethod),
571  nonOverlapPatchName_(dict.get<word>("nonOverlapPatch")),
572  nonOverlapPatchID_(-1),
573  srcMask_(),
574  tgtMask_(),
575  srcScalePtr_
576  (
577  dict.found("scale")
578  ? PatchFunction1<scalar>::New(*this, "scale", dict)
579  : nullptr
580  ),
581  AMITime_
582  (
583  IOobject
584  (
585  "AMITime",
586  boundaryMesh().mesh().pointsInstance(),
587  boundaryMesh().mesh(),
588  IOobject::NO_READ,
589  IOobject::NO_WRITE,
590  false
591  ),
592  dimensionedScalar("time", dimTime, -GREAT)
593  ),
594  prevTimeIndex_(-1)
595 {
596  AMIPtr_->setRequireMatch(false);
597 
598  if (nonOverlapPatchName_ == name)
599  {
601  << "Non-overlapping patch name " << nonOverlapPatchName_
602  << " cannot be the same as this patch " << name
603  << exit(FatalIOError);
604  }
605 
606  // Non-overlapping patch might not be valid yet so cannot determine
607  // associated patchID
608 }
609 
610 
612 (
613  const cyclicACMIPolyPatch& pp,
614  const polyBoundaryMesh& bm
615 )
616 :
617  cyclicAMIPolyPatch(pp, bm),
618  nonOverlapPatchName_(pp.nonOverlapPatchName_),
619  nonOverlapPatchID_(-1),
620  srcMask_(),
621  tgtMask_(),
622  srcScalePtr_
623  (
624  pp.srcScalePtr_.valid()
625  ? pp.srcScalePtr_.clone(*this)
626  : nullptr
627  ),
628  AMITime_
629  (
630  IOobject
631  (
632  "AMITime",
633  boundaryMesh().mesh().pointsInstance(),
634  boundaryMesh().mesh(),
637  false
638  ),
639  dimensionedScalar("time", dimTime, -GREAT)
640  ),
641  prevTimeIndex_(-1)
642 {
643  AMIPtr_->setRequireMatch(false);
644 
645  // Non-overlapping patch might not be valid yet so cannot determine
646  // associated patchID
647 }
648 
649 
651 (
652  const cyclicACMIPolyPatch& pp,
653  const polyBoundaryMesh& bm,
654  const label index,
655  const label newSize,
656  const label newStart,
657  const word& nbrPatchName,
658  const word& nonOverlapPatchName
659 )
660 :
661  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
662  nonOverlapPatchName_(nonOverlapPatchName),
663  nonOverlapPatchID_(-1),
664  srcMask_(),
665  tgtMask_(),
666  srcScalePtr_
667  (
668  pp.srcScalePtr_.valid()
669  ? pp.srcScalePtr_.clone(*this)
670  : nullptr
671  ),
672  AMITime_
673  (
674  IOobject
675  (
676  "AMITime",
677  boundaryMesh().mesh().pointsInstance(),
678  boundaryMesh().mesh(),
681  false
682  ),
683  dimensionedScalar("time", dimTime, -GREAT)
684  ),
685  prevTimeIndex_(-1)
686 {
687  AMIPtr_->setRequireMatch(false);
688 
689  if (nonOverlapPatchName_ == name())
690  {
692  << "Non-overlapping patch name " << nonOverlapPatchName_
693  << " cannot be the same as this patch " << name()
694  << exit(FatalError);
695  }
696 
697  // Non-overlapping patch might not be valid yet so cannot determine
698  // associated patchID
699 }
700 
701 
703 (
704  const cyclicACMIPolyPatch& pp,
705  const polyBoundaryMesh& bm,
706  const label index,
707  const labelUList& mapAddressing,
708  const label newStart
709 )
710 :
711  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
712  nonOverlapPatchName_(pp.nonOverlapPatchName_),
713  nonOverlapPatchID_(-1),
714  srcMask_(),
715  tgtMask_(),
716  srcScalePtr_
717  (
718  pp.srcScalePtr_.valid()
719  ? pp.srcScalePtr_.clone(*this)
720  : nullptr
721  ),
722  AMITime_
723  (
724  IOobject
725  (
726  "AMITime",
727  boundaryMesh().mesh().pointsInstance(),
728  boundaryMesh().mesh(),
731  false
732  ),
733  dimensionedScalar("time", dimTime, -GREAT)
734  ),
735  prevTimeIndex_(-1)
736 {
737  AMIPtr_->setRequireMatch(false);
738 }
739 
740 
741 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
742 
744 (
745  label& newFaces,
746  label& newProcFaces
747 ) const
748 {
749  const List<labelList>& addSourceFaces = AMI().srcAddress();
750  const scalarField& fMask = srcMask();
751 
752  // Add new faces as many weights for AMI
753  forAll (addSourceFaces, faceI)
754  {
755  if (fMask[faceI] > tolerance_)
756  {
757  const labelList& nbrFaceIs = addSourceFaces[faceI];
758 
759  forAll (nbrFaceIs, j)
760  {
761  label nbrFaceI = nbrFaceIs[j];
762 
763  if (nbrFaceI < neighbPatch().size())
764  {
765  // local faces
766  newFaces++;
767  }
768  else
769  {
770  // Proc faces
771  newProcFaces++;
772  }
773  }
774  }
775  }
776 }
777 
778 
780 {
781  const scalarField& fMask = srcMask();
782  const labelListList& srcFaces = AMI().srcAddress();
783  labelListList dOverFaces;
784 
785  dOverFaces.setSize(srcFaces.size());
786  forAll (dOverFaces, faceI)
787  {
788  if (fMask[faceI] > tolerance_)
789  {
790  dOverFaces[faceI].setSize(srcFaces[faceI].size());
791 
792  forAll (dOverFaces[faceI], subFaceI)
793  {
794  dOverFaces[faceI][subFaceI] = srcFaces[faceI][subFaceI];
795  }
796  }
797  }
798  return refPtr<labelListList>(new labelListList(dOverFaces));
799 }
800 
801 
803 {
804  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
805 
806  // Bit of checking now we know neighbour patch
807  if (!owner() && srcScalePtr_.valid())
808  {
810  << "Ignoring \"scale\" setting in slave patch " << name()
811  << endl;
812  srcScalePtr_.clear();
813  tgtScalePtr_.clear();
814  }
815 
816  return refCast<const cyclicACMIPolyPatch>(pp);
817 }
818 
819 
821 {
822  if (nonOverlapPatchID_ == -1)
823  {
824  nonOverlapPatchID_ =
825  this->boundaryMesh().findPatchID(nonOverlapPatchName_);
826 
827  if (nonOverlapPatchID_ == -1)
828  {
830  << "Illegal non-overlapping patch name " << nonOverlapPatchName_
831  << nl << "Valid patch names are "
832  << this->boundaryMesh().names()
833  << exit(FatalError);
834  }
835 
836  if (nonOverlapPatchID_ < index())
837  {
839  << "Boundary ordering error: " << type()
840  << " patch must be defined prior to its non-overlapping patch"
841  << nl
842  << type() << " patch: " << name() << ", ID:" << index() << nl
843  << "Non-overlap patch: " << nonOverlapPatchName_
844  << ", ID:" << nonOverlapPatchID_ << nl
845  << exit(FatalError);
846  }
847 
848  const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
849 
850  bool ok = true;
851 
852  if (size() == noPp.size())
853  {
854  const scalarField magSf(mag(faceAreas()));
855  const scalarField noMagSf(mag(noPp.faceAreas()));
856 
857  forAll(magSf, facei)
858  {
859  scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
860 
861  if (ratio - 1 > tolerance_)
862  {
863  ok = false;
864  break;
865  }
866  }
867  }
868  else
869  {
870  ok = false;
871  }
872 
873  if (!ok)
874  {
876  << "Inconsistent ACMI patches " << name() << " and "
877  << noPp.name() << ". Patches should have identical topology"
878  << exit(FatalError);
879  }
880  }
881 
882  return nonOverlapPatchID_;
883 }
884 
885 
887 (
888  PstreamBuffers& pBufs,
889  const primitivePatch& pp
890 ) const
891 {
893 }
894 
895 
897 (
898  PstreamBuffers& pBufs,
899  const primitivePatch& pp,
901  labelList& rotation
902 ) const
903 {
904  return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
905 }
906 
907 
909 {
911 
912  os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
913 
914  if (owner() && srcScalePtr_.valid())
915  {
916  srcScalePtr_->writeData(os);
917  }
918 }
919 
920 
921 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::cyclicAMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicAMIPolyPatch.C:1161
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::cyclicACMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicACMIPolyPatch.C:412
Foam::cyclicACMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicACMIPolyPatch.C:441
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::cyclicAMIPolyPatch::resetAMI
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
Definition: cyclicAMIPolyPatch.C:372
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:964
cyclicACMIPolyPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::regIOobject::upToDate
bool upToDate(const regIOobject &) const
Return true if up-to-date with respect to given object.
Definition: regIOobject.C:337
Foam::cyclicAMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Definition: cyclicAMIPolyPatch.C:1153
Foam::AMIInterpolation::tgtWeightsSum
const scalarField & tgtWeightsSum() const
Definition: AMIInterpolationI.H:219
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::cyclicACMIPolyPatch::upToDate
bool upToDate(const regIOobject &) const
Return true if given object is up to date with *this.
Definition: cyclicACMIPolyPatch.C:142
Foam::cyclicAMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicAMIPolyPatch.C:485
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Foam::primitiveMesh::clearCellGeom
void clearCellGeom()
Clear cell-based geometry only.
Definition: primitiveMeshClear.C:143
Foam::cyclicAMIPolyPatch::weights
const scalarListList & weights() const
Helper function to return the weights.
Definition: cyclicAMIPolyPatchI.H:68
Foam::cyclicACMIPolyPatch::setUpToDate
void setUpToDate(regIOobject &) const
Set object up to date with *this.
Definition: cyclicACMIPolyPatch.C:151
Foam::AMIInterpolation::srcWeightsSum
const scalarField & srcWeightsSum() const
Definition: AMIInterpolationI.H:153
Foam::cyclicAMIPolyPatch::tolerance_
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
Definition: cyclicAMIPolyPatch.H:376
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::regIOobject::eventNo
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:185
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::cyclicACMIPolyPatch::resetAMI
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
Definition: cyclicACMIPolyPatch.C:265
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::cyclicACMIPolyPatch::newInternalProcFaces
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal sub-faces and new proc faces.
Definition: cyclicACMIPolyPatch.C:744
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cyclicAMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicAMIPolyPatch.C:509
Foam::cyclicAMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicAMIPolyPatch.C:581
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:64
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::cyclicACMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicACMIPolyPatch.C:461
Foam::cyclicACMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicACMIPolyPatch.C:908
Foam::polyBoundaryMesh::mesh
const polyMesh & mesh() const noexcept
Return the mesh reference.
Definition: polyBoundaryMesh.H:152
Foam::cyclicAMIPolyPatch::weightsSum
const scalarField & weightsSum() const
Helper function to return the weights sum.
Definition: cyclicAMIPolyPatchI.H:79
Foam::cyclicAMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicAMIPolyPatch.C:502
Foam::Field< scalar >
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1225
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::cyclicACMIPolyPatch::srcMask
virtual const scalarField & srcMask() const
Return the mask/weighting for the source patch.
Definition: cyclicACMIPolyPatch.C:477
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::cyclicAMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicAMIPolyPatch.C:545
Foam::primitiveMesh::hasCellCentres
bool hasCellCentres() const noexcept
Definition: primitiveMeshI.H:189
Foam::primitiveMesh::hasFaceAreas
bool hasFaceAreas() const noexcept
Definition: primitiveMeshI.H:207
Foam::cyclicACMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicACMIPolyPatch.C:469
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::cyclicACMIPolyPatch::cyclicACMIPolyPatch
cyclicACMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const word &defaultAMIMethod=faceAreaWeightAMI::typeName)
Construct from (base coupled patch) components.
Definition: cyclicACMIPolyPatch.C:512
Foam::cyclicACMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicACMIPolyPatch.C:897
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::cyclicACMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
Definition: cyclicACMIPolyPatch.C:887
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicACMIPolyPatch::nonOverlapPatchID
virtual label nonOverlapPatchID() const
Non-overlapping patch ID.
Definition: cyclicACMIPolyPatch.C:820
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:60
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:393
Foam::regIOobject::setUpToDate
void setUpToDate()
Set as up-to-date.
Definition: regIOobject.C:409
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::cyclicACMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicACMIPolyPatch.C:421
Foam::cyclicACMIPolyPatch::scalePatchFaceAreas
virtual void scalePatchFaceAreas()
Scale patch face areas to maintain physical area.
Definition: cyclicACMIPolyPatch.C:353
Time.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:73
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::cyclicACMIPolyPatch::neighbPatch
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Definition: cyclicACMIPolyPatch.C:802
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::cyclicACMIPolyPatch::mapCollocatedFaces
virtual refPtr< labelListList > mapCollocatedFaces() const
Return collocated faces.
Definition: cyclicACMIPolyPatch.C:779
Foam::cyclicAMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicAMIPolyPatch.C:568
Foam::cyclicAMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicAMIPolyPatch.C:590
Foam::List< vector >
Foam::polyPatch::faceAreas
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
Foam::cyclicACMIPolyPatch::updateAreas
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything.
Definition: cyclicACMIPolyPatch.C:46
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::cyclicACMIPolyPatch::nonOverlapPatch
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
Definition: cyclicACMIPolyPatchI.H:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::cyclicAMIPolyPatch::createAMIFaces_
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
Definition: cyclicAMIPolyPatch.H:150
Foam::cyclicACMIPolyPatch::tgtMask
virtual const scalarField & tgtMask() const
Return the mask/weighting for the target patch.
Definition: cyclicACMIPolyPatch.C:493
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::cyclicACMIPolyPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
Definition: cyclicACMIPolyPatch.H:78
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:60
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::refPtr
A class for managing references or pointers (no reference counting)
Definition: PtrList.H:60
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::cyclicACMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicACMIPolyPatch.C:394
Foam::cyclicACMIPolyPatch::reportCoverage
void reportCoverage(const word &name, const scalarField &weightSum) const
Definition: cyclicACMIPolyPatch.C:158
Foam::cyclicACMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicACMIPolyPatch.C:453
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:68