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-------------------------------------------------------------------------------
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"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
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
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
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
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
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,
273 mesh.faceAreas(), // already scaled
274 uniqueSort(acmipp.faceCells()), // unique cells only
275 mesh.cells(),
276 const_cast<vectorField&>(mesh.cellCentres()),
277 const_cast<scalarField&>(mesh.cellVolumes())
278 );
279}
280
281
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
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_ = min(scalar(1), max(scalar(0), AMI.srcWeightsSum()));
330 tgtMask_ = min(scalar(1), max(scalar(0), AMI.tgtWeightsSum()));
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 << " mask size:" << patch.srcMask().size() << 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 << " mask size:" << patch.neighbPatch().tgtMask().size()
348 << endl;
349 }
350 }
351}
352
353
355{
356 if (!owner() || !canResetAMI())
357 {
358 return;
359 }
360
361 const polyPatch& nonOverlapPatch = this->nonOverlapPatch();
362 const cyclicACMIPolyPatch& nbrPatch = this->neighbPatch();
363 const polyPatch& nbrNonOverlapPatch = nbrPatch.nonOverlapPatch();
364
365 if (srcScalePtr_)
366 {
367 // Use primitivePatch::faceAreas() to avoid need for additional copy?
368
369 // Save overlap geometry for later scaling
370 thisSf_ = this->faceAreas();
371 thisNoSf_ = nonOverlapPatch.faceAreas();
372 nbrSf_ = nbrPatch.faceAreas();
373 nbrNoSf_ = nbrNonOverlapPatch.faceAreas();
374 }
375
376 // In-place scale the polyPatch face areas
377
378 // Note
379 // - using primitivePatch face areas since these are based on the raw
380 // point locations (not affected by ACMI scaling)
381 scalePatchFaceAreas
382 (
383 *this,
384 srcMask_, // unscaled mask
386 nonOverlapPatch.primitivePatch::faceAreas()
387 );
388 scalePatchFaceAreas
389 (
390 nbrPatch,
391 tgtMask_, // unscaled mask
392 nbrPatch.primitivePatch::faceAreas(),
393 nbrNonOverlapPatch.primitivePatch::faceAreas()
394 );
395
396 // Mark current AMI as up to date with points
397 boundaryMesh().mesh().setUpToDatePoints(AMITime_);
398
399 if (debug)
400 {
401 const auto& acmipp = *this;
402 const auto& mesh = boundaryMesh().mesh();
403 const auto& vols = mesh.cellVolumes();
404
405 Info<< "cyclicACMI PP: " << acmipp.name()
406 << " V:" << sum(scalarField(vols, acmipp.faceCells()))
407 << nl
408 << "cyclicACMI N-O PP: " << nonOverlapPatch.name()
409 << " V:" << sum(scalarField(vols, nonOverlapPatch.faceCells()))
410 << nl
411 << "cyclicACMI NBR PP: " << nbrPatch.name()
412 << " V:" << sum(scalarField(vols, nbrPatch.faceCells()))
413 << nl
414 << "cyclicACMI NBR N-O PP: " << nbrNonOverlapPatch.name()
415 << " V:" << sum(scalarField(vols, nbrNonOverlapPatch.faceCells()))
416 << nl
417 << "cyclicACMI PP+N-O AREA: "
418 << sum(faceAreas() + nonOverlapPatch.faceAreas()) << nl
419 << "cyclicACMI NBR PP+N-O AREA: "
420 << sum(nbrPatch.faceAreas() + nbrNonOverlapPatch.faceAreas())
421 << endl;
422 }
423}
424
425
427{
428 DebugPout << "cyclicACMIPolyPatch::initGeometry : " << name() << endl;
429
430 // Note: calculates transformation and triggers face centre calculation
432
433 // Initialise the AMI early to make sure we adapt the face areas before the
434 // cell centre calculation gets triggered.
435 if (!createAMIFaces_ && canResetAMI())
436 {
437 resetAMI();
438 }
439
440 scalePatchFaceAreas();
441}
442
443
445{
446 DebugPout << "cyclicACMIPolyPatch::calcGeometry : " << name() << endl;
447
449}
450
451
453(
454 PstreamBuffers& pBufs,
455 const pointField& p
456)
457{
458 DebugPout<< "cyclicACMIPolyPatch::initMovePoints : " << name() << endl;
459
460 // Note: calculates transformation and triggers face centre calculation
462
463 if (!createAMIFaces_ && canResetAMI())
464 {
465 resetAMI();
466 }
467
468 scalePatchFaceAreas();
469}
470
471
473(
474 PstreamBuffers& pBufs,
475 const pointField& p
476)
477{
478 DebugPout << "cyclicACMIPolyPatch::movePoints : " << name() << endl;
479
480 // When topology is changing, this will scale the duplicate AMI faces
482}
483
484
486{
487 DebugPout << "cyclicACMIPolyPatch::initUpdateMesh : " << name() << endl;
488
490}
491
492
494{
495 DebugPout << "cyclicACMIPolyPatch::updateMesh : " << name() << endl;
496
498}
499
500
502{
503 DebugPout << "cyclicACMIPolyPatch::clearGeom : " << name() << endl;
504
506}
507
508
510{
511 if (srcScalePtr_)
512 {
513 // Make sure areas are up-to-date
514 updateAreas();
515
516 return srcScaledMask_;
517 }
518 else
519 {
520 return srcMask_;
521 }
522}
523
524
526{
527 if (tgtScalePtr_)
528 {
529 // Make sure areas are up-to-date
530 updateAreas();
531
532 return tgtScaledMask_;
533 }
534 else
535 {
536 return tgtMask_;
537 }
538}
539
540
541// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
542
544(
545 const word& name,
546 const label size,
547 const label start,
548 const label index,
549 const polyBoundaryMesh& bm,
550 const word& patchType,
552 const word& defaultAMIMethod
553)
554:
556 (
557 name,
558 size,
559 start,
560 index,
561 bm,
562 patchType,
563 transform,
564 defaultAMIMethod
565 ),
566 nonOverlapPatchName_(word::null),
567 nonOverlapPatchID_(-1),
568 srcMask_(),
569 tgtMask_(),
570 AMITime_
571 (
573 (
574 "AMITime",
575 boundaryMesh().mesh().pointsInstance(),
576 boundaryMesh().mesh(),
577 IOobject::NO_READ,
578 IOobject::NO_WRITE,
579 false
580 ),
581 dimensionedScalar("time", dimTime, -GREAT)
582 ),
583 prevTimeIndex_(-1)
584{
585 AMIPtr_->setRequireMatch(false);
586
587 // Non-overlapping patch might not be valid yet so cannot determine
588 // associated patchID
589}
590
591
593(
594 const word& name,
595 const dictionary& dict,
596 const label index,
597 const polyBoundaryMesh& bm,
598 const word& patchType,
599 const word& defaultAMIMethod
600)
601:
602 cyclicAMIPolyPatch(name, dict, index, bm, patchType, defaultAMIMethod),
603 nonOverlapPatchName_(dict.get<word>("nonOverlapPatch")),
604 nonOverlapPatchID_(-1),
605 srcMask_(),
606 tgtMask_(),
607 srcScalePtr_
608 (
609 dict.found("scale")
610 ? PatchFunction1<scalar>::New(*this, "scale", dict)
611 : nullptr
612 ),
613 AMITime_
614 (
616 (
617 "AMITime",
618 boundaryMesh().mesh().pointsInstance(),
619 boundaryMesh().mesh(),
620 IOobject::NO_READ,
621 IOobject::NO_WRITE,
622 false
623 ),
624 dimensionedScalar("time", dimTime, -GREAT)
625 ),
626 prevTimeIndex_(-1)
627{
628 AMIPtr_->setRequireMatch(false);
629
630 if (nonOverlapPatchName_ == name)
631 {
633 << "Non-overlapping patch name " << nonOverlapPatchName_
634 << " cannot be the same as this patch " << name
635 << exit(FatalIOError);
636 }
637
638 // Non-overlapping patch might not be valid yet so cannot determine
639 // associated patchID
640}
641
642
644(
645 const cyclicACMIPolyPatch& pp,
646 const polyBoundaryMesh& bm
647)
648:
649 cyclicAMIPolyPatch(pp, bm),
650 nonOverlapPatchName_(pp.nonOverlapPatchName_),
651 nonOverlapPatchID_(-1),
652 srcMask_(),
653 tgtMask_(),
654 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
655 AMITime_
656 (
658 (
659 "AMITime",
660 boundaryMesh().mesh().pointsInstance(),
661 boundaryMesh().mesh(),
662 IOobject::NO_READ,
663 IOobject::NO_WRITE,
664 false
665 ),
666 dimensionedScalar("time", dimTime, -GREAT)
667 ),
668 prevTimeIndex_(-1)
669{
670 AMIPtr_->setRequireMatch(false);
671
672 // Non-overlapping patch might not be valid yet so cannot determine
673 // associated patchID
674}
675
676
678(
679 const cyclicACMIPolyPatch& pp,
680 const polyBoundaryMesh& bm,
681 const label index,
682 const label newSize,
683 const label newStart,
684 const word& nbrPatchName,
685 const word& nonOverlapPatchName
686)
687:
688 cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
689 nonOverlapPatchName_(nonOverlapPatchName),
690 nonOverlapPatchID_(-1),
691 srcMask_(),
692 tgtMask_(),
693 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
694 AMITime_
695 (
697 (
698 "AMITime",
699 boundaryMesh().mesh().pointsInstance(),
700 boundaryMesh().mesh(),
701 IOobject::NO_READ,
702 IOobject::NO_WRITE,
703 false
704 ),
705 dimensionedScalar("time", dimTime, -GREAT)
706 ),
707 prevTimeIndex_(-1)
708{
709 AMIPtr_->setRequireMatch(false);
710
711 if (nonOverlapPatchName_ == name())
712 {
714 << "Non-overlapping patch name " << nonOverlapPatchName_
715 << " cannot be the same as this patch " << name()
716 << exit(FatalError);
717 }
718
719 // Non-overlapping patch might not be valid yet so cannot determine
720 // associated patchID
721}
722
723
725(
726 const cyclicACMIPolyPatch& pp,
727 const polyBoundaryMesh& bm,
728 const label index,
729 const labelUList& mapAddressing,
730 const label newStart
731)
732:
733 cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
734 nonOverlapPatchName_(pp.nonOverlapPatchName_),
735 nonOverlapPatchID_(-1),
736 srcMask_(),
737 tgtMask_(),
738 srcScalePtr_(pp.srcScalePtr_.clone(*this)),
739 AMITime_
740 (
742 (
743 "AMITime",
744 boundaryMesh().mesh().pointsInstance(),
745 boundaryMesh().mesh(),
746 IOobject::NO_READ,
747 IOobject::NO_WRITE,
748 false
749 ),
750 dimensionedScalar("time", dimTime, -GREAT)
751 ),
752 prevTimeIndex_(-1)
753{
754 AMIPtr_->setRequireMatch(false);
755}
756
757
758// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
759
761(
762 label& newFaces,
763 label& newProcFaces
764) const
765{
766 const List<labelList>& addSourceFaces = AMI().srcAddress();
767 const scalarField& fMask = srcMask();
768
769 // Add new faces as many weights for AMI
770 forAll(addSourceFaces, faceI)
771 {
772 if (fMask[faceI] > tolerance_)
773 {
774 const labelList& nbrFaceIs = addSourceFaces[faceI];
775
776 forAll(nbrFaceIs, j)
777 {
778 label nbrFaceI = nbrFaceIs[j];
779
780 if (nbrFaceI < neighbPatch().size())
781 {
782 // local faces
783 newFaces++;
784 }
785 else
786 {
787 // Proc faces
788 newProcFaces++;
789 }
790 }
791 }
792 }
793}
794
795
798{
799 const scalarField& fMask = srcMask();
800 const labelListList& srcFaces = AMI().srcAddress();
801 labelListList dOverFaces;
802
803 dOverFaces.setSize(srcFaces.size());
804 forAll(dOverFaces, faceI)
805 {
806 if (fMask[faceI] > tolerance_)
807 {
808 dOverFaces[faceI].setSize(srcFaces[faceI].size());
809
810 forAll(dOverFaces[faceI], subFaceI)
811 {
812 dOverFaces[faceI][subFaceI] = srcFaces[faceI][subFaceI];
813 }
814 }
815 }
816 return refPtr<labelListList>(new labelListList(dOverFaces));
817}
818
819
821{
822 const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
823
824 // Bit of checking now we know neighbour patch
825 if (!owner() && srcScalePtr_)
826 {
828 << "Ignoring \"scale\" setting in slave patch " << name()
829 << endl;
830 srcScalePtr_.clear();
831 tgtScalePtr_.clear();
832 }
833
834 return refCast<const cyclicACMIPolyPatch>(pp);
835}
836
837
839{
840 if (nonOverlapPatchID_ == -1)
841 {
842 nonOverlapPatchID_ =
843 this->boundaryMesh().findPatchID(nonOverlapPatchName_);
844
845 if (nonOverlapPatchID_ == -1)
846 {
848 << "Illegal non-overlapping patch name " << nonOverlapPatchName_
849 << nl << "Valid patch names are "
850 << this->boundaryMesh().names()
851 << exit(FatalError);
852 }
853
854 if (nonOverlapPatchID_ < index())
855 {
857 << "Boundary ordering error: " << type()
858 << " patch must be defined prior to its non-overlapping patch"
859 << nl
860 << type() << " patch: " << name() << ", ID:" << index() << nl
861 << "Non-overlap patch: " << nonOverlapPatchName_
862 << ", ID:" << nonOverlapPatchID_ << nl
863 << exit(FatalError);
864 }
865
866 const polyPatch& noPp = this->boundaryMesh()[nonOverlapPatchID_];
867
868 bool ok = true;
869
870 if (size() == noPp.size())
871 {
872 const scalarField magSf(mag(faceAreas()));
873 const scalarField noMagSf(mag(noPp.faceAreas()));
874
875 forAll(magSf, facei)
876 {
877 scalar ratio = mag(magSf[facei]/(noMagSf[facei] + ROOTVSMALL));
878
879 if (ratio - 1 > tolerance_)
880 {
881 ok = false;
882 break;
883 }
884 }
885 }
886 else
887 {
888 ok = false;
889 }
890
891 if (!ok)
892 {
894 << "Inconsistent ACMI patches " << name() << " and "
895 << noPp.name() << ". Patches should have identical topology"
896 << exit(FatalError);
897 }
898 }
899
900 return nonOverlapPatchID_;
901}
902
903
905(
906 PstreamBuffers& pBufs,
907 const primitivePatch& pp
908) const
909{
911}
912
913
915(
916 PstreamBuffers& pBufs,
917 const primitivePatch& pp,
919 labelList& rotation
920) const
921{
922 return cyclicAMIPolyPatch::order(pBufs, pp, faceMap, rotation);
923}
924
925
927{
929
930 os.writeEntry("nonOverlapPatch", nonOverlapPatchName_);
931
932 if (owner() && srcScalePtr_)
933 {
934 srcScalePtr_->writeData(os);
935 }
936}
937
938
939// ************************************************************************* //
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
bool upToDate() const
Access to the up-to-date flag.
const scalarField & srcWeightsSum() const
const scalarField & tgtWeightsSum() const
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of faces which address into the list of points.
const Field< point_type > & faceAreas() const
Return face area vectors for patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:31
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
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:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
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 refPtr< labelListList > mapCollocatedFaces() const
Return collocated faces.
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
virtual void clearGeom()
Clear geometry.
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.
virtual label nonOverlapPatchID() const
Non-overlapping patch ID.
void reportCoverage(const word &name, const scalarField &weightSum) const
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
virtual bool updateAreas() const
Update the AMI and patch areas. Return true if anything.
Cyclic patch for Arbitrary Mesh Interface (AMI)
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
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 clearGeom()
Clear geometry.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
void calcGeometry()
Calculate the geometry for the patches.
virtual bool write()
Write the output fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: fvPatch.C:188
order
Enumeration specifying required accuracy.
Definition: meshToMesh0.H:148
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:860
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1096
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const vectorField::subField faceAreas() const
Return face normals.
Definition: polyPatch.C:327
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
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
bool hasFaceAreas() const noexcept
const scalarField & cellVolumes() const
const vectorField & cellCentres() const
const vectorField & faceAreas() const
const cellList & cells() const
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
void setUpToDate()
Set as up-to-date.
Definition: regIOobject.C:411
label eventNo() const noexcept
Event number at last update.
Definition: regIOobjectI.H:191
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugPout
Report an information message using Foam::Pout.
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< T > uniqueSort(const UList< T > &input)
Return sorted list with removal of duplicates.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333