cyclicAMIPolyPatch.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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 "cyclicAMIPolyPatch.H"
30#include "SubField.H"
31#include "Time.H"
32#include "unitConversion.H"
33#include "OFstream.H"
34#include "meshTools.H"
36#include "cyclicPolyPatch.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40namespace Foam
41{
43
46}
47
48const Foam::scalar Foam::cyclicAMIPolyPatch::tolerance_ = 1e-10;
49
50// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51
52Foam::vector Foam::cyclicAMIPolyPatch::findFaceNormalMaxRadius
53(
54 const pointField& faceCentres
55) const
56{
57 // Determine a face furthest away from the axis
58
60
61 const scalarField magRadSqr(magSqr(n));
62
63 label facei = findMax(magRadSqr);
64
66 << "Patch: " << name() << nl
67 << " rotFace = " << facei << nl
68 << " point = " << faceCentres[facei] << nl
69 << " distance = " << Foam::sqrt(magRadSqr[facei])
70 << endl;
71
72 return n[facei];
73}
74
75
77(
78 const primitivePatch& half0,
79 const pointField& half0Ctrs,
80 const vectorField& half0Areas,
81 const pointField& half1Ctrs,
82 const vectorField& half1Areas
83)
84{
85 if (transform() != neighbPatch().transform())
86 {
88 << "Patch " << name()
89 << " has transform type " << transformTypeNames[transform()]
90 << ", neighbour patch " << neighbPatchName()
91 << " has transform type "
92 << neighbPatch().transformTypeNames[neighbPatch().transform()]
93 << exit(FatalError);
94 }
95
96
97 // Calculate transformation tensors
98
99 switch (transform())
100 {
101 case ROTATIONAL:
102 {
103 tensor revT = Zero;
104
105 if (rotationAngleDefined_)
106 {
107 const tensor T(rotationAxis_*rotationAxis_);
108
109 const tensor S
110 (
111 0, -rotationAxis_.z(), rotationAxis_.y(),
112 rotationAxis_.z(), 0, -rotationAxis_.x(),
113 -rotationAxis_.y(), rotationAxis_.x(), 0
114 );
115
116 const tensor revTPos
117 (
118 T
119 + cos(rotationAngle_)*(tensor::I - T)
120 + sin(rotationAngle_)*S
121 );
122
123 const tensor revTNeg
124 (
125 T
126 + cos(-rotationAngle_)*(tensor::I - T)
127 + sin(-rotationAngle_)*S
128 );
129
130 // Check - assume correct angle when difference in face areas
131 // is the smallest
132 const vector transformedAreaPos = gSum(half1Areas & revTPos);
133 const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
134 const vector area0 = gSum(half0Areas);
135 const scalar magArea0 = mag(area0) + ROOTVSMALL;
136
137 // Areas have opposite sign, so sum should be zero when correct
138 // rotation applied
139 const scalar errorPos = mag(transformedAreaPos + area0);
140 const scalar errorNeg = mag(transformedAreaNeg + area0);
141
142 const scalar normErrorPos = errorPos/magArea0;
143 const scalar normErrorNeg = errorNeg/magArea0;
144
145 if (errorPos > errorNeg && normErrorNeg < matchTolerance())
146 {
147 revT = revTNeg;
148 rotationAngle_ *= -1;
149 }
150 else
151 {
152 revT = revTPos;
153 }
154
155 const scalar areaError = min(normErrorPos, normErrorNeg);
156
157 if (areaError > matchTolerance())
158 {
160 << "Patch areas are not consistent within "
161 << 100*matchTolerance()
162 << " % indicating a possible error in the specified "
163 << "angle of rotation" << nl
164 << " owner patch : " << name() << nl
165 << " neighbour patch : " << neighbPatch().name()
166 << nl
167 << " angle : "
168 << radToDeg(rotationAngle_) << " deg" << nl
169 << " area error : " << 100*areaError << " %"
170 << " match tolerance : " << matchTolerance()
171 << endl;
172 }
173
174 if (debug)
175 {
176 scalar theta = radToDeg(rotationAngle_);
177
178 Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
179 << name()
180 << " Specified rotation:"
181 << " swept angle: " << theta << " [deg]"
182 << " reverse transform: " << revT
183 << endl;
184 }
185 }
186 else
187 {
188 point n0 = Zero;
189 point n1 = Zero;
190 if (half0Ctrs.size())
191 {
192 n0 = findFaceNormalMaxRadius(half0Ctrs);
193 }
194 if (half1Ctrs.size())
195 {
196 n1 = -findFaceNormalMaxRadius(half1Ctrs);
197 }
198
199 reduce(n0, maxMagSqrOp<point>());
200 reduce(n1, maxMagSqrOp<point>());
201
202 n0.normalise();
203 n1.normalise();
204
205 // Extended tensor from two local coordinate systems calculated
206 // using normal and rotation axis
207 const tensor E0
208 (
209 rotationAxis_,
210 (n0 ^ rotationAxis_),
211 n0
212 );
213 const tensor E1
214 (
215 rotationAxis_,
216 (-n1 ^ rotationAxis_),
217 -n1
218 );
219 revT = E1.T() & E0;
220
221 if (debug)
222 {
223 scalar theta = radToDeg(acos(-(n0 & n1)));
224
225 Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
226 << name()
227 << " Specified rotation:"
228 << " n0:" << n0 << " n1:" << n1
229 << " swept angle: " << theta << " [deg]"
230 << " reverse transform: " << revT
231 << endl;
232 }
233 }
234
235 const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
236 const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
237 const_cast<vectorField&>(separation()).setSize(0);
238 const_cast<boolList&>(collocated()) = boolList(1, false);
239
240 break;
241 }
242 case TRANSLATIONAL:
243 {
244 if (debug)
245 {
246 Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
247 << " Specified translation : " << separationVector_
248 << endl;
249 }
250
251 const_cast<tensorField&>(forwardT()).clear();
252 const_cast<tensorField&>(reverseT()).clear();
253 const_cast<vectorField&>(separation()) = vectorField
254 (
255 1,
256 separationVector_
257 );
258 const_cast<boolList&>(collocated()) = boolList(1, false);
259
260 break;
261 }
262 default:
263 {
264 if (debug)
265 {
266 Pout<< "patch:" << name()
267 << " Assuming cyclic AMI pairs are colocated" << endl;
268 }
269
270 const_cast<tensorField&>(forwardT()).clear();
271 const_cast<tensorField&>(reverseT()).clear();
272 const_cast<vectorField&>(separation()).setSize(0);
273 const_cast<boolList&>(collocated()) = boolList(1, true);
274
275 break;
276 }
277 }
278
279 if (debug)
280 {
281 Pout<< "patch: " << name() << nl
282 << " forwardT = " << forwardT() << nl
283 << " reverseT = " << reverseT() << nl
284 << " separation = " << separation() << nl
285 << " collocated = " << collocated() << nl << endl;
286 }
287}
288
289
292{
293 const label periodicID = periodicPatchID();
294 if (periodicID != -1)
295 {
296 // Get the periodic patch
297 const coupledPolyPatch& perPp
298 (
299 refCast<const coupledPolyPatch>
300 (
301 boundaryMesh()[periodicID]
302 )
303 );
304 if (!perPp.parallel())
305 {
306 vector axis(Zero);
307 point axisPoint(Zero);
308 if (isA<cyclicPolyPatch>(perPp))
309 {
310 axis =
311 refCast<const cyclicPolyPatch>(perPp).rotationAxis();
312 axisPoint =
313 refCast<const cyclicPolyPatch>(perPp).rotationCentre();
314 }
315 else if (isA<cyclicAMIPolyPatch>(perPp))
316 {
317 axis =
318 refCast<const cyclicAMIPolyPatch>(perPp).rotationAxis();
319 axisPoint =
320 refCast<const cyclicAMIPolyPatch>(perPp).rotationCentre();
321 }
322 else
323 {
324 FatalErrorInFunction << "On patch " << name()
325 << " have unsupported periodicPatch " << perPp.name()
326 << exit(FatalError);
327 }
328
329 return autoPtr<coordSystem::cylindrical>::New(axisPoint, axis);
330 }
331 }
332 return nullptr;
333}
334
335
336// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
337
340{
341 const word surfType(surfDict_.getOrDefault<word>("type", "none"));
342
343 if (!surfPtr_ && owner() && surfType != "none")
344 {
345 word surfName(surfDict_.getOrDefault("name", name()));
346
347 const polyMesh& mesh = boundaryMesh().mesh();
348
349 surfPtr_ =
351 (
352 surfType,
354 (
355 surfName,
356 mesh.time().constant(),
357 "triSurface",
358 mesh,
361 ),
362 surfDict_
363 );
364 }
365
366 return surfPtr_;
367}
368
369
371{
372 resetAMI(boundaryMesh().mesh().points());
373}
374
375
377{
379
380 if (!owner())
381 {
382 return;
383 }
384
385 const cyclicAMIPolyPatch& nbr = neighbPatch();
386 pointField srcPoints(localPoints());
387 pointField nbrPoints(nbr.localPoints());
388
389 if (debug)
390 {
391 const Time& t = boundaryMesh().mesh().time();
392 OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
393 meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
394 }
395
396 label patchSize0 = size();
397 label nbrPatchSize0 = nbr.size();
398
399 if (createAMIFaces_)
400 {
401 // AMI is created based on the original patch faces (non-extended patch)
402 if (srcFaceIDs_.size())
403 {
404 patchSize0 = srcFaceIDs_.size();
405 }
406 if (tgtFaceIDs_.size())
407 {
408 nbrPatchSize0 = tgtFaceIDs_.size();
409 }
410 }
411
412 // Transform neighbour patch to local system
413 transformPosition(nbrPoints);
414 primitivePatch nbrPatch0
415 (
416 SubList<face>(nbr.localFaces(), nbrPatchSize0),
417 nbrPoints
418 );
419 primitivePatch patch0
420 (
421 SubList<face>(localFaces(), patchSize0),
422 srcPoints
423 );
424
425
426 if (debug)
427 {
428 const Time& t = boundaryMesh().mesh().time();
429 OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
430 meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
431
432 OFstream osO(t.path()/name() + "_ownerPatch.obj");
433 meshTools::writeOBJ(osO, this->localFaces(), localPoints());
434 }
435
436 // Construct/apply AMI interpolation to determine addressing and weights
437 AMIPtr_->upToDate() = false;
438 AMIPtr_->calculate(patch0, nbrPatch0, surfPtr());
439
440 if (debug)
441 {
442 AMIPtr_->checkSymmetricWeights(true);
443 }
444}
445
446
448{
450
451 const cyclicAMIPolyPatch& half0 = *this;
452 vectorField half0Areas(half0.size());
453 forAll(half0, facei)
454 {
455 half0Areas[facei] = half0[facei].areaNormal(half0.points());
456 }
457
458 const cyclicAMIPolyPatch& half1 = neighbPatch();
459 vectorField half1Areas(half1.size());
460 forAll(half1, facei)
461 {
462 half1Areas[facei] = half1[facei].areaNormal(half1.points());
463 }
464
465 calcTransforms
466 (
467 half0,
468 half0.faceCentres(),
469 half0Areas,
470 half1.faceCentres(),
471 half1Areas
472 );
473
475 << "calcTransforms() : patch: " << name() << nl
476 << " forwardT = " << forwardT() << nl
477 << " reverseT = " << reverseT() << nl
478 << " separation = " << separation() << nl
479 << " collocated = " << collocated() << nl << endl;
480}
481
482
484{
486
487 // Flag AMI as needing update
488 AMIPtr_->upToDate() = false;
489
491
492 // Early calculation of transforms so e.g. cyclicACMI can use them.
493 // Note: also triggers primitiveMesh face centre. Note that cell
494 // centres should -not- be calculated
495 // since e.g. cyclicACMI overrides face areas
496 calcTransforms();
497}
498
499
501{
503}
504
505
507(
508 PstreamBuffers& pBufs,
509 const pointField& p
510)
511{
513
514 // See below. Clear out any local geometry
516
517 // Note: processorPolyPatch::initMovePoints calls
518 // processorPolyPatch::initGeometry which will trigger calculation of
519 // patch faceCentres() and cell volumes...
520
521 if (createAMIFaces_)
522 {
523 // Note: AMI should have been updated in setTopology
524
525 // faceAreas() and faceCentres() have been reset and will be
526 // recalculated on-demand using the mesh points and no longer
527 // correspond to the scaled areas!
528 restoreScaledGeometry();
529
530 // deltas need to be recalculated to use new face centres!
531 }
532 else
533 {
534 AMIPtr_->upToDate() = false;
535 }
536
537 // Early calculation of transforms. See above.
538 calcTransforms();
539}
540
541
543(
544 PstreamBuffers& pBufs,
545 const pointField& p
546)
547{
549
550// Note: not calling movePoints since this will undo our manipulations!
551// polyPatch::movePoints(pBufs, p);
552/*
553 polyPatch::movePoints
554 -> primitivePatch::movePoints
555 -> primitivePatch::clearGeom:
556 deleteDemandDrivenData(localPointsPtr_);
557 deleteDemandDrivenData(faceCentresPtr_);
558 deleteDemandDrivenData(faceAreasPtr_);
559 deleteDemandDrivenData(magFaceAreasPtr_);
560 deleteDemandDrivenData(faceNormalsPtr_);
561 deleteDemandDrivenData(pointNormalsPtr_);
562*/
563}
564
565
567{
569
571
572 if (createAMIFaces_ && boundaryMesh().mesh().topoChanging() && owner())
573 {
574 setAMIFaces();
575 }
576}
577
578
580{
582
583 // Note: this clears out cellCentres(), faceCentres() and faceAreas()
585}
586
587
589{
591
592 if (!updatingAMI_)
593 {
594 AMIPtr_->upToDate() = false;
595 }
596
598}
599
600
601// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
602
604(
605 const word& name,
606 const label size,
607 const label start,
608 const label index,
609 const polyBoundaryMesh& bm,
610 const word& patchType,
612 const word& defaultAMIMethod
613)
614:
615 coupledPolyPatch(name, size, start, index, bm, patchType, transform),
616 nbrPatchName_(word::null),
617 nbrPatchID_(-1),
618 fraction_(Zero),
619 rotationAxis_(Zero),
620 rotationCentre_(Zero),
621 rotationAngleDefined_(false),
622 rotationAngle_(0.0),
623 separationVector_(Zero),
624 periodicPatchName_(word::null),
625 periodicPatchID_(-1),
626 AMIPtr_(AMIInterpolation::New(defaultAMIMethod)),
627 surfDict_(fileName("surface")),
628 surfPtr_(nullptr),
629 createAMIFaces_(false),
630 moveFaceCentres_(false),
631 updatingAMI_(true),
632 srcFaceIDs_(),
633 tgtFaceIDs_(),
634 faceAreas0_(),
635 faceCentres0_()
636{
637 // Neighbour patch might not be valid yet so no transformation
638 // calculation possible
639}
640
641
643(
644 const word& name,
645 const dictionary& dict,
646 const label index,
647 const polyBoundaryMesh& bm,
648 const word& patchType,
649 const word& defaultAMIMethod
650)
651:
652 coupledPolyPatch(name, dict, index, bm, patchType),
653 nbrPatchName_(dict.getOrDefault<word>("neighbourPatch", word::null)),
654 coupleGroup_(dict),
655 nbrPatchID_(-1),
656 fraction_(dict.getOrDefault<scalar>("fraction", Zero)),
657 rotationAxis_(Zero),
658 rotationCentre_(Zero),
659 rotationAngleDefined_(false),
660 rotationAngle_(0.0),
661 separationVector_(Zero),
662 periodicPatchName_(dict.getOrDefault<word>("periodicPatch", word::null)),
663 periodicPatchID_(-1),
664 AMIPtr_
665 (
667 (
668 dict.getOrDefault<word>("AMIMethod", defaultAMIMethod),
669 dict,
670 dict.getOrDefault("flipNormals", false)
671 )
672 ),
673 surfDict_(dict.subOrEmptyDict("surface")),
674 surfPtr_(nullptr),
675 createAMIFaces_(dict.getOrDefault("createAMIFaces", false)),
676 moveFaceCentres_(false),
677 updatingAMI_(true),
678 srcFaceIDs_(),
679 tgtFaceIDs_(),
680 faceAreas0_(),
681 faceCentres0_()
682{
684 {
686 << "No \"neighbourPatch\" or \"coupleGroup\" provided."
687 << exit(FatalIOError);
688 }
689
690 if (nbrPatchName_ == name)
691 {
693 << "Neighbour patch name " << nbrPatchName_
694 << " cannot be the same as this patch " << name
695 << exit(FatalIOError);
696 }
697
698 switch (transform())
699 {
700 case ROTATIONAL:
701 {
702 dict.readEntry("rotationAxis", rotationAxis_);
703 dict.readEntry("rotationCentre", rotationCentre_);
704 if (dict.readIfPresent("rotationAngle", rotationAngle_))
705 {
708
709 if (debug)
710 {
711 Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
712 << endl;
713 }
714 }
715
716 scalar magRot = mag(rotationAxis_);
717 if (magRot < SMALL)
718 {
720 << "Illegal rotationAxis " << rotationAxis_ << endl
721 << "Please supply a non-zero vector."
722 << exit(FatalIOError);
723 }
724 rotationAxis_ /= magRot;
725
726 break;
727 }
728 case TRANSLATIONAL:
729 {
730 dict.readEntry("separationVector", separationVector_);
731 break;
732 }
733 default:
734 {
735 // No additional info required
736 }
737 }
738
739 // Neighbour patch might not be valid yet so no transformation
740 // calculation possible
741
742 // If topology change, recover the sizes of the original patches and
743 // read additional controls
744 if (createAMIFaces_)
745 {
746 srcFaceIDs_.setSize(dict.get<label>("srcSize"));
747 tgtFaceIDs_.setSize(dict.get<label>("tgtSize"));
748 moveFaceCentres_ = dict.getOrDefault("moveFaceCentres", true);
749 }
750}
751
752
754(
755 const cyclicAMIPolyPatch& pp,
756 const polyBoundaryMesh& bm
757)
758:
759 coupledPolyPatch(pp, bm),
760 nbrPatchName_(pp.nbrPatchName_),
761 coupleGroup_(pp.coupleGroup_),
762 nbrPatchID_(-1),
763 fraction_(pp.fraction_),
764 rotationAxis_(pp.rotationAxis_),
765 rotationCentre_(pp.rotationCentre_),
766 rotationAngleDefined_(pp.rotationAngleDefined_),
767 rotationAngle_(pp.rotationAngle_),
768 separationVector_(pp.separationVector_),
769 periodicPatchName_(pp.periodicPatchName_),
770 periodicPatchID_(-1),
771 AMIPtr_(pp.AMIPtr_->clone()),
772 surfDict_(pp.surfDict_),
773 surfPtr_(nullptr),
774 createAMIFaces_(pp.createAMIFaces_),
775 moveFaceCentres_(pp.moveFaceCentres_),
776 updatingAMI_(true),
777 srcFaceIDs_(),
778 tgtFaceIDs_(),
779 faceAreas0_(),
780 faceCentres0_()
781{
782 // Neighbour patch might not be valid yet so no transformation
783 // calculation possible
784}
785
786
788(
789 const cyclicAMIPolyPatch& pp,
790 const polyBoundaryMesh& bm,
791 const label index,
792 const label newSize,
793 const label newStart,
794 const word& nbrPatchName
795)
796:
797 coupledPolyPatch(pp, bm, index, newSize, newStart),
798 nbrPatchName_(nbrPatchName),
799 coupleGroup_(pp.coupleGroup_),
800 nbrPatchID_(-1),
801 fraction_(pp.fraction_),
802 rotationAxis_(pp.rotationAxis_),
803 rotationCentre_(pp.rotationCentre_),
804 rotationAngleDefined_(pp.rotationAngleDefined_),
805 rotationAngle_(pp.rotationAngle_),
806 separationVector_(pp.separationVector_),
807 periodicPatchName_(pp.periodicPatchName_),
808 periodicPatchID_(-1),
809 AMIPtr_(pp.AMIPtr_->clone()),
810 surfDict_(pp.surfDict_),
811 surfPtr_(nullptr),
812 createAMIFaces_(pp.createAMIFaces_),
813 moveFaceCentres_(pp.moveFaceCentres_),
814 updatingAMI_(true),
815 srcFaceIDs_(),
816 tgtFaceIDs_(),
817 faceAreas0_(),
818 faceCentres0_()
819{
820 if (nbrPatchName_ == name())
821 {
823 << "Neighbour patch name " << nbrPatchName_
824 << " cannot be the same as this patch " << name()
825 << exit(FatalError);
826 }
827
828 // Neighbour patch might not be valid yet so no transformation
829 // calculation possible
830}
831
832
834(
835 const cyclicAMIPolyPatch& pp,
836 const polyBoundaryMesh& bm,
837 const label index,
838 const labelUList& mapAddressing,
839 const label newStart
840)
841:
842 coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
843 nbrPatchName_(pp.nbrPatchName_),
844 coupleGroup_(pp.coupleGroup_),
845 nbrPatchID_(-1),
846 fraction_(pp.fraction_),
847 rotationAxis_(pp.rotationAxis_),
848 rotationCentre_(pp.rotationCentre_),
849 rotationAngleDefined_(pp.rotationAngleDefined_),
850 rotationAngle_(pp.rotationAngle_),
851 separationVector_(pp.separationVector_),
852 periodicPatchName_(pp.periodicPatchName_),
853 periodicPatchID_(-1),
854 AMIPtr_(pp.AMIPtr_->clone()),
855 surfDict_(pp.surfDict_),
856 surfPtr_(nullptr),
857 createAMIFaces_(pp.createAMIFaces_),
858 moveFaceCentres_(pp.moveFaceCentres_),
859 updatingAMI_(true),
860 srcFaceIDs_(),
861 tgtFaceIDs_(),
862 faceAreas0_(),
863 faceCentres0_()
864{}
865
866// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
867
869(
870 label& newFaces,
871 label& newProcFaces
872) const
873{
874 const labelListList& addSourceFaces = AMI().srcAddress();
875
876 // Add new faces as many weights for AMI
877 forAll (addSourceFaces, faceI)
878 {
879 const labelList& nbrFaceIs = addSourceFaces[faceI];
880
881 forAll (nbrFaceIs, j)
882 {
883 label nbrFaceI = nbrFaceIs[j];
884
885 if (nbrFaceI < neighbPatch().size())
886 {
887 // local faces
888 newFaces++;
889 }
890 else
891 {
892 // Proc faces
893 newProcFaces++;
894 }
895 }
896 }
897}
898
899
901{
902 if (nbrPatchID_ == -1)
903 {
904 nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
905
906 if (nbrPatchID_ == -1)
907 {
909 << "Illegal neighbourPatch name " << neighbPatchName()
910 << nl << "Valid patch names are "
911 << this->boundaryMesh().names()
912 << exit(FatalError);
913 }
914
915 // Check that it is a cyclic AMI patch
916 const cyclicAMIPolyPatch& nbrPatch =
917 refCast<const cyclicAMIPolyPatch>
918 (
919 this->boundaryMesh()[nbrPatchID_]
920 );
921
922 if (nbrPatch.neighbPatchName() != name())
923 {
925 << "Patch " << name()
926 << " specifies neighbour patch " << neighbPatchName()
927 << nl << " but that in return specifies "
928 << nbrPatch.neighbPatchName() << endl;
929 }
930 }
931
932 return nbrPatchID_;
933}
934
935
937{
938 if (periodicPatchName_ == word::null)
939 {
940 return -1;
941 }
942 else
943 {
944 if (periodicPatchID_ == -1)
945 {
946 periodicPatchID_ = boundaryMesh().findPatchID(periodicPatchName_);
947
948 if (periodicPatchID_ == -1)
949 {
951 << "Illegal neighbourPatch name " << periodicPatchName_
952 << nl << "Valid patch names are "
953 << this->boundaryMesh().names()
954 << exit(FatalError);
955 }
956 }
957 return periodicPatchID_;
958 }
959}
960
961
963{
964 return index() < neighbPatchID();
965}
966
967
969{
970 const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
971 return refCast<const cyclicAMIPolyPatch>(pp);
972}
973
974
976{
977 if (!owner())
978 {
980 << "AMI interpolator only available to owner patch"
981 << abort(FatalError);
982 }
983
984 if (!AMIPtr_->upToDate())
985 {
986 resetAMI();
987 }
988
989 return *AMIPtr_;
990}
991
992
994{
995 if (owner())
996 {
997 return AMI().applyLowWeightCorrection();
998 }
999 else
1000 {
1001 return neighbPatch().AMI().applyLowWeightCorrection();
1002 }
1003}
1004
1005
1007{
1008 if (!parallel())
1009 {
1010 if (transform() == ROTATIONAL)
1011 {
1012 l = Foam::transform(forwardT(), l - rotationCentre_)
1013 + rotationCentre_;
1014 }
1015 else
1016 {
1017 l = Foam::transform(forwardT(), l);
1018 }
1019 }
1020 else if (separated())
1021 {
1022 // transformPosition gets called on the receiving side,
1023 // separation gets calculated on the sending side so subtract
1024
1025 const vectorField& s = separation();
1026 if (s.size() == 1)
1027 {
1028 forAll(l, i)
1029 {
1030 l[i] -= s[0];
1031 }
1032 }
1033 else
1034 {
1035 l -= s;
1036 }
1037 }
1038}
1039
1040
1042(
1043 point& l,
1044 const label facei
1045) const
1046{
1047 if (!parallel())
1048 {
1049 const tensor& T =
1050 (
1051 forwardT().size() == 1
1052 ? forwardT()[0]
1053 : forwardT()[facei]
1054 );
1055
1056 if (transform() == ROTATIONAL)
1057 {
1058 l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1059 }
1060 else
1061 {
1062 l = Foam::transform(T, l);
1063 }
1064 }
1065 else if (separated())
1066 {
1067 const vector& s =
1068 (
1069 separation().size() == 1
1070 ? separation()[0]
1071 : separation()[facei]
1072 );
1073
1074 l -= s;
1075 }
1076}
1077
1078
1080(
1081 point& l,
1082 const label facei
1083) const
1084{
1085 if (!parallel())
1086 {
1087 const tensor& T =
1088 (
1089 reverseT().size() == 1
1090 ? reverseT()[0]
1091 : reverseT()[facei]
1092 );
1093
1094 if (transform() == ROTATIONAL)
1095 {
1096 l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1097 }
1098 else
1099 {
1100 l = Foam::transform(T, l);
1101 }
1102 }
1103 else if (separated())
1104 {
1105 const vector& s =
1106 (
1107 separation().size() == 1
1108 ? separation()[0]
1109 : separation()[facei]
1110 );
1111
1112 l += s;
1113 }
1114}
1115
1116
1118(
1119 vector& d,
1120 const label facei
1121) const
1122{
1123 if (!parallel())
1124 {
1125 const tensor& T =
1126 (
1127 reverseT().size() == 1
1128 ? reverseT()[0]
1129 : reverseT()[facei]
1130 );
1131
1132 d = Foam::transform(T, d);
1133 }
1134}
1135
1136
1138(
1139 const primitivePatch& referPatch,
1140 const pointField& thisCtrs,
1141 const vectorField& thisAreas,
1142 const pointField& thisCc,
1143 const pointField& nbrCtrs,
1144 const vectorField& nbrAreas,
1145 const pointField& nbrCc
1146)
1147{}
1148
1149
1151(
1152 PstreamBuffers& pBufs,
1153 const primitivePatch& pp
1154) const
1155{}
1156
1157
1159(
1160 PstreamBuffers& pBufs,
1161 const primitivePatch& pp,
1163 labelList& rotation
1164) const
1165{
1166 faceMap.setSize(pp.size());
1167 faceMap = -1;
1168
1169 rotation.setSize(pp.size());
1170 rotation = 0;
1171
1172 return false;
1173}
1174
1175
1177(
1178 const label facei,
1179 const vector& n,
1180 point& p
1181) const
1182{
1183 point prt(p);
1184 reverseTransformPosition(prt, facei);
1185
1186 vector nrt(n);
1187 reverseTransformDirection(nrt, facei);
1188
1189 label nbrFacei = -1;
1190
1191 if (owner())
1192 {
1193 nbrFacei = AMI().tgtPointFace
1194 (
1195 *this,
1196 neighbPatch(),
1197 nrt,
1198 facei,
1199 prt
1200 );
1201 }
1202 else
1203 {
1204 nbrFacei = neighbPatch().AMI().srcPointFace
1205 (
1206 neighbPatch(),
1207 *this,
1208 nrt,
1209 facei,
1210 prt
1211 );
1212 }
1213
1214 if (nbrFacei >= 0)
1215 {
1216 p = prt;
1217 }
1218
1219 return nbrFacei;
1220}
1221
1222
1224{
1226 if (!nbrPatchName_.empty())
1227 {
1228 os.writeEntry("neighbourPatch", nbrPatchName_);
1229 }
1230 coupleGroup_.write(os);
1231
1232 switch (transform())
1233 {
1234 case ROTATIONAL:
1235 {
1236 os.writeEntry("rotationAxis", rotationAxis_);
1237 os.writeEntry("rotationCentre", rotationCentre_);
1238
1239 if (rotationAngleDefined_)
1240 {
1241 os.writeEntry("rotationAngle", radToDeg(rotationAngle_));
1242 }
1243
1244 break;
1245 }
1246 case TRANSLATIONAL:
1247 {
1248 os.writeEntry("separationVector", separationVector_);
1249 break;
1250 }
1251 case NOORDERING:
1252 {
1253 break;
1254 }
1255 default:
1256 {
1257 // No additional info to write
1258 }
1259 }
1260
1261 if (periodicPatchName_ != word::null)
1262 {
1263 os.writeEntry("periodicPatch", periodicPatchName_);
1264 }
1265
1266 AMIPtr_->write(os);
1267
1268 if (!surfDict_.empty())
1269 {
1270 surfDict_.writeEntry(surfDict_.dictName(), os);
1271 }
1272
1273 if (createAMIFaces_)
1274 {
1275 os.writeEntry("createAMIFaces", createAMIFaces_);
1276 os.writeEntry("srcSize", srcFaceIDs_.size());
1277 os.writeEntry("tgtSize", tgtFaceIDs_.size());
1278 os.writeEntry("moveFaceCentres", moveFaceCentres_);
1279 }
1280
1281 os.writeEntryIfDifferent<scalar>("fraction", Zero, fraction_);
1282}
1283
1284
1285// ************************************************************************* //
label n
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...
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
Output to file stream, using an OSstream.
Definition: OFstream.H:57
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
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of faces which address into the list of points.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A List obtained as a section of another List.
Definition: SubList.H:70
static const Tensor I
Definition: Tensor.H:81
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
fileName path() const
Return path.
Definition: Time.H:358
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
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
bool valid() const noexcept
Is a valid patchGroup (non-empty) name.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual transformType transform() const
Type of transform.
virtual bool parallel() const
Are the cyclic planes parallel.
Cyclic patch for Arbitrary Mesh Interface (AMI)
scalar rotationAngle_
Rotation angle.
vector separationVector_
Translation vector.
word nbrPatchName_
Name of other half.
autoPtr< coordSystem::cylindrical > cylindricalCS() const
Create a coordinate system from the periodic patch (or nullptr)
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
label pointFace(const label facei, const vector &n, point &p) const
const word & neighbPatchName() const
Neighbour patch name.
virtual bool owner() const
Does this side own the patch?
virtual void clearGeom()
Clear geometry.
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
const autoPtr< searchableSurface > & surfPtr() const
Create and return pointer to the projection surface.
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
vector rotationAxis_
Axis of rotation for rotational cyclics.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
bool moveFaceCentres_
Move face centres (default = no)
label periodicPatchID() const
Periodic patch ID (or -1)
virtual void calcTransforms()
Recalculate the transformation tensors.
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:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
void calcGeometry()
Calculate the geometry for the patches.
A class for handling file names.
Definition: fileName.H:76
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
Default transformation behaviour for position.
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:102
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:74
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:117
Tensor of scalars, i.e. Tensor<scalar>.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
const volScalarField & T
patchWriters clear()
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
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
#define DebugPout
Report an information message using Foam::Pout.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type gSum(const FieldField< Field, Type > &f)
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
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
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
label findMax(const ListType &input, label start=0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
points setSize(newPointi)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.