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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "cyclicAMIPolyPatch.H"
27 #include "SubField.H"
28 #include "Time.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(cyclicAMIPolyPatch, 0);
36 
37  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, word);
38  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, dictionary);
39 }
40 
41 
42 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
43 
44 Foam::vector Foam::cyclicAMIPolyPatch::findFaceNormalMaxRadius
45 (
46  const pointField& faceCentres
47 ) const
48 {
49  // Determine a face furthest away from the axis
50 
51  const vectorField n((faceCentres - rotationCentre_) ^ rotationAxis_);
52 
53  const scalarField magRadSqr(magSqr(n));
54 
55  label facei = findMax(magRadSqr);
56 
57  if (debug)
58  {
59  Info<< "findFaceMaxRadius(const pointField&) : patch: " << name() << nl
60  << " rotFace = " << facei << nl
61  << " point = " << faceCentres[facei] << nl
62  << " distance = " << Foam::sqrt(magRadSqr[facei])
63  << endl;
64  }
65 
66  return n[facei];
67 }
68 
69 
71 (
72  const primitivePatch& half0,
73  const pointField& half0Ctrs,
74  const vectorField& half0Areas,
75  const pointField& half1Ctrs,
76  const vectorField& half1Areas
77 )
78 {
79  if (transform() != neighbPatch().transform())
80  {
82  << "Patch " << name()
83  << " has transform type " << transformTypeNames[transform()]
84  << ", neighbour patch " << neighbPatchName()
85  << " has transform type "
86  << neighbPatch().transformTypeNames[neighbPatch().transform()]
87  << exit(FatalError);
88  }
89 
90 
91  // Calculate transformation tensors
92 
93  switch (transform())
94  {
95  case ROTATIONAL:
96  {
97  tensor revT = Zero;
98 
99  if (rotationAngleDefined_)
100  {
101  const tensor T(rotationAxis_*rotationAxis_);
102 
103  const tensor S
104  (
105  0, -rotationAxis_.z(), rotationAxis_.y(),
106  rotationAxis_.z(), 0, -rotationAxis_.x(),
107  -rotationAxis_.y(), rotationAxis_.x(), 0
108  );
109 
110  const tensor revTPos
111  (
112  T
113  + cos(rotationAngle_)*(tensor::I - T)
114  + sin(rotationAngle_)*S
115  );
116 
117  const tensor revTNeg
118  (
119  T
120  + cos(-rotationAngle_)*(tensor::I - T)
121  + sin(-rotationAngle_)*S
122  );
123 
124  // Check - assume correct angle when difference in face areas
125  // is the smallest
126  const vector transformedAreaPos = gSum(half1Areas & revTPos);
127  const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
128  const vector area0 = gSum(half0Areas);
129  const scalar magArea0 = mag(area0) + ROOTVSMALL;
130 
131  // Areas have opposite sign, so sum should be zero when correct
132  // rotation applied
133  const scalar errorPos = mag(transformedAreaPos + area0);
134  const scalar errorNeg = mag(transformedAreaNeg + area0);
135 
136  const scalar normErrorPos = errorPos/magArea0;
137  const scalar normErrorNeg = errorNeg/magArea0;
138 
139  if (errorPos > errorNeg && normErrorNeg < matchTolerance())
140  {
141  revT = revTNeg;
142  rotationAngle_ *= -1;
143  }
144  else
145  {
146  revT = revTPos;
147  }
148 
149  const scalar areaError = min(normErrorPos, normErrorNeg);
150 
151  if (areaError > matchTolerance())
152  {
154  << "Patch areas are not consistent within "
155  << 100*matchTolerance()
156  << " % indicating a possible error in the specified "
157  << "angle of rotation" << nl
158  << " owner patch : " << name() << nl
159  << " neighbour patch : " << neighbPatch().name()
160  << nl
161  << " angle : "
162  << radToDeg(rotationAngle_) << " deg" << nl
163  << " area error : " << 100*areaError << " %"
164  << " match tolerance : " << matchTolerance()
165  << endl;
166  }
167 
168  if (debug)
169  {
170  scalar theta = radToDeg(rotationAngle_);
171 
172  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
173  << name()
174  << " Specified rotation:"
175  << " swept angle: " << theta << " [deg]"
176  << " reverse transform: " << revT
177  << endl;
178  }
179  }
180  else
181  {
182  point n0 = Zero;
183  point n1 = Zero;
184  if (half0Ctrs.size())
185  {
186  n0 = findFaceNormalMaxRadius(half0Ctrs);
187  }
188  if (half1Ctrs.size())
189  {
190  n1 = -findFaceNormalMaxRadius(half1Ctrs);
191  }
192 
193  reduce(n0, maxMagSqrOp<point>());
194  reduce(n1, maxMagSqrOp<point>());
195 
196  n0 /= mag(n0) + VSMALL;
197  n1 /= mag(n1) + VSMALL;
198 
199  // Extended tensor from two local coordinate systems calculated
200  // using normal and rotation axis
201  const tensor E0
202  (
203  rotationAxis_,
204  (n0 ^ rotationAxis_),
205  n0
206  );
207  const tensor E1
208  (
209  rotationAxis_,
210  (-n1 ^ rotationAxis_),
211  -n1
212  );
213  revT = E1.T() & E0;
214 
215  if (debug)
216  {
217  scalar theta = radToDeg(acos(-(n0 & n1)));
218 
219  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
220  << name()
221  << " Specified rotation:"
222  << " n0:" << n0 << " n1:" << n1
223  << " swept angle: " << theta << " [deg]"
224  << " reverse transform: " << revT
225  << endl;
226  }
227  }
228 
229  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
230  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
231  const_cast<vectorField&>(separation()).setSize(0);
232  const_cast<boolList&>(collocated()) = boolList(1, false);
233 
234  break;
235  }
236  case TRANSLATIONAL:
237  {
238  if (debug)
239  {
240  Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
241  << " Specified translation : " << separationVector_
242  << endl;
243  }
244 
245  const_cast<tensorField&>(forwardT()).clear();
246  const_cast<tensorField&>(reverseT()).clear();
247  const_cast<vectorField&>(separation()) = vectorField
248  (
249  1,
250  separationVector_
251  );
252  const_cast<boolList&>(collocated()) = boolList(1, false);
253 
254  break;
255  }
256  default:
257  {
258  if (debug)
259  {
260  Pout<< "patch:" << name()
261  << " Assuming cyclic AMI pairs are colocated" << endl;
262  }
263 
264  const_cast<tensorField&>(forwardT()).clear();
265  const_cast<tensorField&>(reverseT()).clear();
266  const_cast<vectorField&>(separation()).setSize(0);
267  const_cast<boolList&>(collocated()) = boolList(1, true);
268 
269  break;
270  }
271  }
272 
273  if (debug)
274  {
275  Pout<< "patch: " << name() << nl
276  << " forwardT = " << forwardT() << nl
277  << " reverseT = " << reverseT() << nl
278  << " separation = " << separation() << nl
279  << " collocated = " << collocated() << nl << endl;
280  }
281 }
282 
283 
284 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
285 
287 (
289 ) const
290 {
291  if (owner())
292  {
293  AMIPtr_.clear();
294 
295  const polyPatch& nbr = neighbPatch();
296  pointField nbrPoints
297  (
298  neighbPatch().boundaryMesh().mesh().points(),
299  neighbPatch().meshPoints()
300  );
301 
302  if (debug)
303  {
304  const Time& t = boundaryMesh().mesh().time();
305  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
306  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
307  }
308 
309  // Transform neighbour patch to local system
310  transformPosition(nbrPoints);
311  primitivePatch nbrPatch0
312  (
314  (
315  nbr.localFaces(),
316  nbr.size()
317  ),
318  nbrPoints
319  );
320 
321  if (debug)
322  {
323  const Time& t = boundaryMesh().mesh().time();
324  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
325  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
326 
327  OFstream osO(t.path()/name() + "_ownerPatch.obj");
328  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
329  }
330 
331  // Construct/apply AMI interpolation to determine addressing and weights
332  AMIPtr_.reset
333  (
335  (
336  *this,
337  nbrPatch0,
338  surfPtr(),
340  AMIRequireMatch_,
341  AMIMethod,
342  AMILowWeightCorrection_,
343  AMIReverse_
344  )
345  );
346 
347  if (debug)
348  {
349  Pout<< "cyclicAMIPolyPatch : " << name()
350  << " constructed AMI with " << nl
351  << " " << "srcAddress:" << AMIPtr_().srcAddress().size()
352  << nl
353  << " " << "tgAddress :" << AMIPtr_().tgtAddress().size()
354  << nl << endl;
355  }
356  }
357 }
358 
359 
361 {
362  const cyclicAMIPolyPatch& half0 = *this;
363  vectorField half0Areas(half0.size());
364  forAll(half0, facei)
365  {
366  half0Areas[facei] = half0[facei].normal(half0.points());
367  }
368 
369  const cyclicAMIPolyPatch& half1 = neighbPatch();
370  vectorField half1Areas(half1.size());
371  forAll(half1, facei)
372  {
373  half1Areas[facei] = half1[facei].normal(half1.points());
374  }
375 
376  calcTransforms
377  (
378  half0,
379  half0.faceCentres(),
380  half0Areas,
381  half1.faceCentres(),
382  half1Areas
383  );
384 
385  if (debug)
386  {
387  Pout<< "calcTransforms() : patch: " << name() << nl
388  << " forwardT = " << forwardT() << nl
389  << " reverseT = " << reverseT() << nl
390  << " separation = " << separation() << nl
391  << " collocated = " << collocated() << nl << endl;
392  }
393 }
394 
395 
397 {
398  // Clear the invalid AMI
399  AMIPtr_.clear();
400 
402 }
403 
404 
406 {
407  calcGeometry
408  (
409  *this,
410  faceCentres(),
411  faceAreas(),
412  faceCellCentres(),
413  neighbPatch().faceCentres(),
414  neighbPatch().faceAreas(),
415  neighbPatch().faceCellCentres()
416  );
417 }
418 
419 
421 (
422  PstreamBuffers& pBufs,
423  const pointField& p
424 )
425 {
426  // Clear the invalid AMI
427  AMIPtr_.clear();
428 
429  polyPatch::initMovePoints(pBufs, p);
430 
431  // See below. Clear out any local geometry
433 }
434 
435 
437 (
438  PstreamBuffers& pBufs,
439  const pointField& p
440 )
441 {
442  polyPatch::movePoints(pBufs, p);
443 
444  calcTransforms();
445 }
446 
447 
449 {
450  // Clear the invalid AMI
451  AMIPtr_.clear();
452 
454 }
455 
456 
458 {
459  polyPatch::updateMesh(pBufs);
460 }
461 
462 
464 {
465  AMIPtr_.clear();
467 }
468 
469 
470 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
471 
473 (
474  const word& name,
475  const label size,
476  const label start,
477  const label index,
478  const polyBoundaryMesh& bm,
479  const word& patchType,
480  const transformType transform
481 )
482 :
483  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
484  nbrPatchName_(word::null),
485  nbrPatchID_(-1),
486  rotationAxis_(Zero),
487  rotationCentre_(point::zero),
488  rotationAngleDefined_(false),
489  rotationAngle_(0.0),
490  separationVector_(Zero),
491  AMIPtr_(nullptr),
492  AMIReverse_(false),
493  AMIRequireMatch_(true),
494  AMILowWeightCorrection_(-1.0),
495  surfPtr_(nullptr),
496  surfDict_(fileName("surface"))
497 {
498  // Neighbour patch might not be valid yet so no transformation
499  // calculation possible
500 }
501 
502 
504 (
505  const word& name,
506  const dictionary& dict,
507  const label index,
508  const polyBoundaryMesh& bm,
509  const word& patchType
510 )
511 :
512  coupledPolyPatch(name, dict, index, bm, patchType),
513  nbrPatchName_(dict.lookupOrDefault<word>("neighbourPatch", "")),
514  coupleGroup_(dict),
515  nbrPatchID_(-1),
516  rotationAxis_(Zero),
517  rotationCentre_(point::zero),
518  rotationAngleDefined_(false),
519  rotationAngle_(0.0),
520  separationVector_(Zero),
521  AMIPtr_(nullptr),
522  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
523  AMIRequireMatch_(true),
524  AMILowWeightCorrection_(dict.lookupOrDefault("lowWeightCorrection", -1.0)),
525  surfPtr_(nullptr),
526  surfDict_(dict.subOrEmptyDict("surface"))
527 {
528  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
529  {
531  (
532  dict
533  ) << "No \"neighbourPatch\" or \"coupleGroup\" provided."
534  << exit(FatalIOError);
535  }
536 
537  if (nbrPatchName_ == name)
538  {
540  (
541  dict
542  ) << "Neighbour patch name " << nbrPatchName_
543  << " cannot be the same as this patch " << name
544  << exit(FatalIOError);
545  }
546 
547  switch (transform())
548  {
549  case ROTATIONAL:
550  {
551  dict.lookup("rotationAxis") >> rotationAxis_;
552  dict.lookup("rotationCentre") >> rotationCentre_;
553  if (dict.readIfPresent("rotationAngle", rotationAngle_))
554  {
555  rotationAngleDefined_ = true;
556  rotationAngle_ = degToRad(rotationAngle_);
557 
558  if (debug)
559  {
560  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
561  << endl;
562  }
563  }
564 
565  scalar magRot = mag(rotationAxis_);
566  if (magRot < SMALL)
567  {
569  (
570  dict
571  ) << "Illegal rotationAxis " << rotationAxis_ << endl
572  << "Please supply a non-zero vector."
573  << exit(FatalIOError);
574  }
575  rotationAxis_ /= magRot;
576 
577  break;
578  }
579  case TRANSLATIONAL:
580  {
581  dict.lookup("separationVector") >> separationVector_;
582  break;
583  }
584  default:
585  {
586  // No additional info required
587  }
588  }
589 
590  // Neighbour patch might not be valid yet so no transformation
591  // calculation possible
592 }
593 
594 
596 (
597  const cyclicAMIPolyPatch& pp,
598  const polyBoundaryMesh& bm
599 )
600 :
601  coupledPolyPatch(pp, bm),
602  nbrPatchName_(pp.nbrPatchName_),
603  coupleGroup_(pp.coupleGroup_),
604  nbrPatchID_(-1),
605  rotationAxis_(pp.rotationAxis_),
606  rotationCentre_(pp.rotationCentre_),
607  rotationAngleDefined_(pp.rotationAngleDefined_),
608  rotationAngle_(pp.rotationAngle_),
609  separationVector_(pp.separationVector_),
610  AMIPtr_(nullptr),
611  AMIReverse_(pp.AMIReverse_),
612  AMIRequireMatch_(pp.AMIRequireMatch_),
613  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
614  surfPtr_(nullptr),
615  surfDict_(pp.surfDict_)
616 {
617  // Neighbour patch might not be valid yet so no transformation
618  // calculation possible
619 }
620 
621 
623 (
624  const cyclicAMIPolyPatch& pp,
625  const polyBoundaryMesh& bm,
626  const label index,
627  const label newSize,
628  const label newStart,
629  const word& nbrPatchName
630 )
631 :
632  coupledPolyPatch(pp, bm, index, newSize, newStart),
633  nbrPatchName_(nbrPatchName),
634  coupleGroup_(pp.coupleGroup_),
635  nbrPatchID_(-1),
636  rotationAxis_(pp.rotationAxis_),
637  rotationCentre_(pp.rotationCentre_),
638  rotationAngleDefined_(pp.rotationAngleDefined_),
639  rotationAngle_(pp.rotationAngle_),
640  separationVector_(pp.separationVector_),
641  AMIPtr_(nullptr),
642  AMIReverse_(pp.AMIReverse_),
643  AMIRequireMatch_(pp.AMIRequireMatch_),
644  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
645  surfPtr_(nullptr),
646  surfDict_(pp.surfDict_)
647 {
648  if (nbrPatchName_ == name())
649  {
651  << "Neighbour patch name " << nbrPatchName_
652  << " cannot be the same as this patch " << name()
653  << exit(FatalError);
654  }
655 
656  // Neighbour patch might not be valid yet so no transformation
657  // calculation possible
658 }
659 
660 
662 (
663  const cyclicAMIPolyPatch& pp,
664  const polyBoundaryMesh& bm,
665  const label index,
666  const labelUList& mapAddressing,
667  const label newStart
668 )
669 :
670  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
671  nbrPatchName_(pp.nbrPatchName_),
672  coupleGroup_(pp.coupleGroup_),
673  nbrPatchID_(-1),
674  rotationAxis_(pp.rotationAxis_),
675  rotationCentre_(pp.rotationCentre_),
676  rotationAngleDefined_(pp.rotationAngleDefined_),
677  rotationAngle_(pp.rotationAngle_),
678  separationVector_(pp.separationVector_),
679  AMIPtr_(nullptr),
680  AMIReverse_(pp.AMIReverse_),
681  AMIRequireMatch_(pp.AMIRequireMatch_),
682  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
683  surfPtr_(nullptr),
684  surfDict_(pp.surfDict_)
685 {}
686 
687 
688 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
689 
691 {}
692 
693 
694 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
695 
697 {
698  if (nbrPatchID_ == -1)
699  {
700  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
701 
702  if (nbrPatchID_ == -1)
703  {
705  << "Illegal neighbourPatch name " << neighbPatchName()
706  << nl << "Valid patch names are "
707  << this->boundaryMesh().names()
708  << exit(FatalError);
709  }
710 
711  // Check that it is a cyclic AMI patch
712  const cyclicAMIPolyPatch& nbrPatch =
713  refCast<const cyclicAMIPolyPatch>
714  (
715  this->boundaryMesh()[nbrPatchID_]
716  );
717 
718  if (nbrPatch.neighbPatchName() != name())
719  {
721  << "Patch " << name()
722  << " specifies neighbour patch " << neighbPatchName()
723  << nl << " but that in return specifies "
724  << nbrPatch.neighbPatchName() << endl;
725  }
726  }
727 
728  return nbrPatchID_;
729 }
730 
731 
733 {
734  return index() < neighbPatchID();
735 }
736 
737 
739 {
740  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
741  return refCast<const cyclicAMIPolyPatch>(pp);
742 }
743 
744 
747 {
748  const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
749 
750  if (!surfPtr_.valid() && owner() && surfType != "none")
751  {
752  word surfName(surfDict_.lookupOrDefault("name", name()));
753 
754  const polyMesh& mesh = boundaryMesh().mesh();
755 
756  surfPtr_ =
758  (
759  surfType,
760  IOobject
761  (
762  surfName,
763  mesh.time().constant(),
764  "triSurface",
765  mesh,
768  ),
769  surfDict_
770  );
771  }
772 
773  return surfPtr_;
774 }
775 
776 
778 {
779  if (!owner())
780  {
782  << "AMI interpolator only available to owner patch"
783  << abort(FatalError);
784  }
785 
786  if (!AMIPtr_.valid())
787  {
788  resetAMI();
789  }
790 
791  return AMIPtr_();
792 }
793 
794 
796 {
797  if (owner())
798  {
799  return AMI().applyLowWeightCorrection();
800  }
801  else
802  {
803  return neighbPatch().AMI().applyLowWeightCorrection();
804  }
805 }
806 
807 
809 {
810  if (!parallel())
811  {
812  if (transform() == ROTATIONAL)
813  {
814  l = Foam::transform(forwardT(), l - rotationCentre_)
815  + rotationCentre_;
816  }
817  else
818  {
819  l = Foam::transform(forwardT(), l);
820  }
821  }
822  else if (separated())
823  {
824  // transformPosition gets called on the receiving side,
825  // separation gets calculated on the sending side so subtract
826 
827  const vectorField& s = separation();
828  if (s.size() == 1)
829  {
830  forAll(l, i)
831  {
832  l[i] -= s[0];
833  }
834  }
835  else
836  {
837  l -= s;
838  }
839  }
840 }
841 
842 
844 (
845  point& l,
846  const label facei
847 ) const
848 {
849  if (!parallel())
850  {
851  const tensor& T =
852  (
853  forwardT().size() == 1
854  ? forwardT()[0]
855  : forwardT()[facei]
856  );
857 
858  if (transform() == ROTATIONAL)
859  {
860  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
861  }
862  else
863  {
864  l = Foam::transform(T, l);
865  }
866  }
867  else if (separated())
868  {
869  const vector& s =
870  (
871  separation().size() == 1
872  ? separation()[0]
873  : separation()[facei]
874  );
875 
876  l -= s;
877  }
878 }
879 
880 
882 (
883  point& l,
884  const label facei
885 ) const
886 {
887  if (!parallel())
888  {
889  const tensor& T =
890  (
891  reverseT().size() == 1
892  ? reverseT()[0]
893  : reverseT()[facei]
894  );
895 
896  if (transform() == ROTATIONAL)
897  {
898  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
899  }
900  else
901  {
902  l = Foam::transform(T, l);
903  }
904  }
905  else if (separated())
906  {
907  const vector& s =
908  (
909  separation().size() == 1
910  ? separation()[0]
911  : separation()[facei]
912  );
913 
914  l += s;
915  }
916 }
917 
918 
920 (
921  vector& d,
922  const label facei
923 ) const
924 {
925  if (!parallel())
926  {
927  const tensor& T =
928  (
929  reverseT().size() == 1
930  ? reverseT()[0]
931  : reverseT()[facei]
932  );
933 
934  d = Foam::transform(T, d);
935  }
936 }
937 
938 
940 (
941  const primitivePatch& referPatch,
942  const pointField& thisCtrs,
943  const vectorField& thisAreas,
944  const pointField& thisCc,
945  const pointField& nbrCtrs,
946  const vectorField& nbrAreas,
947  const pointField& nbrCc
948 )
949 {
950  calcTransforms
951  (
952  referPatch,
953  thisCtrs,
954  thisAreas,
955  nbrCtrs,
956  nbrAreas
957  );
958 }
959 
960 
962 (
963  PstreamBuffers& pBufs,
964  const primitivePatch& pp
965 ) const
966 {}
967 
968 
970 (
971  PstreamBuffers& pBufs,
972  const primitivePatch& pp,
973  labelList& faceMap,
974  labelList& rotation
975 ) const
976 {
977  faceMap.setSize(pp.size());
978  faceMap = -1;
979 
980  rotation.setSize(pp.size());
981  rotation = 0;
982 
983  return false;
984 }
985 
986 
988 (
989  const label facei,
990  const vector& n,
991  point& p
992 ) const
993 {
994  point prt(p);
995  reverseTransformPosition(prt, facei);
996 
997  vector nrt(n);
998  reverseTransformDirection(nrt, facei);
999 
1000  label nbrFacei = -1;
1001 
1002  if (owner())
1003  {
1004  nbrFacei = AMI().tgtPointFace
1005  (
1006  *this,
1007  neighbPatch(),
1008  nrt,
1009  facei,
1010  prt
1011  );
1012  }
1013  else
1014  {
1015  nbrFacei = neighbPatch().AMI().srcPointFace
1016  (
1017  neighbPatch(),
1018  *this,
1019  nrt,
1020  facei,
1021  prt
1022  );
1023  }
1024 
1025  if (nbrFacei >= 0)
1026  {
1027  p = prt;
1028  }
1029 
1030  return nbrFacei;
1031 }
1032 
1033 
1035 {
1037  if (!nbrPatchName_.empty())
1038  {
1039  os.writeKeyword("neighbourPatch") << nbrPatchName_
1040  << token::END_STATEMENT << nl;
1041  }
1042  coupleGroup_.write(os);
1043 
1044  switch (transform())
1045  {
1046  case ROTATIONAL:
1047  {
1048  os.writeKeyword("rotationAxis") << rotationAxis_
1049  << token::END_STATEMENT << nl;
1050  os.writeKeyword("rotationCentre") << rotationCentre_
1051  << token::END_STATEMENT << nl;
1052 
1053  if (rotationAngleDefined_)
1054  {
1055  os.writeKeyword("rotationAngle") << radToDeg(rotationAngle_)
1056  << token::END_STATEMENT << nl;
1057  }
1058 
1059  break;
1060  }
1061  case TRANSLATIONAL:
1062  {
1063  os.writeKeyword("separationVector") << separationVector_
1064  << token::END_STATEMENT << nl;
1065  break;
1066  }
1067  case NOORDERING:
1068  {
1069  break;
1070  }
1071  default:
1072  {
1073  // No additional info to write
1074  }
1075  }
1076 
1077  if (AMIReverse_)
1078  {
1079  os.writeKeyword("flipNormals") << AMIReverse_
1080  << token::END_STATEMENT << nl;
1081  }
1082 
1083  if (AMILowWeightCorrection_ > 0)
1084  {
1085  os.writeKeyword("lowWeightCorrection") << AMILowWeightCorrection_
1086  << token::END_STATEMENT << nl;
1087  }
1088 
1089  if (!surfDict_.empty())
1090  {
1091  os.writeKeyword(surfDict_.dictName());
1092  os << surfDict_;
1093  }
1094 }
1095 
1096 
1097 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:115
vector separationVector_
Translation vector.
dimensionedScalar acos(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
fileName path() const
Return path.
Definition: Time.H:269
virtual label neighbPatchID() const
Neighbour patch ID.
virtual void movePoints(const Field< PointType > &)
Correct patch after moving points.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
scalar radToDeg(const scalar rad)
Conversion from radians to degrees.
A class for handling file names.
Definition: fileName.H:69
bool AMIRequireMatch_
Flag to indicate that patches should match/overlap.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:115
tUEqn clear()
Output to file stream.
Definition: OFstream.H:82
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual void calcTransforms()
Recalculate the transformation tensors.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
vector rotationAxis_
Axis of rotation for rotational cyclics.
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
virtual void initMovePoints(PstreamBuffers &, const pointField &)
Initialise the patches for moving points.
Definition: polyPatch.H:108
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:69
interpolationMethod
Enumeration specifying interpolation method.
virtual void movePoints(PstreamBuffers &, const pointField &p)
Correct patches after moving points.
Definition: polyPatch.C:57
const dictionary surfDict_
Dictionary used during projection surface construction.
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
points setSize(newPointi)
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from (base couped patch) components.
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
const word & neighbPatchName() const
Neighbour patch name.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar cos(const dimensionedScalar &ds)
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
const autoPtr< searchableSurface > & surfPtr() const
Return a reference to the projection surface.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const scalar AMILowWeightCorrection_
Low weight correction threshold for AMI.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
static const zero Zero
Definition: zero.H:91
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:62
Cyclic patch for Arbitrary Mesh Interface (AMI)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
virtual ~cyclicAMIPolyPatch()
Destructor.
Foam::polyBoundaryMesh.
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
static const Tensor I
Definition: Tensor.H:80
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:100
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
void setSize(const label)
Reset size of List.
Definition: List.C:281
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
vector point
Point is a vector.
Definition: point.H:41
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual bool owner() const
Does this side own the patch?
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const bMesh & mesh() const
Definition: boundaryMesh.H:203
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
word nbrPatchName_
Name of other half.
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
scalar rotationAngle_
Rotation angle.
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p.
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:727
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
IOerror FatalIOError
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284