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 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  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  for (scalar& w : wghts)
255  {
256  w /= sum;
257  }
258  sum = 1.0;
259  }
260  }
261  }
262 
263  const polyMesh& mesh = boundaryMesh().mesh();
264 
265  // Recompute the cell volumes adjacent to the patches since the cells with
266  // duplicate faces are only closed after the duplicate faces have been
267  // scaled.
268 
270  (
271  mesh,
272  mesh.faceCentres(),
273  mesh.faceAreas(), // already scaled
274  uniqueSort(acmipp.faceCells()), // unique cells only
276  const_cast<vectorField&>(mesh.cellCentres()),
277  const_cast<scalarField&>(mesh.cellVolumes())
278  );
279 }
280 
283 {
284  resetAMI(boundaryMesh().mesh().points());
285 }
286 
287 
289 {
290  if (!owner())
291  {
292  return;
293  }
294 
295  const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
296 
297  DebugPout
298  << "cyclicACMIPolyPatch::resetAMI : recalculating weights"
299  << " for " << name() << " and " << nonOverlapPatch.name()
300  << endl;
301 
302  const polyMesh& mesh = boundaryMesh().mesh();
303 
304  // At this point we want face geometry but not cell geometry since we want
305  // correct the face area on duplicate baffles before calculating the cell
306  // centres and volumes.
307  if (!mesh.hasFaceAreas())
308  {
310  << "primitiveMesh must already have face geometry"
311  << abort(FatalError);
312  }
313 
314  // Calculate the AMI using partial face-area-weighted. This leaves
315  // the weights as fractions of local areas (sum(weights) = 1 means
316  // face is fully covered)
318 
319  const AMIPatchToPatchInterpolation& AMI = this->AMI();
320 
321  // Output some statistics
322  reportCoverage("source", AMI.srcWeightsSum());
323  reportCoverage("target", AMI.tgtWeightsSum());
324 
325  // Set the mask fields
326  // Note:
327  // - assumes that the non-overlap patches are decomposed using the same
328  // decomposition as the coupled patches (per side)
329  srcMask_ = clamp(AMI.srcWeightsSum(), zero_one{});
330  tgtMask_ = clamp(AMI.tgtWeightsSum(), zero_one{});
331 
332  if (debug)
333  {
334  Pout<< "resetAMI" << endl;
335  {
336  const cyclicACMIPolyPatch& patch = *this;
337  Pout<< "patch:" << patch.name() << " size:" << patch.size()
338  << " non-overlap patch: " << patch.nonOverlapPatch().name()
339  << " size:" << patch.nonOverlapPatch().size()
340  << endl;
341  }
342  {
343  const cyclicACMIPolyPatch& patch = this->neighbPatch();
344  Pout<< "patch:" << patch.name() << " size:" << patch.size()
345  << " non-overlap patch: " << patch.nonOverlapPatch().name()
346  << " size:" << patch.nonOverlapPatch().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_)
365  {
366  // Use primitivePatch::faceAreas() to avoid need for additional copy?
367 
368  // Save overlap geometry for later scaling
369  thisSf_ = this->faceAreas();
370  thisNoSf_ = nonOverlapPatch.faceAreas();
371  nbrSf_ = nbrPatch.faceAreas();
372  nbrNoSf_ = nbrNonOverlapPatch.faceAreas();
373  }
374 
375  // In-place scale the polyPatch face areas
376 
377  // Note
378  // - using primitivePatch face areas since these are based on the raw
379  // point locations (not affected by ACMI scaling)
380  scalePatchFaceAreas
381  (
382  *this,
383  srcMask_, // unscaled mask
385  nonOverlapPatch.primitivePatch::faceAreas()
386  );
387  scalePatchFaceAreas
388  (
389  nbrPatch,
390  tgtMask_, // unscaled mask
391  nbrPatch.primitivePatch::faceAreas(),
392  nbrNonOverlapPatch.primitivePatch::faceAreas()
393  );
394 
395  // Mark current AMI as up to date with points
396  boundaryMesh().mesh().setUpToDatePoints(AMITime_);
397 
398  if (debug)
399  {
400  const auto& acmipp = *this;
401  const auto& mesh = boundaryMesh().mesh();
402  const auto& vols = mesh.cellVolumes();
403 
404  Info<< "cyclicACMI PP: " << acmipp.name()
405  << " V:" << sum(scalarField(vols, acmipp.faceCells()))
406  << nl
407  << "cyclicACMI N-O PP: " << nonOverlapPatch.name()
408  << " V:" << sum(scalarField(vols, nonOverlapPatch.faceCells()))
409  << nl
410  << "cyclicACMI NBR PP: " << nbrPatch.name()
411  << " V:" << sum(scalarField(vols, nbrPatch.faceCells()))
412  << nl
413  << "cyclicACMI NBR N-O PP: " << nbrNonOverlapPatch.name()
414  << " V:" << sum(scalarField(vols, nbrNonOverlapPatch.faceCells()))
415  << nl
416  << "cyclicACMI PP+N-O AREA: "
417  << sum(faceAreas() + nonOverlapPatch.faceAreas()) << nl
418  << "cyclicACMI NBR PP+N-O AREA: "
419  << sum(nbrPatch.faceAreas() + nbrNonOverlapPatch.faceAreas())
420  << endl;
421  }
422 }
423 
424 
425 void Foam::cyclicACMIPolyPatch::initGeometry(PstreamBuffers& pBufs)
426 {
427  DebugPout << "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
428 
429  // Note: calculates transformation and triggers face centre calculation
431 
432  // Initialise the AMI early to make sure we adapt the face areas before the
433  // cell centre calculation gets triggered.
434  if (!createAMIFaces_ && canResetAMI())
435  {
436  resetAMI();
437  }
438 
439  scalePatchFaceAreas();
440 }
441 
442 
443 void Foam::cyclicACMIPolyPatch::calcGeometry(PstreamBuffers& pBufs)
444 {
445  DebugPout << "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
446 
448 }
449 
450 
452 (
453  PstreamBuffers& pBufs,
454  const pointField& p
455 )
456 {
457  DebugPout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
458 
459  // Note: calculates transformation and triggers face centre calculation
461 
462  if (!createAMIFaces_ && canResetAMI())
463  {
464  resetAMI();
465  }
466 
467  scalePatchFaceAreas();
468 }
469 
470 
472 (
473  PstreamBuffers& pBufs,
474  const pointField& p
475 )
476 {
477  DebugPout << "cyclicACMIPolyPatch::movePoints : " << name() << endl;
478 
479  // When topology is changing, this will scale the duplicate AMI faces
481 }
482 
483 
485 {
486  DebugPout << "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
487 
489 }
490 
491 
493 {
494  DebugPout << "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
495 
497 }
498 
499 
501 {
502  DebugPout << "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
503 
505 }
506 
507 
509 {
510  if (srcScalePtr_)
511  {
512  // Make sure areas are up-to-date
513  updateAreas();
514 
515  return srcScaledMask_;
516  }
517  else
518  {
519  return srcMask_;
520  }
521 }
522 
523 
525 {
526  if (tgtScalePtr_)
527  {
528  // Make sure areas are up-to-date
529  updateAreas();
530 
531  return tgtScaledMask_;
532  }
533  else
534  {
535  return tgtMask_;
536  }
537 }
538 
539 
540 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
541 
543 (
544  const word& name,
545  const label size,
546  const label start,
547  const label index,
548  const polyBoundaryMesh& bm,
549  const word& patchType,
550  const transformType transform,
551  const word& defaultAMIMethod
552 )
553 :
554  cyclicAMIPolyPatch
555  (
556  name,
557  size,
558  start,
559  index,
560  bm,
561  patchType,
562  transform,
563  defaultAMIMethod
564  ),
565  nonOverlapPatchName_(),
566  nonOverlapPatchID_(-1),
567  srcMask_(),
568  tgtMask_(),
569  AMITime_
570  (
571  IOobject
572  (
573  "AMITime",
574  boundaryMesh().mesh().pointsInstance(),
575  boundaryMesh().mesh(),
576  IOobject::NO_READ,
577  IOobject::NO_WRITE,
578  IOobject::NO_REGISTER
579  ),
580  dimensionedScalar("time", dimTime, -GREAT)
581  ),
582  prevTimeIndex_(-1)
583 {
584  AMIPtr_->setRequireMatch(false);
586  // Non-overlapping patch might not be valid yet so cannot determine
587  // associated patchID
588 }
589 
590 
592 (
593  const word& name,
594  const dictionary& dict,
595  const label index,
596  const polyBoundaryMesh& bm,
597  const word& patchType,
598  const word& defaultAMIMethod
599 )
600 :
601  cyclicAMIPolyPatch(name, dict, index, bm, patchType, defaultAMIMethod),
602  nonOverlapPatchName_(dict.get<word>("nonOverlapPatch")),
603  nonOverlapPatchID_(-1),
604  srcMask_(),
605  tgtMask_(),
606  srcScalePtr_(PatchFunction1<scalar>::NewIfPresent(*this, "scale", dict)),
607  AMITime_
608  (
609  IOobject
610  (
611  "AMITime",
612  boundaryMesh().mesh().pointsInstance(),
613  boundaryMesh().mesh(),
614  IOobject::NO_READ,
615  IOobject::NO_WRITE,
616  IOobject::NO_REGISTER
617  ),
618  dimensionedScalar("time", dimTime, -GREAT)
619  ),
620  prevTimeIndex_(-1)
621 {
622  AMIPtr_->setRequireMatch(false);
623 
624  if (nonOverlapPatchName_ == name)
625  {
627  << "Non-overlapping patch name " << nonOverlapPatchName_
628  << " cannot be the same as this patch " << name
629  << exit(FatalIOError);
630  }
632  // Non-overlapping patch might not be valid yet so cannot determine
633  // associated patchID
634 }
635 
636 
638 (
639  const cyclicACMIPolyPatch& pp,
640  const polyBoundaryMesh& bm
641 )
642 :
643  cyclicAMIPolyPatch(pp, bm),
644  nonOverlapPatchName_(pp.nonOverlapPatchName_),
645  nonOverlapPatchID_(-1),
646  srcMask_(),
647  tgtMask_(),
648  srcScalePtr_(pp.srcScalePtr_.clone(*this)),
649  AMITime_
650  (
651  IOobject
652  (
653  "AMITime",
654  boundaryMesh().mesh().pointsInstance(),
655  boundaryMesh().mesh(),
656  IOobject::NO_READ,
657  IOobject::NO_WRITE,
658  IOobject::NO_REGISTER
659  ),
660  dimensionedScalar("time", dimTime, -GREAT)
661  ),
662  prevTimeIndex_(-1)
663 {
664  AMIPtr_->setRequireMatch(false);
666  // Non-overlapping patch might not be valid yet so cannot determine
667  // associated patchID
668 }
669 
670 
672 (
673  const cyclicACMIPolyPatch& pp,
674  const polyBoundaryMesh& bm,
675  const label index,
676  const label newSize,
677  const label newStart,
678  const word& nbrPatchName,
679  const word& nonOverlapPatchName
680 )
681 :
682  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
683  nonOverlapPatchName_(nonOverlapPatchName),
684  nonOverlapPatchID_(-1),
685  srcMask_(),
686  tgtMask_(),
687  srcScalePtr_(pp.srcScalePtr_.clone(*this)),
688  AMITime_
689  (
690  IOobject
691  (
692  "AMITime",
693  boundaryMesh().mesh().pointsInstance(),
694  boundaryMesh().mesh(),
695  IOobject::NO_READ,
696  IOobject::NO_WRITE,
697  IOobject::NO_REGISTER
698  ),
699  dimensionedScalar("time", dimTime, -GREAT)
700  ),
701  prevTimeIndex_(-1)
702 {
703  AMIPtr_->setRequireMatch(false);
704 
705  if (nonOverlapPatchName_ == name())
706  {
708  << "Non-overlapping patch name " << nonOverlapPatchName_
709  << " cannot be the same as this patch " << name()
710  << exit(FatalError);
711  }
713  // Non-overlapping patch might not be valid yet so cannot determine
714  // associated patchID
715 }
716 
717 
719 (
720  const cyclicACMIPolyPatch& pp,
721  const polyBoundaryMesh& bm,
722  const label index,
723  const labelUList& mapAddressing,
724  const label newStart
725 )
726 :
727  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
728  nonOverlapPatchName_(pp.nonOverlapPatchName_),
729  nonOverlapPatchID_(-1),
730  srcMask_(),
731  tgtMask_(),
732  srcScalePtr_(pp.srcScalePtr_.clone(*this)),
733  AMITime_
734  (
735  IOobject
736  (
737  "AMITime",
738  boundaryMesh().mesh().pointsInstance(),
739  boundaryMesh().mesh(),
740  IOobject::NO_READ,
741  IOobject::NO_WRITE,
742  IOobject::NO_REGISTER
743  ),
744  dimensionedScalar("time", dimTime, -GREAT)
745  ),
746  prevTimeIndex_(-1)
747 {
748  AMIPtr_->setRequireMatch(false);
749 }
750 
751 
752 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
753 
755 (
756  label& newFaces,
757  label& newProcFaces
758 ) const
759 {
760  const List<labelList>& addSourceFaces = AMI().srcAddress();
761  const scalarField& fMask = srcMask();
762 
763  // Add new faces as many weights for AMI
764  forAll(addSourceFaces, faceI)
765  {
766  if (fMask[faceI] > tolerance_)
767  {
768  const labelList& nbrFaceIs = addSourceFaces[faceI];
769 
770  forAll(nbrFaceIs, j)
771  {
772  label nbrFaceI = nbrFaceIs[j];
773 
774  if (nbrFaceI < neighbPatch().size())
775  {
776  // local faces
777  newFaces++;
778  }
779  else
780  {
781  // Proc faces
782  newProcFaces++;
783  }
784  }
785  }
786  }
787 }
788 
789 
792 {
793  const scalarField& fMask = srcMask();
794  const labelListList& srcFaces = AMI().srcAddress();
795  labelListList dOverFaces;
796 
797  dOverFaces.setSize(srcFaces.size());
798  forAll(dOverFaces, faceI)
799  {
800  if (fMask[faceI] > tolerance_)
801  {
802  dOverFaces[faceI].setSize(srcFaces[faceI].size());
803 
804  forAll(dOverFaces[faceI], subFaceI)
805  {
806  dOverFaces[faceI][subFaceI] = srcFaces[faceI][subFaceI];
807  }
808  }
809  }
810  return refPtr<labelListList>(new labelListList(dOverFaces));
811 }
812 
813 
815 {
816  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
817 
818  // Bit of checking now we know neighbour patch
819  if (!owner() && srcScalePtr_)
820  {
822  << "Ignoring \"scale\" setting in slave patch " << name()
823  << endl;
824  srcScalePtr_.clear();
825  tgtScalePtr_.clear();
826  }
827 
828  return refCast<const cyclicACMIPolyPatch>(pp);
829 }
830 
831 
833 {
834  if (nonOverlapPatchID_ == -1)
835  {
836  nonOverlapPatchID_ =
837  this->boundaryMesh().findPatchID(nonOverlapPatchName_);
838 
839  if (nonOverlapPatchID_ == -1)
840  {
842  << "Illegal non-overlapping patch name " << nonOverlapPatchName_
843  << nl << "Valid patch names are "
844  << this->boundaryMesh().names()
845  << exit(FatalError);
846  }
847 
848  if (nonOverlapPatchID_ < index())
849  {
851  << "Boundary ordering error: " << type()
852  << " patch must be defined prior to its non-overlapping patch"
853  << nl
854  << type() << " patch: " << name() << ", ID:" << index() << nl
855  << "Non-overlap patch: " << nonOverlapPatchName_
856  << ", ID:" << nonOverlapPatchID_ << nl
857  << exit(FatalError);
858  }
859 
860  const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
861 
862  bool ok = true;
863 
864  if (size() == noPp.size())
865  {
866  const scalarField magSf(mag(faceAreas()));
867  const scalarField noMagSf(mag(noPp.faceAreas()));
868 
869  forAll(magSf, facei)
870  {
871  scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
872 
873  if (ratio - 1 > tolerance_)
874  {
875  ok = false;
876  break;
877  }
878  }
879  }
880  else
881  {
882  ok = false;
883  }
884 
885  if (!ok)
886  {
888  << "Inconsistent ACMI patches " << name() << " and "
889  << noPp.name() << ". Patches should have identical topology"
890  << exit(FatalError);
891  }
892  }
893 
894  return nonOverlapPatchID_;
895 }
896 
897 
899 (
900  PstreamBuffers& pBufs,
901  const primitivePatch& pp
902 ) const
903 {
905 }
906 
907 
909 (
910  PstreamBuffers& pBufs,
911  const primitivePatch& pp,
913  labelList& rotation
914 ) const
915 {
916  return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
917 }
918 
919 
921 {
923 
924  os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
925 
926  if (owner() && srcScalePtr_)
927  {
928  srcScalePtr_->writeData(os);
929  }
930 }
931 
932 
933 // ************************************************************************* //
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:598
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:186
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:316
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:309
virtual void scalePatchFaceAreas()
Scale patch face areas to maintain physical area.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:365
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:321
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:627
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
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
#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:66
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:74
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:172
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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 ...