Loading...
Searching...
No Matches
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-2023 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"
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
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()).clear();
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()).clear();
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 refCast<const coupledPolyPatch>(boundaryMesh()[periodicID]);
299
300 if (!perPp.parallel())
301 {
302 vector axis(Zero);
303 point axisPoint(Zero);
304
305 if
306 (
307 const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(perPp)
308 )
309 {
310 axis = cpp->rotationAxis();
311 axisPoint = cpp->rotationCentre();
312 }
313 else if
314 (
315 const cyclicAMIPolyPatch* cpp = isA<cyclicAMIPolyPatch>(perPp)
316 )
317 {
318 axis = cpp->rotationAxis();
319 axisPoint = cpp->rotationCentre();
320 }
321 else
322 {
323 FatalErrorInFunction << "On patch " << name()
324 << " have unsupported periodicPatch " << perPp.name()
325 << exit(FatalError);
326 }
327
328 return autoPtr<coordSystem::cylindrical>::New(axisPoint, axis);
329 }
330 }
331 return nullptr;
332}
333
334
335// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
336
339{
340 const word surfType(surfDict_.getOrDefault<word>("type", "none"));
341
342 if (!surfPtr_ && owner() && surfType != "none")
343 {
344 word surfName(surfDict_.getOrDefault("name", name()));
345
346 const polyMesh& mesh = boundaryMesh().mesh();
347
348 surfPtr_ =
350 (
351 surfType,
353 (
354 surfName,
355 mesh.time().constant(),
356 "triSurface",
357 mesh,
360 ),
361 surfDict_
362 );
363 }
364
365 return surfPtr_;
366}
367
370{
372}
373
374
376{
378
379 if (!owner())
380 {
381 return;
382 }
383
384 point refPt = point::max;
385 bool restoredFromCache = false;
386 if (AMIPtr_->cacheActive())
387 {
388 const auto& mesh = boundaryMesh().mesh();
389 const label comm = mesh.comm();
390
391 if (UPstream::parRun())
392 {
393 label refProci = -1;
394 if (size() > 0)
395 {
396 refProci = UPstream::myProcNo(comm);
397 }
398 reduce(refProci, maxOp<label>(), UPstream::msgType(), comm);
399
400 if (refProci == -1)
401 {
402 // No AMI faces
403 Info<< "AMI: Patch " << name() << " has no faces "
404 << "- deactivating cache" << endl;
405 AMIPtr_->cache().setSize(-1);
406 }
407 else
408 {
409 if (refProci == UPstream::myProcNo(comm))
410 {
411 const label meshPointi = meshPoints()[0];
412 refPt = points[meshPointi];
413
414 // Need to ensure that this patch is on the moving side!
415 const scalar d(magSqr(refPt - mesh.oldPoints()[meshPointi]));
416
417 if (d < tolerance_)
418 {
420 << "Attempted to use AMI caching on a static patch "
421 << name()
422 << ". Potential patch ordering issue - "
423 << "flip owner and neighbour patches?"
424 << endl;
425 }
426 }
427 reduce(refPt, minOp<point>(), UPstream::msgType(), comm);
428 }
429 }
430 else
431 {
432 if (size())
433 {
434 refPt = points[meshPoints()[0]];
435 }
436 else
437 {
438 // No AMI faces
439 Info<< "AMI: Patch " << name() << " has no faces "
440 << "- deactivating cache" << endl;
441 AMIPtr_->cache().setSize(-1);
442 }
443 }
444
445 // Sets cache indices to use and time interpolation weight
446 restoredFromCache = AMIPtr_->restoreCache(refPt);
447
448 if (returnReduceOr(restoredFromCache, comm))
449 {
450 // Restored AMI weight and addressing from cache - all done
451 Info<< "AMI: Cached weights and addresses restored" << endl;
452 return;
453 }
454
455 Info<< "AMI: Cached weights and addresses unavailable" << endl;
456 }
457
458 const cyclicAMIPolyPatch& nbr = neighbPatch();
459 const pointField srcPoints(points, meshPoints());
460 pointField nbrPoints(points, nbr.meshPoints());
461
462 Info<< "AMI: Creating AMI for source:" << name()
463 << " and target:" << nbr.name() << endl;
464
465 if (debug)
466 {
467 const Time& t = boundaryMesh().mesh().time();
468 OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
469 meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
470 }
471
472 label patchSize0 = size();
473 label nbrPatchSize0 = nbr.size();
474
475 if (createAMIFaces_)
476 {
477 // AMI is created based on the original patch faces (non-extended patch)
478 if (srcFaceIDs_.size())
479 {
480 patchSize0 = srcFaceIDs_.size();
481 }
482 if (tgtFaceIDs_.size())
483 {
484 nbrPatchSize0 = tgtFaceIDs_.size();
485 }
486 }
487
488 // Transform neighbour patch to local system
489 transformPosition(nbrPoints);
490 primitivePatch nbrPatch0
491 (
492 SubList<face>(nbr.localFaces(), nbrPatchSize0),
493 nbrPoints
494 );
495 primitivePatch patch0
496 (
497 SubList<face>(localFaces(), patchSize0),
498 srcPoints
499 );
500
501 if (debug)
502 {
503 const Time& t = boundaryMesh().mesh().time();
504 OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
505 meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
506
507 OFstream osO(t.path()/name() + "_ownerPatch.obj");
508 meshTools::writeOBJ(osO, this->localFaces(), srcPoints);
509 }
510
511 // Construct/apply AMI interpolation to determine addressing and weights
512 AMIPtr_->upToDate(false);
513
514 // Note: e.g. redistributePar might construct the AMI with a different
515 // worldComm so reset it to the mesh.comm.
516 AMIPtr_->comm(boundaryMesh().mesh().comm());
517 AMIPtr_->calculate(patch0, nbrPatch0, surfPtr());
518
519 if (debug)
520 {
521 AMIPtr_->checkSymmetricWeights(true);
522 }
523
524 AMIPtr_->addToCache(refPt);
525}
526
527
529{
531
532 const cyclicAMIPolyPatch& half0 = *this;
533 vectorField half0Areas(half0.size());
534 forAll(half0, facei)
535 {
536 half0Areas[facei] = half0[facei].areaNormal(half0.points());
537 }
538
539 const cyclicAMIPolyPatch& half1 = neighbPatch();
540 vectorField half1Areas(half1.size());
541 forAll(half1, facei)
542 {
543 half1Areas[facei] = half1[facei].areaNormal(half1.points());
544 }
545
546 calcTransforms
547 (
548 half0,
549 half0.faceCentres(),
550 half0Areas,
551 half1.faceCentres(),
552 half1Areas
553 );
554
556 << "calcTransforms() : patch: " << name() << nl
557 << " forwardT = " << forwardT() << nl
558 << " reverseT = " << reverseT() << nl
559 << " separation = " << separation() << nl
560 << " collocated = " << collocated() << nl << endl;
561}
562
563
565{
567
568 // Flag AMI as needing update
569 AMIPtr_->upToDate(false);
570
572
573 // Early calculation of transforms so e.g. cyclicACMI can use them.
574 // Note: also triggers primitiveMesh face centre. Note that cell
575 // centres should -not- be calculated
576 // since e.g. cyclicACMI overrides face areas
577 calcTransforms();
578}
579
582{
584}
585
586
588(
589 PstreamBuffers& pBufs,
590 const pointField& p
591)
592{
594
595 // See below. Clear out any local geometry
597
598 // Note: processorPolyPatch::initMovePoints calls
599 // processorPolyPatch::initGeometry which will trigger calculation of
600 // patch faceCentres() and cell volumes...
601
602 if (createAMIFaces_)
603 {
604 // Note: AMI should have been updated in setTopology
605
606 // faceAreas() and faceCentres() have been reset and will be
607 // recalculated on-demand using the mesh points and no longer
608 // correspond to the scaled areas!
609 restoreScaledGeometry();
610
611 // deltas need to be recalculated to use new face centres!
612 }
613 else
614 {
615 AMIPtr_->upToDate(false);
617
618 // Early calculation of transforms. See above.
619 calcTransforms();
620}
621
622
624(
625 PstreamBuffers& pBufs,
626 const pointField& p
627)
628{
630
631// Note: not calling movePoints since this will undo our manipulations!
632// polyPatch::movePoints(pBufs, p);
633/*
634 polyPatch::movePoints
635 -> primitivePatch::movePoints
636 -> primitivePatch::clearGeom:
637 localPointsPtr_.reset(nullptr);
638 faceCentresPtr_.reset(nullptr);
639 faceAreasPtr_.reset(nullptr);
640 magFaceAreasPtr_.reset(nullptr);
641 faceNormalsPtr_.reset(nullptr);
642 pointNormalsPtr_.reset(nullptr);
643*/
644}
645
646
648{
650
652
653 if (createAMIFaces_ && boundaryMesh().mesh().topoChanging() && owner())
654 {
655 setAMIFaces();
656 }
657}
658
659
660void Foam::cyclicAMIPolyPatch::updateMesh(PstreamBuffers& pBufs)
661{
663
664 // Note: this clears out cellCentres(), faceCentres() and faceAreas()
666}
667
668
670{
672
673 if (!updatingAMI_)
674 {
675 AMIPtr_->upToDate(false);
676 }
679}
680
681
682// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
683
685(
686 const word& name,
687 const label size,
688 const label start,
689 const label index,
690 const polyBoundaryMesh& bm,
691 const word& patchType,
692 const transformType transform,
693 const word& defaultAMIMethod
694)
695:
696 coupledPolyPatch(name, size, start, index, bm, patchType, transform),
697 nbrPatchName_(),
698 nbrPatchID_(-1),
699 fraction_(Zero),
700 rotationAxis_(Zero),
701 rotationCentre_(Zero),
702 rotationAngleDefined_(false),
703 rotationAngle_(0.0),
704 separationVector_(Zero),
705 periodicPatchName_(),
706 periodicPatchID_(-1),
707 AMIPtr_(AMIInterpolation::New(defaultAMIMethod)),
708 surfDict_(fileName("surface")),
709 surfPtr_(nullptr),
710 createAMIFaces_(false),
711 moveFaceCentres_(false),
712 updatingAMI_(true),
713 srcFaceIDs_(),
714 tgtFaceIDs_(),
715 faceAreas0_(),
717{
718 // Neighbour patch might not be valid yet so no transformation
719 // calculation possible
720}
721
722
724(
725 const word& name,
726 const dictionary& dict,
727 const label index,
728 const polyBoundaryMesh& bm,
729 const word& patchType,
730 const word& defaultAMIMethod
731)
732:
733 coupledPolyPatch(name, dict, index, bm, patchType),
734 nbrPatchName_(dict.getOrDefault<word>("neighbourPatch", word::null)),
735 coupleGroup_(dict),
736 nbrPatchID_(-1),
737 fraction_(dict.getOrDefault<scalar>("fraction", Zero)),
738 rotationAxis_(Zero),
739 rotationCentre_(Zero),
740 rotationAngleDefined_(false),
741 rotationAngle_(0.0),
742 separationVector_(Zero),
743 periodicPatchName_(dict.getOrDefault<word>("periodicPatch", word::null)),
744 periodicPatchID_(-1),
745 AMIPtr_
746 (
748 (
749 dict.getOrDefault<word>("AMIMethod", defaultAMIMethod),
750 dict,
751 dict.getOrDefault("flipNormals", false)
752 )
753 ),
754 surfDict_(dict.subOrEmptyDict("surface")),
755 surfPtr_(nullptr),
756 createAMIFaces_(dict.getOrDefault("createAMIFaces", false)),
757 moveFaceCentres_(false),
758 updatingAMI_(true),
759 srcFaceIDs_(),
760 tgtFaceIDs_(),
761 faceAreas0_(),
762 faceCentres0_()
763{
764 if (nbrPatchName_.empty() && !coupleGroup_.good())
765 {
766 FatalIOErrorInFunction(dict)
767 << "No \"neighbourPatch\" or \"coupleGroup\" provided."
768 << exit(FatalIOError);
769 }
770
772 {
774 << "Neighbour patch name " << nbrPatchName_
775 << " cannot be the same as this patch " << name
776 << exit(FatalIOError);
777 }
778
779 switch (transform())
780 {
781 case ROTATIONAL:
782 {
783 dict.readEntry("rotationAxis", rotationAxis_);
784 dict.readEntry("rotationCentre", rotationCentre_);
785 if (dict.readIfPresent("rotationAngle", rotationAngle_))
786 {
789
790 if (debug)
791 {
792 Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
793 << endl;
794 }
795 }
796
797 scalar magRot = mag(rotationAxis_);
798 if (magRot < SMALL)
799 {
801 << "Illegal rotationAxis " << rotationAxis_ << endl
802 << "Please supply a non-zero vector."
803 << exit(FatalIOError);
804 }
805 rotationAxis_ /= magRot;
806
807 break;
808 }
809 case TRANSLATIONAL:
810 {
811 dict.readEntry("separationVector", separationVector_);
812 break;
813 }
814 default:
815 {
816 // No additional info required
817 }
818 }
819
820 // Neighbour patch might not be valid yet so no transformation
821 // calculation possible
822
823 // If topology change, recover the sizes of the original patches and
824 // read additional controls
825 if (createAMIFaces_)
826 {
827 if (AMIPtr_->cacheActive())
828 {
830 << "createAMIFaces and AMI caching cannot be used together."
831 << exit(FatalIOError);
832 }
833
834 srcFaceIDs_.setSize(dict.get<label>("srcSize"));
835 tgtFaceIDs_.setSize(dict.get<label>("tgtSize"));
836 moveFaceCentres_ = dict.getOrDefault("moveFaceCentres", true);
837 }
838}
839
840
842(
843 const cyclicAMIPolyPatch& pp,
844 const polyBoundaryMesh& bm
845)
846:
847 coupledPolyPatch(pp, bm),
848 nbrPatchName_(pp.nbrPatchName_),
849 coupleGroup_(pp.coupleGroup_),
850 nbrPatchID_(-1),
851 fraction_(pp.fraction_),
852 rotationAxis_(pp.rotationAxis_),
853 rotationCentre_(pp.rotationCentre_),
854 rotationAngleDefined_(pp.rotationAngleDefined_),
855 rotationAngle_(pp.rotationAngle_),
856 separationVector_(pp.separationVector_),
857 periodicPatchName_(pp.periodicPatchName_),
858 periodicPatchID_(-1),
859 AMIPtr_(pp.AMIPtr_->clone()),
860 surfDict_(pp.surfDict_),
861 surfPtr_(nullptr),
862 createAMIFaces_(pp.createAMIFaces_),
863 moveFaceCentres_(pp.moveFaceCentres_),
864 updatingAMI_(true),
865 srcFaceIDs_(),
866 tgtFaceIDs_(),
867 faceAreas0_(),
869{
870 // Neighbour patch might not be valid yet so no transformation
871 // calculation possible
872}
873
874
876(
877 const cyclicAMIPolyPatch& pp,
878 const polyBoundaryMesh& bm,
879 const label index,
880 const label newSize,
881 const label newStart,
882 const word& nbrPatchName
883)
884:
885 coupledPolyPatch(pp, bm, index, newSize, newStart),
886 nbrPatchName_(nbrPatchName),
887 coupleGroup_(pp.coupleGroup_),
888 nbrPatchID_(-1),
889 fraction_(pp.fraction_),
890 rotationAxis_(pp.rotationAxis_),
891 rotationCentre_(pp.rotationCentre_),
892 rotationAngleDefined_(pp.rotationAngleDefined_),
893 rotationAngle_(pp.rotationAngle_),
894 separationVector_(pp.separationVector_),
895 periodicPatchName_(pp.periodicPatchName_),
896 periodicPatchID_(-1),
897 AMIPtr_(pp.AMIPtr_->clone()),
898 surfDict_(pp.surfDict_),
899 surfPtr_(nullptr),
900 createAMIFaces_(pp.createAMIFaces_),
901 moveFaceCentres_(pp.moveFaceCentres_),
902 updatingAMI_(true),
903 srcFaceIDs_(),
904 tgtFaceIDs_(),
905 faceAreas0_(),
906 faceCentres0_()
907{
908 if (nbrPatchName_ == name())
909 {
911 << "Neighbour patch name " << nbrPatchName_
912 << " cannot be the same as this patch " << name()
913 << exit(FatalError);
915
916 // Neighbour patch might not be valid yet so no transformation
917 // calculation possible
918}
919
920
922(
923 const cyclicAMIPolyPatch& pp,
924 const polyBoundaryMesh& bm,
925 const label index,
926 const labelUList& mapAddressing,
927 const label newStart
928)
929:
930 coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
931 nbrPatchName_(pp.nbrPatchName_),
932 coupleGroup_(pp.coupleGroup_),
933 nbrPatchID_(-1),
934 fraction_(pp.fraction_),
935 rotationAxis_(pp.rotationAxis_),
936 rotationCentre_(pp.rotationCentre_),
937 rotationAngleDefined_(pp.rotationAngleDefined_),
938 rotationAngle_(pp.rotationAngle_),
939 separationVector_(pp.separationVector_),
940 periodicPatchName_(pp.periodicPatchName_),
941 periodicPatchID_(-1),
942 AMIPtr_(pp.AMIPtr_->clone()),
943 surfDict_(pp.surfDict_),
944 surfPtr_(nullptr),
945 createAMIFaces_(pp.createAMIFaces_),
946 moveFaceCentres_(pp.moveFaceCentres_),
947 updatingAMI_(true),
948 srcFaceIDs_(),
950 faceAreas0_(),
952{}
953
954// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
955
957(
958 label& newFaces,
959 label& newProcFaces
960) const
961{
962 const labelListList& addSourceFaces = AMI().srcAddress();
963
964 // Add new faces as many weights for AMI
965 forAll (addSourceFaces, faceI)
966 {
967 const labelList& nbrFaceIs = addSourceFaces[faceI];
968
969 forAll (nbrFaceIs, j)
970 {
971 label nbrFaceI = nbrFaceIs[j];
972
973 if (nbrFaceI < neighbPatch().size())
974 {
975 // local faces
976 newFaces++;
977 }
978 else
979 {
980 // Proc faces
981 newProcFaces++;
982 }
983 }
984 }
985}
986
987
989{
990 if (nbrPatchID_ == -1)
991 {
992 nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
993
994 if (nbrPatchID_ == -1)
995 {
997 << "Illegal neighbourPatch name " << neighbPatchName()
998 << nl << "Valid patch names are "
999 << this->boundaryMesh().names()
1000 << exit(FatalError);
1001 }
1002
1003 // Check that it is a cyclic AMI patch
1004 const cyclicAMIPolyPatch& nbrPatch =
1006 (
1007 this->boundaryMesh()[nbrPatchID_]
1008 );
1009
1010 if (nbrPatch.neighbPatchName() != name())
1011 {
1013 << "Patch " << name()
1014 << " specifies neighbour patch " << neighbPatchName()
1015 << nl << " but that in return specifies "
1016 << nbrPatch.neighbPatchName() << endl;
1018 }
1019
1020 return nbrPatchID_;
1021}
1022
1023
1025{
1026 if (periodicPatchName_.empty())
1027 {
1028 return -1;
1029 }
1030 else
1031 {
1032 if (periodicPatchID_ == -1)
1033 {
1034 periodicPatchID_ = boundaryMesh().findPatchID(periodicPatchName_);
1035
1036 if (periodicPatchID_ == -1)
1037 {
1039 << "Illegal neighbourPatch name " << periodicPatchName_
1040 << nl << "Valid patch names are "
1041 << this->boundaryMesh().names()
1042 << exit(FatalError);
1044 }
1045 return periodicPatchID_;
1046 }
1047}
1048
1051{
1052 return index() < neighbPatchID();
1053}
1054
1055
1057{
1058 const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
1060}
1061
1062
1064{
1065 if (!owner())
1066 {
1068 << "AMI interpolator only available to owner patch"
1069 << abort(FatalError);
1070 }
1071
1072 if (!AMIPtr_->upToDate())
1073 {
1075 }
1076
1077 return *AMIPtr_;
1078}
1079
1080
1082{
1083 if (owner())
1084 {
1085 return AMI().applyLowWeightCorrection();
1086 }
1087 else
1088 {
1089 return neighbPatch().AMI().applyLowWeightCorrection();
1090 }
1091}
1092
1093
1095{
1096 if (!parallel())
1097 {
1098 if (transform() == ROTATIONAL)
1099 {
1100 l = Foam::transform(forwardT(), l - rotationCentre_)
1101 + rotationCentre_;
1102 }
1103 else
1104 {
1105 l = Foam::transform(forwardT(), l);
1106 }
1107 }
1108 else if (separated())
1109 {
1110 // transformPosition gets called on the receiving side,
1111 // separation gets calculated on the sending side so subtract
1112
1113 const vectorField& s = separation();
1114 if (s.size() == 1)
1115 {
1116 forAll(l, i)
1117 {
1118 l[i] -= s[0];
1119 }
1120 }
1121 else
1123 l -= s;
1124 }
1125 }
1126}
1127
1128
1130(
1131 point& l,
1132 const label facei
1133) const
1134{
1135 if (!parallel())
1136 {
1137 const tensor& T =
1138 (
1139 forwardT().size() == 1
1140 ? forwardT()[0]
1141 : forwardT()[facei]
1142 );
1143
1144 if (transform() == ROTATIONAL)
1145 {
1146 l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1147 }
1148 else
1149 {
1150 l = Foam::transform(T, l);
1151 }
1152 }
1153 else if (separated())
1154 {
1155 const vector& s =
1156 (
1157 separation().size() == 1
1158 ? separation()[0]
1159 : separation()[facei]
1161
1162 l -= s;
1163 }
1164}
1165
1166
1168(
1169 point& l,
1170 const label facei
1171) const
1172{
1173 if (!parallel())
1174 {
1175 const tensor& T =
1176 (
1177 reverseT().size() == 1
1178 ? reverseT()[0]
1179 : reverseT()[facei]
1180 );
1181
1182 if (transform() == ROTATIONAL)
1183 {
1184 l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1185 }
1186 else
1187 {
1188 l = Foam::transform(T, l);
1189 }
1190 }
1191 else if (separated())
1192 {
1193 const vector& s =
1194 (
1195 separation().size() == 1
1196 ? separation()[0]
1197 : separation()[facei]
1199
1200 l += s;
1201 }
1202}
1203
1204
1206(
1207 vector& d,
1208 const label facei
1209) const
1210{
1211 if (!parallel())
1212 {
1213 const tensor& T =
1214 (
1215 reverseT().size() == 1
1216 ? reverseT()[0]
1217 : reverseT()[facei]
1219
1220 d = Foam::transform(T, d);
1221 }
1222}
1223
1224
1226(
1227 const primitivePatch& referPatch,
1228 const pointField& thisCtrs,
1229 const vectorField& thisAreas,
1230 const pointField& thisCc,
1231 const pointField& nbrCtrs,
1232 const vectorField& nbrAreas,
1233 const pointField& nbrCc
1234)
1235{}
1236
1237
1240 PstreamBuffers& pBufs,
1241 const primitivePatch& pp
1242) const
1243{}
1244
1245
1247(
1248 PstreamBuffers& pBufs,
1249 const primitivePatch& pp,
1251 labelList& rotation
1252) const
1253{
1254 faceMap.setSize(pp.size());
1255 faceMap = -1;
1256
1257 rotation.setSize(pp.size());
1258 rotation = 0;
1259
1260 return false;
1261}
1262
1263
1265(
1266 const label facei,
1267 const vector& n,
1268 point& p
1269) const
1270{
1271 point prt(p);
1272 reverseTransformPosition(prt, facei);
1273
1274 vector nrt(n);
1275 reverseTransformDirection(nrt, facei);
1276
1277 label nbrFacei = -1;
1278
1279 if (owner())
1280 {
1281 nbrFacei = AMI().tgtPointFace
1282 (
1283 *this,
1284 neighbPatch(),
1285 nrt,
1286 facei,
1287 prt
1288 );
1289 }
1290 else
1291 {
1292 nbrFacei = neighbPatch().AMI().srcPointFace
1293 (
1294 neighbPatch(),
1295 *this,
1296 nrt,
1297 facei,
1298 prt
1299 );
1300 }
1301
1302 if (nbrFacei >= 0)
1303 {
1304 p = prt;
1305 }
1306
1307 return nbrFacei;
1308}
1309
1310
1311void Foam::cyclicAMIPolyPatch::write(Ostream& os) const
1312{
1314 if (!nbrPatchName_.empty())
1315 {
1316 os.writeEntry("neighbourPatch", nbrPatchName_);
1317 }
1318 coupleGroup_.write(os);
1319
1320 switch (transform())
1321 {
1322 case ROTATIONAL:
1323 {
1324 os.writeEntry("rotationAxis", rotationAxis_);
1325 os.writeEntry("rotationCentre", rotationCentre_);
1326
1327 if (rotationAngleDefined_)
1328 {
1329 os.writeEntry("rotationAngle", radToDeg(rotationAngle_));
1330 }
1331
1332 break;
1333 }
1334 case TRANSLATIONAL:
1335 {
1336 os.writeEntry("separationVector", separationVector_);
1337 break;
1338 }
1339 case NOORDERING:
1340 {
1341 break;
1342 }
1343 default:
1344 {
1345 // No additional info to write
1346 }
1347 }
1348
1349 if (!periodicPatchName_.empty())
1350 {
1351 os.writeEntry("periodicPatch", periodicPatchName_);
1352 }
1353
1354 AMIPtr_->write(os);
1355
1356 if (!surfDict_.empty())
1357 {
1358 surfDict_.writeEntry(surfDict_.dictName(), os);
1359 }
1360
1361 if (createAMIFaces_)
1362 {
1363 os.writeEntry("createAMIFaces", createAMIFaces_);
1364 os.writeEntry("srcSize", srcFaceIDs_.size());
1365 os.writeEntry("tgtSize", tgtFaceIDs_.size());
1366 os.writeEntry("moveFaceCentres", moveFaceCentres_);
1367 }
1368
1369 os.writeEntryIfDifferent<scalar>("fraction", Zero, fraction_);
1370}
1371
1372
1373// ************************************************************************* //
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.
if(patchID !=-1)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
void setSize(label n)
Alias for resize().
Definition List.H:536
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
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:346
const Field< point_type > & points() const noexcept
virtual void movePoints(const Field< point_type > &)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
static const SubList< face > & null() noexcept
Definition SubList.H:70
static const Tensor I
Definition Tensor.H:81
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
fileName path() const
The path for the case = rootPath/caseName.
Definition TimePathsI.H:102
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
void size(const label n)
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
static const Form max
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
bool good() const noexcept
The patchGroup has a 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 const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
virtual const boolList & collocated() const
Are faces collocated. Either size 0,1 or length of patch.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
coupledPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform)
Construct from components.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrary Mesh Interface (AMI).
scalar rotationAngle_
Rotation angle.
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.
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.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p following trajectory vector n.
const word & neighbPatchName() const
Neighbour patch name.
virtual bool owner() const
Does this side own the patch?
autoPtr< searchableSurface > surfPtr_
Projection surface.
const scalar fraction_
Particle displacement fraction across AMI.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
vectorField faceCentres0_
Temporary storage for AMI face centres.
virtual void restoreScaledGeometry()
Helper to re-apply the geometric scaling lost during mesh updates.
label periodicPatchID_
Periodic patch.
virtual void clearGeom()
Clear geometry.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
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.
vectorField faceAreas0_
Temporary storage for AMI face areas.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not refer to *this (except for name() and type() etc....
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
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.
virtual void setAMIFaces()
Set properties of newly inserted faces after topological changes.
word periodicPatchName_
Periodic patch name.
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
label nbrPatchID_
Index of other half.
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 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.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
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.
virtual void resetAMI(const UList< point > &points) const
Reset the AMI interpolator, supply patch points.
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.
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
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.
Cyclic plane patch.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
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,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const pointField & oldPoints() const
Return old points (mesh motion).
Definition polyMesh.C:1113
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition polyPatch.C:59
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition polyPatch.H:117
virtual void clearGeom()
Clear geometry.
Definition polyPatch.C:66
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition polyPatch.H:140
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
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:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const volScalarField & T
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)
auto & name
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))
surface1 clear()
#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.
Namespace for handling debugging switches.
Definition debug.C:45
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedScalar sin(const dimensionedScalar &ds)
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.
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr scalar radToDeg() noexcept
Multiplication factor for radians to degrees conversion.
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 ...
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
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.
Vector< scalar > vector
Definition vector.H:57
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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:50
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.