Loading...
Searching...
No Matches
cyclicPolyPatch.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) 2015-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 "cyclicPolyPatch.H"
31#include "polyBoundaryMesh.H"
32#include "polyMesh.H"
33#include "OFstream.H"
34#include "matchPoints.H"
35#include "edgeHashes.H"
36#include "Time.H"
37#include "transformField.H"
38#include "SubField.H"
39#include "unitConversion.H"
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
44{
46
49}
50
51
52// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53
54Foam::label Foam::cyclicPolyPatch::findMaxArea
55(
56 const UList<point>& points,
57 const UList<face>& faces
58)
59{
60 label maxI = -1;
61 scalar maxAreaSqr = -GREAT;
62
63 forAll(faces, facei)
64 {
65 scalar areaSqr = faces[facei].magSqr(points);
66
67 if (maxAreaSqr < areaSqr)
68 {
69 maxAreaSqr = areaSqr;
70 maxI = facei;
71 }
72 }
73 return maxI;
74}
75
76
78{
79 if (size())
80 {
81 // Half0
82 const cyclicPolyPatch& half0 = *this;
83 vectorField half0Areas(half0.size());
84 forAll(half0, facei)
85 {
86 half0Areas[facei] = half0[facei].areaNormal(half0.points());
87 }
88
89 // Half1
90 const cyclicPolyPatch& half1 = neighbPatch();
91 vectorField half1Areas(half1.size());
92 forAll(half1, facei)
93 {
94 half1Areas[facei] = half1[facei].areaNormal(half1.points());
95 }
96
97 calcTransforms
98 (
99 half0,
100 half0.faceCentres(),
101 half0Areas,
102 half1.faceCentres(),
103 half1Areas
104 );
105 }
106}
107
108
110(
111 const primitivePatch& half0,
112 const pointField& half0Ctrs,
113 const vectorField& half0Areas,
114 const pointField& half1Ctrs,
115 const vectorField& half1Areas
116)
117{
118 if (debug && owner())
119 {
120 fileName casePath(boundaryMesh().mesh().time().path());
121 {
122 fileName nm0(casePath/name()+"_faces.obj");
123 Pout<< "cyclicPolyPatch::calcTransforms : Writing " << name()
124 << " faces to OBJ file " << nm0 << endl;
125 writeOBJ(nm0, half0, half0.points());
126 }
127 const cyclicPolyPatch& half1 = neighbPatch();
128 {
129 fileName nm1(casePath/half1.name()+"_faces.obj");
130 Pout<< "cyclicPolyPatch::calcTransforms : Writing " << half1.name()
131 << " faces to OBJ file " << nm1 << endl;
132 writeOBJ(nm1, half1, half1.points());
133 }
134 {
135 OFstream str(casePath/name()+"_to_" + half1.name() + ".obj");
136 label vertI = 0;
137 Pout<< "cyclicPolyPatch::calcTransforms :"
138 << " Writing coupled face centres as lines to " << str.name()
139 << endl;
140 forAll(half0Ctrs, i)
141 {
142 const point& p0 = half0Ctrs[i];
143 str << "v " << p0.x() << ' ' << p0.y() << ' ' << p0.z() << nl;
144 vertI++;
145 const point& p1 = half1Ctrs[i];
146 str << "v " << p1.x() << ' ' << p1.y() << ' ' << p1.z() << nl;
147 vertI++;
148 str << "l " << vertI-1 << ' ' << vertI << nl;
149 }
150 }
151 }
152
153
154 // Some sanity checks
155
156 if (half0Ctrs.size() != half1Ctrs.size())
157 {
159 << "For patch " << name()
160 << " there are " << half0Ctrs.size()
161 << " face centres, for the neighbour patch " << neighbPatch().name()
162 << " there are " << half1Ctrs.size()
163 << exit(FatalError);
164 }
165
166 if (transform() != neighbPatch().transform())
167 {
169 << "Patch " << name()
170 << " has transform type " << transformTypeNames[transform()]
171 << ", neighbour patch " << neighbPatchName()
172 << " has transform type "
173 << neighbPatch().transformTypeNames[neighbPatch().transform()]
174 << exit(FatalError);
175 }
176
177
178 // Calculate transformation tensors
179
180 if (half0Ctrs.size() > 0)
181 {
182 vectorField half0Normals(half0Areas.size());
183 vectorField half1Normals(half1Areas.size());
184
185 scalar maxAreaDiff = -GREAT;
186 label maxAreaFacei = -1;
187
188 forAll(half0, facei)
189 {
190 scalar magSf = mag(half0Areas[facei]);
191 scalar nbrMagSf = mag(half1Areas[facei]);
192 scalar avSf = (magSf + nbrMagSf)/2.0;
193
194 if (magSf < ROOTVSMALL && nbrMagSf < ROOTVSMALL)
195 {
196 // Undetermined normal. Use dummy normal to force separation
197 // check. (note use of sqrt(VSMALL) since that is how mag
198 // scales)
199 half0Normals[facei] = point(1, 0, 0);
200 half1Normals[facei] = half0Normals[facei];
201 }
202 else
203 {
204 scalar areaDiff = mag(magSf - nbrMagSf)/avSf;
205
206 if (areaDiff > maxAreaDiff)
207 {
208 maxAreaDiff = areaDiff;
209 maxAreaFacei = facei;
210 }
211
212 if (areaDiff > matchTolerance())
213 {
215 << "face " << facei
216 << " area does not match neighbour by "
217 << 100*areaDiff
218 << "% -- possible face ordering problem." << endl
219 << "patch:" << name()
220 << " my area:" << magSf
221 << " neighbour area:" << nbrMagSf
222 << " matching tolerance:" << matchTolerance()
223 << endl
224 << "Mesh face:" << start()+facei
225 << " fc:" << half0Ctrs[facei]
226 << endl
227 << "Neighbour fc:" << half1Ctrs[facei]
228 << endl
229 << "If you are certain your matching is correct"
230 << " you can increase the 'matchTolerance' setting"
231 << " in the patch dictionary in the boundary file."
232 << endl
233 << "Rerun with cyclic debug flag set"
234 << " for more information." << exit(FatalError);
235 }
236 else
237 {
238 half0Normals[facei] = half0Areas[facei] / magSf;
239 half1Normals[facei] = half1Areas[facei] / nbrMagSf;
240 }
241 }
242 }
243
244
245 // Print area match
246 if (debug)
247 {
248 Pout<< "cyclicPolyPatch::calcTransforms :"
249 << " patch:" << name()
250 << " Max area error:" << 100*maxAreaDiff << "% at face:"
251 << maxAreaFacei << " at:" << half0Ctrs[maxAreaFacei]
252 << " coupled face at:" << half1Ctrs[maxAreaFacei]
253 << endl;
254 }
255
256
257 // Calculate transformation tensors
258
259 if (transform() == ROTATIONAL)
260 {
261 // Calculate using the given rotation axis and centre. Do not
262 // use calculated normals.
263 vector n0 = findFaceMaxRadius(half0Ctrs);
264 vector n1 = -findFaceMaxRadius(half1Ctrs);
265 n0.normalise();
266 n1.normalise();
267
268 if (debug)
269 {
270 Pout<< "cyclicPolyPatch::calcTransforms :"
271 << " patch:" << name()
272 << " Specified rotation :"
273 << " n0:" << n0 << " n1:" << n1
274 << " swept angle: " << radToDeg(acos(n0 & n1))
275 << " [deg]" << endl;
276 }
277
278 // Extended tensor from two local coordinate systems calculated
279 // using normal and rotation axis
280 const tensor E0
281 (
282 rotationAxis_,
283 (n0 ^ rotationAxis_),
284 n0
285 );
286 const tensor E1
287 (
288 rotationAxis_,
289 (-n1 ^ rotationAxis_),
290 -n1
291 );
292 const tensor revT(E1.T() & E0);
293
294 const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
295 const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
296 const_cast<vectorField&>(separation()).clear();
297 const_cast<boolList&>(collocated()) = boolList(1, false);
298 }
299 else
300 {
301 scalarField half0Tols
302 (
303 matchTolerance()
304 *calcFaceTol
305 (
306 half0,
307 half0.points(),
308 static_cast<const pointField&>(half0Ctrs)
309 )
310 );
311
312 calcTransformTensors
313 (
314 static_cast<const pointField&>(half0Ctrs),
315 static_cast<const pointField&>(half1Ctrs),
316 half0Normals,
317 half1Normals,
318 half0Tols,
319 matchTolerance(),
320 transform()
321 );
322
323
324 if (transform() == TRANSLATIONAL)
325 {
326 if (debug)
327 {
328 Pout<< "cyclicPolyPatch::calcTransforms :"
329 << " patch:" << name()
330 << " Specified separation vector : "
331 << separationVector_ << endl;
332 }
333
334 // Check that separation vectors are same.
335 const scalar avgTol = average(half0Tols);
336 if
337 (
338 mag(separationVector_ + neighbPatch().separationVector_)
339 > avgTol
340 )
341 {
343 << "Specified separation vector " << separationVector_
344 << " differs by that of neighbouring patch "
345 << neighbPatch().separationVector_
346 << " by more than tolerance " << avgTol << endl
347 << "patch:" << name()
348 << " neighbour:" << neighbPatchName()
349 << endl;
350 }
351
352
353 // Override computed transform with specified.
354 if
355 (
356 separation().size() != 1
357 || mag(separation()[0] - separationVector_) > avgTol
358 )
359 {
361 << "Specified separationVector " << separationVector_
362 << " differs from computed separation vector "
363 << separation() << endl
364 << "This probably means your geometry is not consistent"
365 << " with the specified separation and might lead"
366 << " to problems." << endl
367 << "Continuing with specified separation vector "
368 << separationVector_ << endl
369 << "patch:" << name()
370 << " neighbour:" << neighbPatchName()
371 << endl;
372 }
373
374 // Set tensors
375 const_cast<tensorField&>(forwardT()).clear();
376 const_cast<tensorField&>(reverseT()).clear();
377 const_cast<vectorField&>(separation()) = vectorField
378 (
379 1,
380 separationVector_
381 );
382 const_cast<boolList&>(collocated()) = boolList(1, false);
383 }
384 }
385 }
386}
387
388
389void Foam::cyclicPolyPatch::getCentresAndAnchors
390(
391 const primitivePatch& pp0,
392 const primitivePatch& pp1,
393
394 pointField& half0Ctrs,
395 pointField& half1Ctrs,
396 pointField& anchors0,
397 scalarField& tols
398) const
399{
400 // Get geometric data on both halves.
401 half0Ctrs = pp0.faceCentres();
402 anchors0 = getAnchorPoints(pp0, pp0.points(), transform());
403 half1Ctrs = pp1.faceCentres();
404
405 if (debug)
406 {
407 Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
408 << " patch:" << name() << nl
409 << "half0 untransformed faceCentres (avg) : "
410 << gAverage(half0Ctrs) << nl
411 << "half1 untransformed faceCentres (avg) : "
412 << gAverage(half1Ctrs) << endl;
413 }
414
415 if (half0Ctrs.size())
416 {
417 switch (transform())
418 {
419 case ROTATIONAL:
420 {
421 vector n0 = findFaceMaxRadius(half0Ctrs);
422 vector n1 = -findFaceMaxRadius(half1Ctrs);
423 n0.normalise();
424 n1.normalise();
425
426 if (debug)
427 {
428 Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
429 << " patch:" << name()
430 << " Specified rotation :"
431 << " n0:" << n0 << " n1:" << n1
432 << " swept angle: "
433 << radToDeg(acos(n0 & n1)) << " [deg]"
434 << endl;
435 }
436
437 // Extended tensor from two local coordinate systems calculated
438 // using normal and rotation axis
439 const tensor E0
440 (
441 rotationAxis_,
442 (n0 ^ rotationAxis_),
443 n0
444 );
445 const tensor E1
446 (
447 rotationAxis_,
448 (-n1 ^ rotationAxis_),
449 -n1
450 );
451 const tensor revT(E1.T() & E0);
452
453 // Rotation
454 forAll(half0Ctrs, facei)
455 {
456 half0Ctrs[facei] =
458 (
459 revT,
460 half0Ctrs[facei] - rotationCentre_
461 )
462 + rotationCentre_;
463 anchors0[facei] =
465 (
466 revT,
467 anchors0[facei] - rotationCentre_
468 )
469 + rotationCentre_;
470 }
471
472 break;
473 }
474 case TRANSLATIONAL:
475 {
476 // Transform 0 points.
477
478 if (debug)
479 {
480 Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
481 << " patch:" << name()
482 << "Specified translation : " << separationVector_
483 << endl;
484 }
485
486 // Note: getCentresAndAnchors gets called on the neighbour side
487 // so separationVector is owner-neighbour points.
488
489 half0Ctrs -= separationVector_;
490 anchors0 -= separationVector_;
491 break;
492 }
493 default:
494 {
495 // Assumes that cyclic is rotational. This is also the initial
496 // condition for patches without faces.
497
498 // Determine the face with max area on both halves. These
499 // two faces are used to determine the transformation tensors
500 label max0I = findMaxArea(pp0.points(), pp0);
501 const vector n0 = pp0[max0I].unitNormal(pp0.points());
502
503 label max1I = findMaxArea(pp1.points(), pp1);
504 const vector n1 = pp1[max1I].unitNormal(pp1.points());
505
506 if (mag(n0 & n1) < 1-matchTolerance())
507 {
508 if (debug)
509 {
510 Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
511 << " patch:" << name()
512 << " Detected rotation :"
513 << " n0:" << n0 << " n1:" << n1 << endl;
514 }
515
516 // Rotation (around origin)
517 const tensor revT(rotationTensor(n0, -n1));
518
519 // Rotation
520 forAll(half0Ctrs, facei)
521 {
522 half0Ctrs[facei] = Foam::transform
523 (
524 revT,
525 half0Ctrs[facei]
526 );
527 anchors0[facei] = Foam::transform
528 (
529 revT,
530 anchors0[facei]
531 );
532 }
533 }
534 else
535 {
536 // Parallel translation. Get average of all used points.
537
538 const point ctr0(sum(pp0.localPoints())/pp0.nPoints());
539 const point ctr1(sum(pp1.localPoints())/pp1.nPoints());
540
541 if (debug)
542 {
543 Pout<< "cyclicPolyPatch::getCentresAndAnchors :"
544 << " patch:" << name()
545 << " Detected translation :"
546 << " n0:" << n0 << " n1:" << n1
547 << " ctr0:" << ctr0 << " ctr1:" << ctr1 << endl;
548 }
549
550 half0Ctrs += ctr1 - ctr0;
551 anchors0 += ctr1 - ctr0;
552 }
553 break;
554 }
555 }
556 }
557
558 // Calculate typical distance per face
559 tols = matchTolerance()*calcFaceTol(pp1, pp1.points(), half1Ctrs);
560}
561
562
563Foam::vector Foam::cyclicPolyPatch::findFaceMaxRadius
564(
565 const pointField& faceCentres
566) const
567{
568 // Determine a face furthest away from the axis
569
570 const vectorField n((faceCentres - rotationCentre_) ^ rotationAxis_);
571
572 const scalarField magRadSqr(magSqr(n));
573
574 label facei = findMax(magRadSqr);
575
576 if (debug)
577 {
578 Info<< "findFaceMaxRadius(const pointField&) : patch: " << name() << nl
579 << " rotFace = " << facei << nl
580 << " point = " << faceCentres[facei] << nl
581 << " distance = " << Foam::sqrt(magRadSqr[facei])
582 << endl;
583 }
585 return n[facei];
586}
587
588
589// * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
590
592(
593 const word& name,
594 const label size,
595 const label start,
596 const label index,
597 const polyBoundaryMesh& bm,
598 const word& patchType,
599 const transformType transform
600)
601:
602 coupledPolyPatch(name, size, start, index, bm, patchType, transform),
603 neighbPatchName_(),
604 neighbPatchID_(-1),
605 rotationAxis_(Zero),
606 rotationCentre_(Zero),
607 separationVector_(Zero),
608 coupledPointsPtr_(nullptr),
609 coupledEdgesPtr_(nullptr)
610{
611 // Neighbour patch might not be valid yet so no transformation
612 // calculation possible.
613}
614
615
617(
618 const word& name,
619 const label size,
620 const label start,
621 const label index,
622 const polyBoundaryMesh& bm,
623 const word& neighbPatchName,
624 const transformType transform,
625 const vector& rotationAxis,
626 const point& rotationCentre,
627 const vector& separationVector
628)
629:
630 coupledPolyPatch(name, size, start, index, bm, typeName, transform),
631 neighbPatchName_(neighbPatchName),
632 neighbPatchID_(-1),
633 rotationAxis_(rotationAxis),
634 rotationCentre_(rotationCentre),
635 separationVector_(separationVector),
636 coupledPointsPtr_(nullptr),
637 coupledEdgesPtr_(nullptr)
638{
639 // Neighbour patch might not be valid yet so no transformation
640 // calculation possible.
641}
642
643
645(
646 const word& name,
647 const dictionary& dict,
648 const label index,
649 const polyBoundaryMesh& bm,
650 const word& patchType
651)
652:
653 coupledPolyPatch(name, dict, index, bm, patchType),
654 neighbPatchName_(dict.getOrDefault("neighbourPatch", word::null)),
655 coupleGroup_(dict),
656 neighbPatchID_(-1),
657 rotationAxis_(Zero),
658 rotationCentre_(Zero),
659 separationVector_(Zero),
660 coupledPointsPtr_(nullptr),
661 coupledEdgesPtr_(nullptr)
662{
663 if (neighbPatchName_.empty() && !coupleGroup_.good())
664 {
665 FatalIOErrorInFunction(dict)
666 << "No \"neighbourPatch\" provided." << endl
667 << "Is your mesh uptodate with split cyclics?" << endl
668 << "Run foamUpgradeCyclics to convert mesh and fields"
669 << " to split cyclics." << exit(FatalIOError);
670 }
671
672 if (neighbPatchName_ == name)
673 {
675 << "Neighbour patch name " << neighbPatchName_
676 << " cannot be the same as this patch " << name
677 << exit(FatalIOError);
678 }
679
680 switch (transform())
681 {
682 case ROTATIONAL:
683 {
684 dict.readEntry("rotationAxis", rotationAxis_);
685 dict.readEntry("rotationCentre", rotationCentre_);
686
687 scalar magRot = mag(rotationAxis_);
688 if (magRot < SMALL)
689 {
691 << "Illegal rotationAxis " << rotationAxis_ << endl
692 << "Please supply a non-zero vector."
693 << exit(FatalIOError);
694 }
695 rotationAxis_ /= magRot;
696
697 break;
698 }
699 case TRANSLATIONAL:
700 {
701 dict.readEntry("separationVector", separationVector_);
702 break;
703 }
704 default:
705 {
706 // no additional info required
707 }
709
710 // Neighbour patch might not be valid yet so no transformation
711 // calculation possible.
712}
713
714
716(
717 const cyclicPolyPatch& pp,
718 const polyBoundaryMesh& bm
719)
720:
721 coupledPolyPatch(pp, bm),
722 neighbPatchName_(pp.neighbPatchName_),
723 coupleGroup_(pp.coupleGroup_),
724 neighbPatchID_(-1),
725 rotationAxis_(pp.rotationAxis_),
726 rotationCentre_(pp.rotationCentre_),
727 separationVector_(pp.separationVector_),
728 coupledPointsPtr_(nullptr),
729 coupledEdgesPtr_(nullptr)
730{
731 // Neighbour patch might not be valid yet so no transformation
732 // calculation possible.
733}
734
735
737(
738 const cyclicPolyPatch& pp,
739 const label nrbPatchID,
740 const labelList& faceCells
741)
742:
744 neighbPatchName_(pp.neighbPatchName_),
745 coupleGroup_(pp.coupleGroup_),
746 neighbPatchID_(nrbPatchID),
747 rotationAxis_(pp.rotationAxis_),
748 rotationCentre_(pp.rotationCentre_),
749 separationVector_(pp.separationVector_),
750 coupledPointsPtr_(nullptr),
751 coupledEdgesPtr_(nullptr)
752{
753 // Neighbour patch might not be valid yet so no transformation
754 // calculation possible.
755}
756
757
759(
760 const cyclicPolyPatch& pp,
761 const polyBoundaryMesh& bm,
762 const label index,
763 const label newSize,
764 const label newStart,
765 const word& neighbName
766)
767:
768 coupledPolyPatch(pp, bm, index, newSize, newStart),
769 neighbPatchName_(neighbName),
770 coupleGroup_(pp.coupleGroup_),
771 neighbPatchID_(-1),
772 rotationAxis_(pp.rotationAxis_),
773 rotationCentre_(pp.rotationCentre_),
774 separationVector_(pp.separationVector_),
775 coupledPointsPtr_(nullptr),
776 coupledEdgesPtr_(nullptr)
777{
778 if (neighbName == name())
779 {
781 << "Neighbour patch name " << neighbName
782 << " cannot be the same as this patch " << name()
783 << exit(FatalError);
785
786 // Neighbour patch might not be valid yet so no transformation
787 // calculation possible.
788}
789
790
792(
793 const cyclicPolyPatch& pp,
794 const polyBoundaryMesh& bm,
795 const label index,
796 const labelUList& mapAddressing,
797 const label newStart
798)
799:
800 coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
801 neighbPatchName_(pp.neighbPatchName_),
802 coupleGroup_(pp.coupleGroup_),
803 neighbPatchID_(-1),
804 rotationAxis_(pp.rotationAxis_),
805 rotationCentre_(pp.rotationCentre_),
806 separationVector_(pp.separationVector_),
807 coupledPointsPtr_(nullptr),
808 coupledEdgesPtr_(nullptr)
809{}
810
811
812// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
815{}
816
817
818// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
819
821{
822 if (neighbPatchName_.empty())
823 {
824 // Try and use patchGroup to find samplePatch and sampleRegion
825 label patchID = coupleGroup_.findOtherPatchID(*this);
827 neighbPatchName_ = boundaryMesh()[patchID].name();
828 }
829 return neighbPatchName_;
830}
831
832
833Foam::label Foam::cyclicPolyPatch::neighbPatchID() const
834{
835 if (neighbPatchID_ == -1)
836 {
837 neighbPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
838
839 if (neighbPatchID_ == -1)
840 {
842 << "Illegal neighbourPatch name " << neighbPatchName()
843 << endl << "Valid patch names are "
844 << this->boundaryMesh().names()
845 << exit(FatalError);
846 }
847
848 // Check that it is a cyclic
850 (
851 this->boundaryMesh()[neighbPatchID_]
852 );
853
854 if (nbrPatch.neighbPatchName() != name())
855 {
857 << "Patch " << name()
858 << " specifies neighbour patch " << neighbPatchName()
859 << endl << " but that in return specifies "
860 << nbrPatch.neighbPatchName()
862 }
863 }
864 return neighbPatchID_;
865}
866
867
869{
870 if (!parallel())
871 {
872 if (transform() == ROTATIONAL)
873 {
874 l =
875 Foam::transform(forwardT(), l-rotationCentre_)
876 + rotationCentre_;
877 }
878 else
879 {
880 l = Foam::transform(forwardT(), l);
881 }
882 }
883 else if (separated())
884 {
885 // transformPosition gets called on the receiving side,
886 // separation gets calculated on the sending side so subtract.
887
888 const vectorField& s = separation();
889 if (s.size() == 1)
890 {
891 forAll(l, i)
892 {
893 l[i] -= s[0];
894 }
895 }
896 else
898 l -= s;
899 }
900 }
901}
902
903
904void Foam::cyclicPolyPatch::transformPosition(point& l, const label facei) const
905{
906 if (!parallel())
907 {
908 const tensor& T =
909 (
910 forwardT().size() == 1
911 ? forwardT()[0]
912 : forwardT()[facei]
913 );
914
915 if (transform() == ROTATIONAL)
916 {
917 l = Foam::transform(T, l-rotationCentre_) + rotationCentre_;
918 }
919 else
920 {
921 l = Foam::transform(T, l);
922 }
923 }
924 else if (separated())
925 {
926 const vector& s =
927 (
928 separation().size() == 1
929 ? separation()[0]
930 : separation()[facei]
931 );
932
933 l -= s;
934 }
935}
936
945(
946 const primitivePatch& referPatch,
947 pointField& nbrCtrs,
948 vectorField& nbrAreas,
949 pointField& nbrCc
950)
951{}
952
953
955(
956 const primitivePatch& referPatch,
957 const pointField& thisCtrs,
958 const vectorField& thisAreas,
959 const pointField& thisCc,
960 const pointField& nbrCtrs,
961 const vectorField& nbrAreas,
962 const pointField& nbrCc
963)
964{
965 calcTransforms
966 (
967 referPatch,
968 thisCtrs,
969 thisAreas,
970 nbrCtrs,
971 nbrAreas
972 );
973}
974
975
977{
978 calcGeometry
979 (
980 *this,
981 faceCentres(),
982 faceAreas(),
983 faceCellCentres(),
992(
993 PstreamBuffers& pBufs,
995)
996{
998}
999
1000
1002(
1003 PstreamBuffers& pBufs,
1004 const pointField& p
1006{
1007 polyPatch::movePoints(pBufs, p);
1008 calcTransforms();
1009}
1010
1017
1020 polyPatch::updateMesh(pBufs);
1021 coupledPointsPtr_.reset(nullptr);
1022 coupledEdgesPtr_.reset(nullptr);
1023}
1024
1025
1027{
1028 if (!coupledPointsPtr_)
1029 {
1030 const faceList& nbrLocalFaces = neighbPatch().localFaces();
1031 const labelList& nbrMeshPoints = neighbPatch().meshPoints();
1032
1033 // Now all we know is that relative face index in *this is same
1034 // as coupled face in nbrPatch and also that the 0th vertex
1035 // corresponds.
1036
1037 // From local point to nbrPatch or -1.
1038 labelList coupledPoint(nPoints(), -1);
1039
1040 forAll(*this, patchFacei)
1041 {
1042 const face& fA = localFaces()[patchFacei];
1043 const face& fB = nbrLocalFaces[patchFacei];
1044
1045 forAll(fA, indexA)
1046 {
1047 label patchPointA = fA[indexA];
1048
1049 if (coupledPoint[patchPointA] == -1)
1050 {
1051 label indexB = (fB.size() - indexA) % fB.size();
1052
1053 // Filter out points on wedge axis
1054 if (meshPoints()[patchPointA] != nbrMeshPoints[fB[indexB]])
1055 {
1056 coupledPoint[patchPointA] = fB[indexB];
1057 }
1058 }
1059 }
1060 }
1061
1062 coupledPointsPtr_.reset(new edgeList(nPoints()));
1063 auto& connected = *coupledPointsPtr_;
1064
1065 // Extract coupled points.
1066 label connectedI = 0;
1067
1068 forAll(coupledPoint, i)
1069 {
1070 if (coupledPoint[i] != -1)
1071 {
1072 connected[connectedI++] = edge(i, coupledPoint[i]);
1073 }
1074 }
1075
1076 connected.setSize(connectedI);
1077
1078 if (debug)
1079 {
1080 OFstream str
1081 (
1082 boundaryMesh().mesh().time().path()
1083 /name() + "_coupledPoints.obj"
1084 );
1085 label vertI = 0;
1086
1087 Pout<< "Writing file " << str.name() << " with coordinates of "
1088 << "coupled points" << endl;
1089
1090 forAll(connected, i)
1091 {
1092 const point& a = points()[meshPoints()[connected[i][0]]];
1093 const point& b = points()[nbrMeshPoints[connected[i][1]]];
1094
1095 str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
1096 str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
1097 vertI += 2;
1098
1099 str<< "l " << vertI-1 << ' ' << vertI << nl;
1101 }
1102 }
1103 return *coupledPointsPtr_;
1104}
1105
1106
1108{
1109 if (!coupledEdgesPtr_)
1110 {
1111 const edgeList& pointCouples = coupledPoints();
1112
1113 // Build map from points on *this (A) to points on neighbourpatch (B)
1114 Map<label> aToB(2*pointCouples.size());
1115
1116 for (const edge& e : pointCouples)
1117 {
1118 aToB.insert(e[0], e[1]);
1119 }
1120
1121 // Map from edge on A to points (in B indices)
1122 EdgeMap<label> edgeMap(nEdges());
1123
1124 forAll(*this, patchFacei)
1125 {
1126 const labelList& fEdges = faceEdges()[patchFacei];
1127
1128 for (const label edgeI : fEdges)
1129 {
1130 const edge& e = edges()[edgeI];
1131
1132 // Convert edge end points to corresponding points on B side.
1133 const auto fnd0 = aToB.cfind(e[0]);
1134 if (fnd0.good())
1135 {
1136 const auto fnd1 = aToB.cfind(e[1]);
1137 if (fnd1.good())
1138 {
1139 edgeMap.insert(edge(fnd0(), fnd1()), edgeI);
1140 }
1141 }
1142 }
1143 }
1144
1145 // Use the edgeMap to get the edges on the B side.
1146
1147 const cyclicPolyPatch& neighbPatch = this->neighbPatch();
1148 const labelList& nbrMp = neighbPatch.meshPoints();
1149 const labelList& mp = meshPoints();
1150
1151
1152 coupledEdgesPtr_.reset(new edgeList(edgeMap.size()));
1153 auto& coupledEdges = *coupledEdgesPtr_;
1154 label coupleI = 0;
1155
1156 forAll(neighbPatch, patchFacei)
1157 {
1158 const labelList& fEdges = neighbPatch.faceEdges()[patchFacei];
1159
1160 for (const label edgeI : fEdges)
1161 {
1162 const edge& e = neighbPatch.edges()[edgeI];
1163
1164 // Look up A edge from HashTable.
1165 auto iter = edgeMap.find(e);
1166
1167 if (iter.good())
1168 {
1169 const label edgeA = iter.val();
1170 const edge& eA = edges()[edgeA];
1171
1172 // Store correspondence. Filter out edges on wedge axis.
1173 if
1174 (
1175 edge(mp[eA[0]], mp[eA[1]])
1176 != edge(nbrMp[e[0]], nbrMp[e[1]])
1177 )
1178 {
1179 coupledEdges[coupleI++] = edge(edgeA, edgeI);
1180 }
1181
1182 // Remove so we build unique list
1183 edgeMap.erase(iter);
1184 }
1185 }
1186 }
1187 coupledEdges.setSize(coupleI);
1188
1189
1190 // Some checks
1191
1192 forAll(coupledEdges, i)
1193 {
1194 const edge& e = coupledEdges[i];
1195
1196 if (e[0] < 0 || e[1] < 0)
1197 {
1199 << "Problem : at position " << i
1200 << " illegal couple:" << e
1201 << abort(FatalError);
1202 }
1203 }
1204
1205 if (debug)
1206 {
1207 OFstream str
1208 (
1209 boundaryMesh().mesh().time().path()
1210 /name() + "_coupledEdges.obj"
1211 );
1212 label vertI = 0;
1213
1214 Pout<< "Writing file " << str.name() << " with centres of "
1215 << "coupled edges" << endl;
1216
1217 for (const edge& e : coupledEdges)
1218 {
1219 const point& a = edges()[e[0]].centre(localPoints());
1220 const point& b = neighbPatch.edges()[e[1]].centre
1221 (
1222 neighbPatch.localPoints()
1223 );
1224
1225 str<< "v " << a.x() << ' ' << a.y() << ' ' << a.z() << nl;
1226 str<< "v " << b.x() << ' ' << b.y() << ' ' << b.z() << nl;
1227 vertI += 2;
1228
1229 str<< "l " << vertI-1 << ' ' << vertI << nl;
1231 }
1232 }
1233 return *coupledEdgesPtr_;
1234}
1235
1236
1238(
1240 const primitivePatch& pp
1241) const
1242{
1243 if (owner())
1244 {
1245 // Save patch for use in non-owner side ordering. Equivalent to
1246 // processorPolyPatch using OPstream.
1247 ownerPatchPtr_.reset
1248 (
1249 new primitivePatch
1250 (
1251 pp,
1252 pp.points()
1253 )
1254 );
1255 }
1256}
1257
1258
1260(
1261 PstreamBuffers& pBufs,
1262 const primitivePatch& pp,
1264 labelList& rotation
1265) const
1266{
1267 if (debug)
1268 {
1269 Pout<< "order : of " << pp.size()
1270 << " faces of patch:" << name()
1271 << " neighbour:" << neighbPatchName()
1272 << endl;
1273 }
1274 faceMap.resize_nocopy(pp.size());
1275 faceMap = -1;
1276
1277 rotation.resize_nocopy(pp.size());
1278 rotation = 0;
1279
1280 if (transform() == NOORDERING)
1281 {
1282 // No faces, nothing to change.
1283 return false;
1284 }
1285
1286 if (owner())
1287 {
1288 // Do nothing (i.e. identical mapping, zero rotation).
1289 // See comment at top.
1290 forAll(faceMap, patchFacei)
1291 {
1292 faceMap[patchFacei] = patchFacei;
1293 }
1294
1295 return false;
1296 }
1297 else
1298 {
1299 // Get stored geometry from initOrder invocation of owner.
1300 const primitivePatch& pp0 = neighbPatch().ownerPatchPtr_();
1301
1302 // Get geometric quantities
1303 pointField half0Ctrs, half1Ctrs, anchors0;
1304 scalarField tols;
1305 getCentresAndAnchors
1306 (
1307 pp0,
1308 pp,
1309
1310 half0Ctrs,
1311 half1Ctrs,
1312 anchors0,
1313 tols
1314 );
1315
1316 if (debug)
1317 {
1318 Pout<< "half0 transformed faceCentres (avg) : "
1319 << gAverage(half0Ctrs) << nl
1320 << "half1 untransformed faceCentres (avg) : "
1321 << gAverage(half1Ctrs) << endl;
1322 }
1323
1324 // Geometric match of face centre vectors
1325 bool matchedAll = matchPoints
1326 (
1327 half1Ctrs,
1328 half0Ctrs,
1329 tols,
1330 true,
1331 faceMap
1332 );
1333
1334 if (!matchedAll || debug)
1335 {
1336 // Dump halves
1337 fileName nm0
1338 (
1339 boundaryMesh().mesh().time().path()
1340 /neighbPatch().name()+"_faces.obj"
1341 );
1342 Pout<< "cyclicPolyPatch::order : Writing neighbour"
1343 << " faces to OBJ file " << nm0 << endl;
1344 writeOBJ(nm0, pp0, pp0.points());
1345
1346 fileName nm1
1347 (
1348 boundaryMesh().mesh().time().path()
1349 /name()+"_faces.obj"
1350 );
1351 Pout<< "cyclicPolyPatch::order : Writing my"
1352 << " faces to OBJ file " << nm1 << endl;
1353 writeOBJ(nm1, pp, pp.points());
1354
1355 OFstream ccStr
1356 (
1357 boundaryMesh().mesh().time().path()
1358 /name() + "_faceCentres.obj"
1359 );
1360 Pout<< "cyclicPolyPatch::order : "
1361 << "Dumping currently found cyclic match as lines between"
1362 << " corresponding face centres to file " << ccStr.name()
1363 << endl;
1364
1365 // Recalculate untransformed face centres
1366 //pointField rawHalf0Ctrs =
1367 // calcFaceCentres(half0Faces, pp.points());
1368 label vertI = 0;
1369
1370 forAll(half1Ctrs, i)
1371 {
1372 if (faceMap[i] != -1)
1373 {
1374 // Write edge between c1 and c0
1375 const point& c0 = half0Ctrs[faceMap[i]];
1376 const point& c1 = half1Ctrs[i];
1377 writeOBJ(ccStr, c0, c1, vertI);
1378 }
1379 }
1380 }
1381
1382 if (!matchedAll)
1383 {
1385 << "Patch:" << name() << " : "
1386 << "Cannot match vectors to faces on both sides of patch"
1387 << endl
1388 << " Perhaps your faces do not match?"
1389 << " The obj files written contain the current match." << endl
1390 << " Continuing with incorrect face ordering from now on!"
1391 << endl;
1392
1393 return false;
1394 }
1395
1396
1397 // Set rotation.
1398 forAll(faceMap, oldFacei)
1399 {
1400 // The face f will be at newFacei (after morphing) and we want its
1401 // anchorPoint (= f[0]) to align with the anchorpoint for the
1402 // corresponding face on the other side.
1403
1404 label newFacei = faceMap[oldFacei];
1405
1406 const point& wantedAnchor = anchors0[newFacei];
1407
1408 rotation[newFacei] = getRotation
1409 (
1410 pp.points(),
1411 pp[oldFacei],
1412 wantedAnchor,
1413 tols[oldFacei]
1414 );
1415
1416 if (rotation[newFacei] == -1)
1417 {
1419 << "in patch " << name()
1420 << " : "
1421 << "Cannot find point on face " << pp[oldFacei]
1422 << " with vertices "
1423 << UIndirectList<point>(pp.points(), pp[oldFacei])
1424 << " that matches point " << wantedAnchor
1425 << " when matching the halves of processor patch " << name()
1426 << "Continuing with incorrect face ordering from now on!"
1427 << endl;
1428
1429 return false;
1430 }
1431 }
1432
1433 ownerPatchPtr_.reset(nullptr);
1434
1435 // Return false if no change necessary, true otherwise.
1436
1437 forAll(faceMap, facei)
1438 {
1439 if (faceMap[facei] != facei || rotation[facei] != 0)
1440 {
1441 return true;
1442 }
1444
1445 return false;
1446 }
1447}
1448
1449
1451{
1453 if (!neighbPatchName_.empty())
1454 {
1455 os.writeEntry("neighbourPatch", neighbPatchName_);
1456 }
1457 coupleGroup_.write(os);
1458 switch (transform())
1459 {
1460 case ROTATIONAL:
1461 {
1462 os.writeEntry("rotationAxis", rotationAxis_);
1463 os.writeEntry("rotationCentre", rotationCentre_);
1464 break;
1465 }
1466 case TRANSLATIONAL:
1467 {
1468 os.writeEntry("separationVector", separationVector_);
1469 break;
1470 }
1471 case NOORDERING:
1472 {
1473 break;
1474 }
1475 default:
1476 {
1477 // no additional info to write
1478 }
1479 }
1480}
1481
1482
1483// ************************************************************************* //
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())
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
Foam::label erase(InputIter first, InputIter last)
Definition HashTable.C:515
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
label size() const noexcept
The number of elements in the list.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void setSize(label n)
Alias for resize().
Definition List.H:536
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
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
const Field< point_type > & points() const noexcept
Return reference to global points.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static const SubList< face > & null() noexcept
Definition SubList.H:70
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
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 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 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.
static label getRotation(const pointField &points, const face &f, const point &anchor, const scalar tol)
Get the number of vertices face f needs to be rotated such that.
virtual const tensorField & forwardT() const
Return face transformation tensor.
static void writeOBJ(Ostream &os, const point &pt)
Write point in OBJ format.
Cyclic plane patch.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
const word & neighbPatchName() const
Neighbour patch name.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
virtual bool owner() const
Does this side own the patch ?
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
cyclicPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from components.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local).
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
const vector & separationVector() const noexcept
Translation vector for translational cyclics.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
virtual void movePoints(PstreamBuffers &, const pointField &)
Correct patches after moving points.
const cyclicPolyPatch & neighbPatch() const
virtual void transformPosition(pointField &l) const
Transform a patch-based position from other side to this side.
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
virtual ~cyclicPolyPatch()
Destructor.
const point & rotationCentre() const noexcept
Point on axis of rotation for rotational cyclics.
virtual void calcTransforms()
Recalculate the transformation tensors.
const vector & rotationAxis() const noexcept
Axis of rotation for rotational cyclics.
virtual label neighbPatchID() const
Neighbour patchID.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
label erase(const label vertex)
Remove an existing vertex from the edge and set its location to '-1'. A negative vertex label never r...
Definition edgeI.H:320
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
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,...
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
const vectorField::subField faceAreas() const
Return face normals.
Definition polyPatch.C:326
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition polyPatch.H:117
const faceList::subList faces() const
Return mesh faces for the patch.
Definition polyPatch.C:301
const vectorField::subField faceCentres() const
Return face centres.
Definition polyPatch.C:320
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition polyPatch.C:53
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
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition polyPatch.C:332
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition polyPatch.H:129
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
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
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
label nPoints
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))
Determine correspondence between points. See below.
surface1 clear()
#define WarningInFunction
Report a warning using Foam::Warning.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
const dimensionedScalar mp
Proton mass.
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
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.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition transform.H:47
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.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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.
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
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
Determine correspondence between pointFields. Gets passed.
Definition matchPoints.C:29
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
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.
Unit conversion functions.