Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
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
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
143bool 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
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>());
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
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
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,
292 mesh.faceAreas(), // already scaled
293 uniqueSort(acmipp.faceCells()), // unique cells only
294 mesh.cells(),
295 const_cast<vectorField&>(mesh.cellCentres()),
296 const_cast<scalarField&>(mesh.cellVolumes())
297 );
298}
299
302{
304}
305
306
308{
309 if (!owner())
310 {
311 return;
312 }
313
314 const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
315
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 const auto& cache = AMI.cache();
340
341 if (cache.index0() == -1 && cache.index1() == -1)
342 {
343 // No caching
344
345 // Output some statistics
346 reportCoverage("source", AMI.srcWeightsSum());
347 reportCoverage("target", AMI.tgtWeightsSum());
348
349 // Set the mask fields
350 // Note:
351 // - assumes that the non-overlap patches are decomposed using the same
352 // decomposition as the coupled patches (per side)
353 srcMask_ = clamp(AMI.srcWeightsSum(), zero_one{});
354 tgtMask_ = clamp(AMI.tgtWeightsSum(), zero_one{});
355 }
356 else
357 {
358 cache.setDirection(owner());
359
360 if (cache.applyLower())
361 {
362 // Output some statistics
363 reportCoverage("source", cache.cSrcWeightsSum0());
364 reportCoverage("target", cache.cTgtWeightsSum0());
365
366 srcMask_ = clamp(cache.cSrcWeightsSum0(), zero_one{});
367 tgtMask_ = clamp(cache.cTgtWeightsSum0(), zero_one{});
368 }
369 else if (cache.applyUpper())
370 {
371 // Output some statistics
372 reportCoverage("source", cache.cSrcWeightsSum1());
373 reportCoverage("target", cache.cTgtWeightsSum1());
374
375 srcMask_ = clamp(cache.cSrcWeightsSum1(), zero_one{});
376 tgtMask_ = clamp(cache.cTgtWeightsSum1(), zero_one{});
377 }
378 else if (cache.applyInterpolate())
379 {
380 srcMask_ = cache.cSrcWeightsSum0();
381 tgtMask_ = cache.cTgtWeightsSum0();
382
383 const auto srcMask1(cache.cSrcWeightsSum1());
384 const auto tgtMask1(cache.cTgtWeightsSum1());
385
386 srcMask_ += (srcMask1 - srcMask_)*cache.weight();
387 tgtMask_ += (tgtMask1 - tgtMask_)*cache.weight();
388
389 // Output some statistics
390 reportCoverage("source", srcMask_);
391 reportCoverage("target", tgtMask_);
392
393 srcMask_ = clamp(srcMask_, zero_one{});
394 tgtMask_ = clamp(tgtMask_, zero_one{});
395 }
396 }
397
398 if (debug)
399 {
400 Pout<< "resetAMI" << endl;
401 {
402 const cyclicACMIPolyPatch& patch = *this;
403 Pout<< "patch:" << patch.name() << " size:" << patch.size()
404 << " non-overlap patch: " << patch.nonOverlapPatch().name()
405 << " size:" << patch.nonOverlapPatch().size()
406 << endl;
407 }
408 {
409 const cyclicACMIPolyPatch& patch = this->neighbPatch();
410 Pout<< "patch:" << patch.name() << " size:" << patch.size()
411 << " non-overlap patch: " << patch.nonOverlapPatch().name()
412 << " size:" << patch.nonOverlapPatch().size()
413 << endl;
414 }
415 }
416}
417
418
420{
421 if (!owner() || !canResetAMI())
422 {
423 return;
424 }
425
426 const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
427 const cyclicACMIPolyPatch& nbrPatch = this->neighbPatch();
428 const polyPatch& nbrNonOverlapPatch = nbrPatch.nonOverlapPatch();
429
430 if (srcScalePtr_)
431 {
432 // Use primitivePatch::faceAreas() to avoid need for additional copy?
433
434 // Save overlap geometry for later scaling
435 thisSf_ = this->faceAreas();
436 thisNoSf_ = nonOverlapPatch.faceAreas();
437 nbrSf_ = nbrPatch.faceAreas();
438 nbrNoSf_ = nbrNonOverlapPatch.faceAreas();
439 }
440
441 // In-place scale the polyPatch face areas
442
443 // Note
444 // - using primitivePatch face areas since these are based on the raw
445 // point locations (not affected by ACMI scaling)
446 scalePatchFaceAreas
447 (
448 *this,
449 srcMask_, // unscaled mask
451 nonOverlapPatch.primitivePatch::faceAreas()
452 );
453 scalePatchFaceAreas
454 (
455 nbrPatch,
456 tgtMask_, // unscaled mask
457 nbrPatch.primitivePatch::faceAreas(),
458 nbrNonOverlapPatch.primitivePatch::faceAreas()
459 );
460
461 // Mark current AMI as up to date with points
462 boundaryMesh().mesh().setUpToDatePoints(AMITime_);
463
464 if (debug)
465 {
466 const auto& acmipp = *this;
467 const auto& mesh = boundaryMesh().mesh();
468 const auto& vols = mesh.cellVolumes();
469
470 Info<< "cyclicACMI PP: " << acmipp.name()
471 << " V:" << sum(scalarField(vols, acmipp.faceCells()))
472 << nl
473 << "cyclicACMI N-O PP: " << nonOverlapPatch.name()
474 << " V:" << sum(scalarField(vols, nonOverlapPatch.faceCells()))
475 << nl
476 << "cyclicACMI NBR PP: " << nbrPatch.name()
477 << " V:" << sum(scalarField(vols, nbrPatch.faceCells()))
478 << nl
479 << "cyclicACMI NBR N-O PP: " << nbrNonOverlapPatch.name()
480 << " V:" << sum(scalarField(vols, nbrNonOverlapPatch.faceCells()))
481 << nl
482 << "cyclicACMI PP+N-O AREA: "
483 << sum(faceAreas() + nonOverlapPatch.faceAreas()) << nl
484 << "cyclicACMI NBR PP+N-O AREA: "
485 << sum(nbrPatch.faceAreas() + nbrNonOverlapPatch.faceAreas())
486 << endl;
487 }
488}
489
490
491void Foam::cyclicACMIPolyPatch::initGeometry(PstreamBuffers& pBufs)
492{
493 DebugPout << "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
494
495 // Note: calculates transformation and triggers face centre calculation
497
498 // Initialise the AMI early to make sure we adapt the face areas before the
499 // cell centre calculation gets triggered.
500 if (!createAMIFaces_ && canResetAMI())
501 {
503 }
504
505 scalePatchFaceAreas();
506}
507
508
509void Foam::cyclicACMIPolyPatch::calcGeometry(PstreamBuffers& pBufs)
511 DebugPout << "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
512
514}
515
516
518(
519 PstreamBuffers& pBufs,
520 const pointField& p
521)
522{
523 DebugPout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
524
525 // Note: calculates transformation and triggers face centre calculation
527
528 if (!createAMIFaces_ && canResetAMI())
529 {
531 }
532
533 scalePatchFaceAreas();
534}
535
536
538(
539 PstreamBuffers& pBufs,
540 const pointField& p
541)
542{
543 DebugPout << "cyclicACMIPolyPatch::movePoints : " << name() << endl;
544
545 // When topology is changing, this will scale the duplicate AMI faces
547}
548
549
552 DebugPout << "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
553
555}
556
557
560 DebugPout << "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
561
563}
564
565
568 DebugPout << "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
569
571}
572
573
575{
576 if (srcScalePtr_)
577 {
578 // Make sure areas are up-to-date
579 updateAreas();
580
581 return srcScaledMask_;
582 }
583 else
584 {
585 return srcMask_;
586 }
587}
588
589
591{
592 if (tgtScalePtr_)
593 {
594 // Make sure areas are up-to-date
595 updateAreas();
596
597 return tgtScaledMask_;
598 }
599 else
600 {
601 return tgtMask_;
602 }
603}
604
605
606// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
607
609(
610 const word& name,
611 const label size,
612 const label start,
613 const label index,
614 const polyBoundaryMesh& bm,
615 const word& patchType,
616 const transformType transform,
617 const word& defaultAMIMethod
618)
619:
620 cyclicAMIPolyPatch
621 (
622 name,
623 size,
624 start,
625 index,
626 bm,
627 patchType,
628 transform,
629 defaultAMIMethod
630 ),
631 nonOverlapPatchName_(),
632 nonOverlapPatchID_(-1),
633 srcMask_(),
634 tgtMask_(),
635 AMITime_
636 (
637 IOobject
638 (
639 "AMITime",
640 boundaryMesh().mesh().pointsInstance(),
641 boundaryMesh().mesh(),
642 IOobject::NO_READ,
643 IOobject::NO_WRITE,
644 IOobject::NO_REGISTER
645 ),
646 dimensionedScalar("time", dimTime, -GREAT)
647 ),
648 prevTimeIndex_(-1)
649{
650 AMIPtr_->setRequireMatch(false);
651
652 // Non-overlapping patch might not be valid yet so cannot determine
653 // associated patchID
654}
655
656
658(
659 const word& name,
660 const dictionary& dict,
661 const label index,
662 const polyBoundaryMesh& bm,
663 const word& patchType,
664 const word& defaultAMIMethod
665)
666:
667 cyclicAMIPolyPatch(name, dict, index, bm, patchType, defaultAMIMethod),
668 nonOverlapPatchName_(dict.get<word>("nonOverlapPatch")),
669 nonOverlapPatchID_(-1),
670 srcMask_(),
671 tgtMask_(),
672 srcScalePtr_(PatchFunction1<scalar>::NewIfPresent(*this, "scale", dict)),
673 AMITime_
674 (
676 (
677 "AMITime",
678 boundaryMesh().mesh().pointsInstance(),
679 boundaryMesh().mesh(),
680 IOobject::NO_READ,
681 IOobject::NO_WRITE,
682 IOobject::NO_REGISTER
683 ),
684 dimensionedScalar("time", dimTime, -GREAT)
685 ),
686 prevTimeIndex_(-1)
687{
688 AMIPtr_->setRequireMatch(false);
689
690 if (nonOverlapPatchName_ == name)
691 {
693 << "Non-overlapping patch name " << nonOverlapPatchName_
694 << " cannot be the same as this patch " << name
695 << exit(FatalIOError);
697
698 // Non-overlapping patch might not be valid yet so cannot determine
699 // associated patchID
700}
701
702
704(
705 const cyclicACMIPolyPatch& pp,
706 const polyBoundaryMesh& bm
707)
708:
709 cyclicAMIPolyPatch(pp, bm),
710 nonOverlapPatchName_(pp.nonOverlapPatchName_),
711 nonOverlapPatchID_(-1),
712 srcMask_(),
713 tgtMask_(),
714 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
715 AMITime_
716 (
717 IOobject
718 (
719 "AMITime",
720 boundaryMesh().mesh().pointsInstance(),
721 boundaryMesh().mesh(),
722 IOobject::NO_READ,
723 IOobject::NO_WRITE,
724 IOobject::NO_REGISTER
725 ),
726 dimensionedScalar("time", dimTime, -GREAT)
727 ),
728 prevTimeIndex_(-1)
729{
730 AMIPtr_->setRequireMatch(false);
731
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 label newSize,
743 const label newStart,
744 const word& nbrPatchName,
745 const word& nonOverlapPatchName
746)
747:
748 cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
749 nonOverlapPatchName_(nonOverlapPatchName),
750 nonOverlapPatchID_(-1),
751 srcMask_(),
752 tgtMask_(),
753 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
754 AMITime_
755 (
757 (
758 "AMITime",
759 boundaryMesh().mesh().pointsInstance(),
760 boundaryMesh().mesh(),
761 IOobject::NO_READ,
762 IOobject::NO_WRITE,
763 IOobject::NO_REGISTER
764 ),
765 dimensionedScalar("time", dimTime, -GREAT)
766 ),
767 prevTimeIndex_(-1)
768{
769 AMIPtr_->setRequireMatch(false);
770
771 if (nonOverlapPatchName_ == name())
772 {
774 << "Non-overlapping patch name " << nonOverlapPatchName_
775 << " cannot be the same as this patch " << name()
776 << exit(FatalError);
778
779 // Non-overlapping patch might not be valid yet so cannot determine
780 // associated patchID
781}
782
783
785(
786 const cyclicACMIPolyPatch& pp,
787 const polyBoundaryMesh& bm,
788 const label index,
789 const labelUList& mapAddressing,
790 const label newStart
791)
792:
793 cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
794 nonOverlapPatchName_(pp.nonOverlapPatchName_),
795 nonOverlapPatchID_(-1),
796 srcMask_(),
797 tgtMask_(),
798 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
799 AMITime_
800 (
801 IOobject
802 (
803 "AMITime",
804 boundaryMesh().mesh().pointsInstance(),
805 boundaryMesh().mesh(),
806 IOobject::NO_READ,
807 IOobject::NO_WRITE,
808 IOobject::NO_REGISTER
809 ),
810 dimensionedScalar("time", dimTime, -GREAT)
811 ),
812 prevTimeIndex_(-1)
814 AMIPtr_->setRequireMatch(false);
815}
816
817
818// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
819
821(
822 label& newFaces,
823 label& newProcFaces
824) const
825{
826 const List<labelList>& addSourceFaces = AMI().srcAddress();
827 const scalarField& fMask = srcMask();
828
829 // Add new faces as many weights for AMI
830 forAll(addSourceFaces, faceI)
831 {
832 if (fMask[faceI] > tolerance_)
833 {
834 const labelList& nbrFaceIs = addSourceFaces[faceI];
835
836 forAll(nbrFaceIs, j)
837 {
838 label nbrFaceI = nbrFaceIs[j];
839
840 if (nbrFaceI < neighbPatch().size())
841 {
842 // local faces
843 newFaces++;
844 }
845 else
846 {
847 // Proc faces
848 newProcFaces++;
849 }
851 }
852 }
853}
854
855
856Foam::refPtr<Foam::labelListList>
858{
859 const scalarField& fMask = srcMask();
860 const labelListList& srcFaces = AMI().srcAddress();
861 labelListList dOverFaces;
862
863 dOverFaces.setSize(srcFaces.size());
864 forAll(dOverFaces, faceI)
865 {
866 if (fMask[faceI] > tolerance_)
867 {
868 dOverFaces[faceI].setSize(srcFaces[faceI].size());
869
870 forAll(dOverFaces[faceI], subFaceI)
871 {
872 dOverFaces[faceI][subFaceI] = srcFaces[faceI][subFaceI];
874 }
875 }
876 return refPtr<labelListList>(new labelListList(dOverFaces));
877}
878
879
880const Foam::cyclicACMIPolyPatch& Foam::cyclicACMIPolyPatch::neighbPatch() const
881{
882 const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
883
884 // Bit of checking now we know neighbour patch
885 if (!owner() && srcScalePtr_)
886 {
888 << "Ignoring \"scale\" setting in slave patch " << name()
889 << endl;
890 srcScalePtr_.clear();
891 tgtScalePtr_.clear();
892 }
893
895}
896
897
899{
900 if (nonOverlapPatchID_ == -1)
901 {
902 nonOverlapPatchID_ =
903 this->boundaryMesh().findPatchID(nonOverlapPatchName_);
904
905 if (nonOverlapPatchID_ == -1)
906 {
908 << "Illegal non-overlapping patch name " << nonOverlapPatchName_
909 << nl << "Valid patch names are "
910 << this->boundaryMesh().names()
911 << exit(FatalError);
912 }
913
914 if (nonOverlapPatchID_ < index())
915 {
917 << "Boundary ordering error: " << type()
918 << " patch must be defined prior to its non-overlapping patch"
919 << nl
920 << type() << " patch: " << name() << ", ID:" << index() << nl
921 << "Non-overlap patch: " << nonOverlapPatchName_
922 << ", ID:" << nonOverlapPatchID_ << nl
923 << exit(FatalError);
924 }
925
926 const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
927
928 bool ok = true;
929
930 if (size() == noPp.size())
931 {
932 const scalarField magSf(mag(faceAreas()));
933 const scalarField noMagSf(mag(noPp.faceAreas()));
934
935 forAll(magSf, facei)
936 {
937 scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
938
939 if (ratio - 1 > tolerance_)
940 {
941 ok = false;
942 break;
943 }
944 }
945 }
946 else
947 {
948 ok = false;
949 }
950
951 if (!ok)
952 {
954 << "Inconsistent ACMI patches " << name() << " and "
955 << noPp.name() << ". Patches should have identical topology"
956 << exit(FatalError);
958 }
959
960 return nonOverlapPatchID_;
961}
962
963
965(
966 PstreamBuffers& pBufs,
968) const
969{
971}
972
973
975(
976 PstreamBuffers& pBufs,
977 const primitivePatch& pp,
979 labelList& rotation
980) const
981{
982 return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
983}
984
985
987{
989
990 os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
991
992 if (owner() && srcScalePtr_)
993 {
994 srcScalePtr_->writeData(os);
995 }
996}
997
998
999// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
SubField< vector > subField
Definition Field.H:183
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
void setUpToDate(regIOobject &) const
Set object up to date with *this.
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
virtual const scalarField & tgtMask() const
Return the mask/weighting for the target patch.
virtual void scalePatchFaceAreas()
Scale patch face areas to maintain physical area.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
virtual refPtr< labelListList > mapCollocatedFaces() const
Return collocated faces.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
virtual void clearGeom()
Clear geometry.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal sub-faces and new proc faces.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
virtual const cyclicACMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual const scalarField & srcMask() const
Return the mask/weighting for the source patch.
bool upToDate(const regIOobject &) const
Return true if given object is up to date with *this.
virtual label nonOverlapPatchID() const
Non-overlapping patch ID.
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.
void reportCoverage(const word &name, const scalarField &weightSum) const
Report coverage statics, e.g. number of uncovered/blended/covered faces.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
const word & nonOverlapPatchName() const
Non-overlapping patch name.
virtual void resetAMI(const UList< point > &points) const
Reset the AMI interpolator, supply patch points.
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything.
Cyclic patch for Arbitrary Mesh Interface (AMI).
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
const scalarField & weightsSum() const
Helper function to return the weights sum.
const scalarListList & weights() const
Helper function to return the weights.
virtual bool owner() const
Does this side own the patch?
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
bool canResetAMI() const
Flag to indicate whether the AMI can be reset.
virtual void clearGeom()
Clear geometry.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not refer to *this (except for name() and type() etc....
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
cyclicAMIPolyPatch(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.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
virtual label neighbPatchID() const
Neighbour patch ID.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyMesh & mesh() const noexcept
Return the mesh reference.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField::subField faceAreas() const
Return face normals.
Definition polyPatch.C:326
tmp< scalarField > areaFraction(const pointField &points) const
Calculate the area fraction as the ratio of the stored face area and the area given by the face point...
Definition polyPatch.C:352
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const labelUList & faceCells() const
Return face-cell addressing.
Definition polyPatch.C:401
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.
const vectorField & faceCentres() const
const scalarField & cellVolumes() const
const vectorField & faceAreas() const
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
label eventNo() const noexcept
Event number at last update.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugPout
Report an information message using Foam::Pout.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
List< T > uniqueSort(const UList< T > &input)
Return sorted list with removal of duplicates.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
List< vector > vectorList
List of vector.
Definition vectorList.H:32
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299