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