cyclicAMIPolyPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 {
288  if (owner())
289  {
290  const polyPatch& nbr = neighbPatch();
291  pointField nbrPoints
292  (
293  neighbPatch().boundaryMesh().mesh().points(),
294  neighbPatch().meshPoints()
295  );
296 
297  if (debug)
298  {
299  const Time& t = boundaryMesh().mesh().time();
300  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
301  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
302  }
303 
304  // Transform neighbour patch to local system
305  transformPosition(nbrPoints);
306  primitivePatch nbrPatch0
307  (
309  (
310  nbr.localFaces(),
311  nbr.size()
312  ),
313  nbrPoints
314  );
315 
316  if (debug)
317  {
318  const Time& t = boundaryMesh().mesh().time();
319  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
320  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
321 
322  OFstream osO(t.path()/name() + "_ownerPatch.obj");
323  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
324  }
325 
326  // Construct/apply AMI interpolation to determine addressing and weights
327  AMIs_.resize(1);
328  AMIs_.set
329  (
330  0,
332  (
333  *this,
334  nbrPatch0,
335  surfPtr(),
337  AMIRequireMatch_,
338  AMIMethod_,
339  AMILowWeightCorrection_,
340  AMIReverse_
341  )
342  );
343 
344  AMITransforms_.resize(1, vectorTensorTransform::I);
345 
346  if (debug)
347  {
348  Pout<< "cyclicAMIPolyPatch : " << name()
349  << " constructed AMI with " << nl
350  << " " << "srcAddress:" << AMIs_[0].srcAddress().size()
351  << nl
352  << " " << "tgAddress :" << AMIs_[0].tgtAddress().size()
353  << nl << endl;
354  }
355  }
356 }
357 
358 
360 {
361  const cyclicAMIPolyPatch& half0 = *this;
362  vectorField half0Areas(half0.size());
363  forAll(half0, facei)
364  {
365  half0Areas[facei] = half0[facei].area(half0.points());
366  }
367 
368  const cyclicAMIPolyPatch& half1 = neighbPatch();
369  vectorField half1Areas(half1.size());
370  forAll(half1, facei)
371  {
372  half1Areas[facei] = half1[facei].area(half1.points());
373  }
374 
375  calcTransforms
376  (
377  half0,
378  half0.faceCentres(),
379  half0Areas,
380  half1.faceCentres(),
381  half1Areas
382  );
383 
384  if (debug)
385  {
386  Pout<< "calcTransforms() : patch: " << name() << nl
387  << " forwardT = " << forwardT() << nl
388  << " reverseT = " << reverseT() << nl
389  << " separation = " << separation() << nl
390  << " collocated = " << collocated() << nl << endl;
391  }
392 }
393 
394 
396 {
397  // Clear the invalid AMIs and transforms
398  AMIs_.clear();
399  AMITransforms_.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 AMIs and transforms
427  AMIs_.clear();
428  AMITransforms_.clear();
429 
430  polyPatch::initMovePoints(pBufs, p);
431 
432  // See below. Clear out any local geometry
434 }
435 
436 
438 (
439  PstreamBuffers& pBufs,
440  const pointField& p
441 )
442 {
443  polyPatch::movePoints(pBufs, p);
444 
445  calcTransforms();
446 }
447 
448 
450 {
451  // Clear the invalid AMIs and transforms
452  AMIs_.clear();
453  AMITransforms_.clear();
454 
456 }
457 
458 
460 {
461  polyPatch::updateMesh(pBufs);
462 }
463 
464 
466 {
467  // Clear the invalid AMIs and transforms
468  AMIs_.clear();
469  AMITransforms_.clear();
470 
472 }
473 
474 
475 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
476 
478 (
479  const word& name,
480  const label size,
481  const label start,
482  const label index,
483  const polyBoundaryMesh& bm,
484  const word& patchType,
485  const transformType transform,
486  const bool AMIRequireMatch,
488 )
489 :
490  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
491  nbrPatchName_(word::null),
492  nbrPatchID_(-1),
493  rotationAxis_(Zero),
494  rotationCentre_(point::zero),
495  rotationAngleDefined_(false),
496  rotationAngle_(0.0),
497  separationVector_(Zero),
498  AMIs_(),
499  AMITransforms_(),
500  AMIReverse_(false),
501  AMIRequireMatch_(AMIRequireMatch),
502  AMILowWeightCorrection_(-1.0),
503  AMIMethod_(AMIMethod),
504  surfPtr_(nullptr),
505  surfDict_(fileName("surface"))
506 {
507  // Neighbour patch might not be valid yet so no transformation
508  // calculation possible
509 }
510 
511 
513 (
514  const word& name,
515  const dictionary& dict,
516  const label index,
517  const polyBoundaryMesh& bm,
518  const word& patchType,
519  const bool AMIRequireMatch,
521 )
522 :
523  coupledPolyPatch(name, dict, index, bm, patchType),
524  nbrPatchName_(dict.lookupOrDefault<word>("neighbourPatch", "")),
525  coupleGroup_(dict),
526  nbrPatchID_(-1),
527  rotationAxis_(Zero),
528  rotationCentre_(point::zero),
529  rotationAngleDefined_(false),
530  rotationAngle_(0.0),
531  separationVector_(Zero),
532  AMIs_(),
533  AMITransforms_(),
534  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
535  AMIRequireMatch_(AMIRequireMatch),
536  AMILowWeightCorrection_(dict.lookupOrDefault("lowWeightCorrection", -1.0)),
537  AMIMethod_
538  (
539  dict.found("method")
541  (
542  dict.lookup("method")
543  )
544  : AMIMethod
545  ),
546  surfPtr_(nullptr),
547  surfDict_(dict.subOrEmptyDict("surface"))
548 {
549  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
550  {
552  (
553  dict
554  ) << "No \"neighbourPatch\" or \"coupleGroup\" provided."
555  << exit(FatalIOError);
556  }
557 
558  if (nbrPatchName_ == name)
559  {
561  (
562  dict
563  ) << "Neighbour patch name " << nbrPatchName_
564  << " cannot be the same as this patch " << name
565  << exit(FatalIOError);
566  }
567 
568  switch (transform())
569  {
570  case ROTATIONAL:
571  {
572  dict.lookup("rotationAxis") >> rotationAxis_;
573  dict.lookup("rotationCentre") >> rotationCentre_;
574  if (dict.readIfPresent("rotationAngle", rotationAngle_))
575  {
576  rotationAngleDefined_ = true;
577  rotationAngle_ = degToRad(rotationAngle_);
578 
579  if (debug)
580  {
581  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
582  << endl;
583  }
584  }
585 
586  scalar magRot = mag(rotationAxis_);
587  if (magRot < small)
588  {
590  (
591  dict
592  ) << "Illegal rotationAxis " << rotationAxis_ << endl
593  << "Please supply a non-zero vector."
594  << exit(FatalIOError);
595  }
596  rotationAxis_ /= magRot;
597 
598  break;
599  }
600  case TRANSLATIONAL:
601  {
602  dict.lookup("separationVector") >> separationVector_;
603  break;
604  }
605  default:
606  {
607  // No additional info required
608  }
609  }
610 
611  // Neighbour patch might not be valid yet so no transformation
612  // calculation possible
613 }
614 
615 
617 (
618  const cyclicAMIPolyPatch& pp,
619  const polyBoundaryMesh& bm
620 )
621 :
622  coupledPolyPatch(pp, bm),
623  nbrPatchName_(pp.nbrPatchName_),
624  coupleGroup_(pp.coupleGroup_),
625  nbrPatchID_(-1),
626  rotationAxis_(pp.rotationAxis_),
627  rotationCentre_(pp.rotationCentre_),
628  rotationAngleDefined_(pp.rotationAngleDefined_),
629  rotationAngle_(pp.rotationAngle_),
630  separationVector_(pp.separationVector_),
631  AMIs_(),
632  AMITransforms_(),
633  AMIReverse_(pp.AMIReverse_),
634  AMIRequireMatch_(pp.AMIRequireMatch_),
635  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
636  AMIMethod_(pp.AMIMethod_),
637  surfPtr_(nullptr),
638  surfDict_(pp.surfDict_)
639 {
640  // Neighbour patch might not be valid yet so no transformation
641  // calculation possible
642 }
643 
644 
646 (
647  const cyclicAMIPolyPatch& pp,
648  const polyBoundaryMesh& bm,
649  const label index,
650  const label newSize,
651  const label newStart,
652  const word& nbrPatchName
653 )
654 :
655  coupledPolyPatch(pp, bm, index, newSize, newStart),
656  nbrPatchName_(nbrPatchName),
657  coupleGroup_(pp.coupleGroup_),
658  nbrPatchID_(-1),
659  rotationAxis_(pp.rotationAxis_),
660  rotationCentre_(pp.rotationCentre_),
661  rotationAngleDefined_(pp.rotationAngleDefined_),
662  rotationAngle_(pp.rotationAngle_),
663  separationVector_(pp.separationVector_),
664  AMIs_(),
665  AMITransforms_(),
666  AMIReverse_(pp.AMIReverse_),
667  AMIRequireMatch_(pp.AMIRequireMatch_),
668  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
669  AMIMethod_(pp.AMIMethod_),
670  surfPtr_(nullptr),
671  surfDict_(pp.surfDict_)
672 {
673  if (nbrPatchName_ == name())
674  {
676  << "Neighbour patch name " << nbrPatchName_
677  << " cannot be the same as this patch " << name()
678  << exit(FatalError);
679  }
680 
681  // Neighbour patch might not be valid yet so no transformation
682  // calculation possible
683 }
684 
685 
687 (
688  const cyclicAMIPolyPatch& pp,
689  const polyBoundaryMesh& bm,
690  const label index,
691  const labelUList& mapAddressing,
692  const label newStart
693 )
694 :
695  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
696  nbrPatchName_(pp.nbrPatchName_),
697  coupleGroup_(pp.coupleGroup_),
698  nbrPatchID_(-1),
699  rotationAxis_(pp.rotationAxis_),
700  rotationCentre_(pp.rotationCentre_),
701  rotationAngleDefined_(pp.rotationAngleDefined_),
702  rotationAngle_(pp.rotationAngle_),
703  separationVector_(pp.separationVector_),
704  AMIs_(),
705  AMITransforms_(),
706  AMIReverse_(pp.AMIReverse_),
707  AMIRequireMatch_(pp.AMIRequireMatch_),
708  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
709  AMIMethod_(pp.AMIMethod_),
710  surfPtr_(nullptr),
711  surfDict_(pp.surfDict_)
712 {}
713 
714 
715 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
716 
718 {}
719 
720 
721 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
722 
724 {
725  if (nbrPatchID_ == -1)
726  {
727  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
728 
729  if (nbrPatchID_ == -1)
730  {
732  << "Illegal neighbourPatch name " << neighbPatchName()
733  << nl << "Valid patch names are "
734  << this->boundaryMesh().names()
735  << exit(FatalError);
736  }
737 
738  // Check that it is a cyclic AMI patch
739  const cyclicAMIPolyPatch& nbrPatch =
740  refCast<const cyclicAMIPolyPatch>
741  (
742  this->boundaryMesh()[nbrPatchID_]
743  );
744 
745  if (nbrPatch.neighbPatchName() != name())
746  {
748  << "Patch " << name()
749  << " specifies neighbour patch " << neighbPatchName()
750  << nl << " but that in return specifies "
751  << nbrPatch.neighbPatchName() << endl;
752  }
753  }
754 
755  return nbrPatchID_;
756 }
757 
758 
760 {
761  return index() < neighbPatchID();
762 }
763 
764 
766 {
767  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
768  return refCast<const cyclicAMIPolyPatch>(pp);
769 }
770 
771 
774 {
775  const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
776 
777  if (!surfPtr_.valid() && owner() && surfType != "none")
778  {
779  word surfName(surfDict_.lookupOrDefault("name", name()));
780 
781  const polyMesh& mesh = boundaryMesh().mesh();
782 
783  surfPtr_ =
785  (
786  surfType,
787  IOobject
788  (
789  surfName,
790  mesh.time().constant(),
791  "triSurface",
792  mesh,
795  ),
796  surfDict_
797  );
798  }
799 
800  return surfPtr_;
801 }
802 
803 
806 {
807  if (!owner())
808  {
810  << "AMI interpolators only available to owner patch"
811  << abort(FatalError);
812  }
813 
814  if (AMIs_.empty())
815  {
816  resetAMI();
817  }
818 
819  return AMIs_;
820 }
821 
822 
825 {
826  if (!owner())
827  {
829  << "AMI transforms only available to owner patch"
830  << abort(FatalError);
831  }
832 
833  if (AMIs_.empty())
834  {
835  resetAMI();
836  }
837 
838  return AMITransforms_;
839 }
840 
841 
843 {
844  if (owner())
845  {
846  return AMILowWeightCorrection_ > 0;
847  }
848  else
849  {
850  return neighbPatch().AMILowWeightCorrection_ > 0;
851  }
852 }
853 
854 
856 {
857  if (owner())
858  {
859  return AMIs()[0].srcWeightsSum();
860  }
861  else
862  {
863  return neighbPatch().AMIs()[0].tgtWeightsSum();
864  }
865 }
866 
867 
869 {
870  if (owner())
871  {
872  return AMIs()[0].tgtWeightsSum();
873  }
874  else
875  {
876  return neighbPatch().AMIs()[0].srcWeightsSum();
877  }
878 }
879 
880 
882 {
883  if (!parallel())
884  {
885  if (transform() == ROTATIONAL)
886  {
887  l = Foam::transform(forwardT(), l - rotationCentre_)
888  + rotationCentre_;
889  }
890  else
891  {
892  l = Foam::transform(forwardT(), l);
893  }
894  }
895  else if (separated())
896  {
897  // transformPosition gets called on the receiving side,
898  // separation gets calculated on the sending side so subtract
899 
900  const vectorField& s = separation();
901  if (s.size() == 1)
902  {
903  forAll(l, i)
904  {
905  l[i] -= s[0];
906  }
907  }
908  else
909  {
910  l -= s;
911  }
912  }
913 }
914 
915 
917 (
918  point& l,
919  const label facei
920 ) const
921 {
922  if (!parallel())
923  {
924  const tensor& T =
925  (
926  forwardT().size() == 1
927  ? forwardT()[0]
928  : forwardT()[facei]
929  );
930 
931  if (transform() == ROTATIONAL)
932  {
933  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
934  }
935  else
936  {
937  l = Foam::transform(T, l);
938  }
939  }
940  else if (separated())
941  {
942  const vector& s =
943  (
944  separation().size() == 1
945  ? separation()[0]
946  : separation()[facei]
947  );
948 
949  l -= s;
950  }
951 }
952 
953 
955 (
956  point& l,
957  const label facei
958 ) const
959 {
960  if (!parallel())
961  {
962  const tensor& T =
963  (
964  reverseT().size() == 1
965  ? reverseT()[0]
966  : reverseT()[facei]
967  );
968 
969  if (transform() == ROTATIONAL)
970  {
971  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
972  }
973  else
974  {
975  l = Foam::transform(T, l);
976  }
977  }
978  else if (separated())
979  {
980  const vector& s =
981  (
982  separation().size() == 1
983  ? separation()[0]
984  : separation()[facei]
985  );
986 
987  l += s;
988  }
989 }
990 
991 
993 (
994  vector& d,
995  const label facei
996 ) const
997 {
998  if (!parallel())
999  {
1000  const tensor& T =
1001  (
1002  reverseT().size() == 1
1003  ? reverseT()[0]
1004  : reverseT()[facei]
1005  );
1006 
1007  d = Foam::transform(T, d);
1008  }
1009 }
1010 
1011 
1014  const scalarField& fld,
1015  const direction cmpt,
1016  const direction rank,
1017  const scalarUList& defaultValues
1018 ) const
1019 {
1020  const cyclicAMIPolyPatch& nei = neighbPatch();
1021 
1022  tmp<scalarField> result(new scalarField(size(), Zero));
1023 
1024  if (owner())
1025  {
1026  forAll(AMIs(), i)
1027  {
1028  const scalar r =
1029  pow(inv(AMITransforms()[i]).R()(cmpt, cmpt), rank);
1030 
1031  result.ref() +=
1032  AMIs()[i].interpolateToSource(r*fld, defaultValues);
1033  }
1034  }
1035  else
1036  {
1037  forAll(nei.AMIs(), i)
1038  {
1039  const scalar r =
1040  pow(nei.AMITransforms()[i].R()(cmpt, cmpt), rank);
1041 
1042  result.ref() +=
1043  nei.AMIs()[i].interpolateToTarget(r*fld, defaultValues);
1044  }
1045  }
1046 
1047  return result;
1048 }
1049 
1050 
1053  const primitivePatch& referPatch,
1054  const pointField& thisCtrs,
1055  const vectorField& thisAreas,
1056  const pointField& thisCc,
1057  const pointField& nbrCtrs,
1058  const vectorField& nbrAreas,
1059  const pointField& nbrCc
1060 )
1061 {
1062  calcTransforms
1063  (
1064  referPatch,
1065  thisCtrs,
1066  thisAreas,
1067  nbrCtrs,
1068  nbrAreas
1069  );
1070 }
1071 
1072 
1075  PstreamBuffers& pBufs,
1076  const primitivePatch& pp
1077 ) const
1078 {}
1079 
1080 
1083  PstreamBuffers& pBufs,
1084  const primitivePatch& pp,
1085  labelList& faceMap,
1086  labelList& rotation
1087 ) const
1088 {
1089  faceMap.setSize(pp.size());
1090  faceMap = -1;
1091 
1092  rotation.setSize(pp.size());
1093  rotation = 0;
1094 
1095  return false;
1096 }
1097 
1098 
1101  const label facei,
1102  const vector& n,
1103  point& p
1104 ) const
1105 {
1106  point pt(p);
1107  reverseTransformPosition(pt, facei);
1108 
1109  vector nt(n);
1110  reverseTransformDirection(nt, facei);
1111 
1112  if (owner())
1113  {
1114  forAll(AMIs(), i)
1115  {
1116  point ptt = AMITransforms()[i].transformPosition(pt);
1117  const vector ntt = AMITransforms()[i].transform(nt);
1118 
1119  const label nbrFacei =
1120  AMIs()[i].tgtPointFace(*this, neighbPatch(), ntt, facei, ptt);
1121 
1122  if (nbrFacei >= 0)
1123  {
1124  p = ptt;
1125  return labelPair(i, nbrFacei);
1126  }
1127  }
1128  }
1129  else
1130  {
1131  forAll(neighbPatch().AMIs(), i)
1132  {
1133  point ptt =
1134  neighbPatch().AMITransforms()[i].invTransformPosition(pt);
1135  const vector ntt =
1136  neighbPatch().AMITransforms()[i].invTransform(nt);
1137 
1138  const label nbrFacei =
1139  neighbPatch().AMIs()[i].srcPointFace
1140  (
1141  neighbPatch(),
1142  *this,
1143  ntt,
1144  facei,
1145  ptt
1146  );
1147 
1148  if (nbrFacei >= 0)
1149  {
1150  p = ptt;
1151  return labelPair(i, nbrFacei);
1152  }
1153  }
1154  }
1155 
1156  return labelPair(-1, -1);
1157 }
1158 
1159 
1161 {
1162  const cyclicAMIPolyPatch& patch = owner() ? *this : neighbPatch();
1163 
1164  const label proc = patch.AMIs()[0].singlePatchProc();
1165 
1166  for (label i = 1; i < patch.AMIs().size(); ++ i)
1167  {
1168  if (patch.AMIs()[i].singlePatchProc() != proc)
1169  {
1170  return -1;
1171  }
1172  }
1173 
1174  return proc;
1175 }
1176 
1177 
1179 {
1181  if (!nbrPatchName_.empty())
1182  {
1183  os.writeKeyword("neighbourPatch") << nbrPatchName_
1184  << token::END_STATEMENT << nl;
1185  }
1186  coupleGroup_.write(os);
1187 
1188  switch (transform())
1189  {
1190  case ROTATIONAL:
1191  {
1192  os.writeKeyword("rotationAxis") << rotationAxis_
1193  << token::END_STATEMENT << nl;
1194  os.writeKeyword("rotationCentre") << rotationCentre_
1195  << token::END_STATEMENT << nl;
1196 
1197  if (rotationAngleDefined_)
1198  {
1199  os.writeKeyword("rotationAngle") << radToDeg(rotationAngle_)
1200  << token::END_STATEMENT << nl;
1201  }
1202 
1203  break;
1204  }
1205  case TRANSLATIONAL:
1206  {
1207  os.writeKeyword("separationVector") << separationVector_
1208  << token::END_STATEMENT << nl;
1209  break;
1210  }
1211  case NOORDERING:
1212  {
1213  break;
1214  }
1215  default:
1216  {
1217  // No additional info to write
1218  }
1219  }
1220 
1221  if (AMIReverse_)
1222  {
1223  os.writeKeyword("flipNormals") << AMIReverse_
1224  << token::END_STATEMENT << nl;
1225  }
1226 
1227  if (AMILowWeightCorrection_ > 0)
1228  {
1229  os.writeKeyword("lowWeightCorrection") << AMILowWeightCorrection_
1230  << token::END_STATEMENT << nl;
1231  }
1232 
1233  os.writeKeyword("method")
1235  << token::END_STATEMENT << nl;
1236 
1237  if (!surfDict_.empty())
1238  {
1239  os.writeKeyword(surfDict_.dictName());
1240  os << surfDict_;
1241  }
1242 }
1243 
1244 
1245 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
static word interpolationMethodToWord(const interpolationMethod &method)
Convert interpolationMethod to word representation.
virtual void resetAMI() const
Reset the AMI interpolator.
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:115
vector separationVector_
Translation vector.
dimensionedScalar acos(const dimensionedScalar &ds)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
cyclicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN, const bool AMIRequireMatch=true, const AMIPatchToPatchInterpolation::interpolationMethod AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight)
Construct from (base couped patch) components.
fileName path() const
Return path.
Definition: Time.H:266
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
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()
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Output to file stream.
Definition: OFstream.H:82
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
uint8_t direction
Definition: direction.H:45
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
static interpolationMethod wordTointerpolationMethod(const word &method)
Convert word to interpolationMethod.
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:256
virtual const scalarField & neighbWeightsSum() const
Return the weights sum for the neighbour patch.
labelPair pointAMIAndFace(const label facei, const vector &n, point &p) const
Return the transform and face indices on neighbour patch which.
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.
const PtrList< AMIPatchToPatchInterpolation > & AMIs() const
Return a reference to the AMI interpolators.
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
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate field.
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 bool AMIRequireMatch_
Flag to indicate that patches should match/overlap.
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.
const List< vectorTensorTransform > & AMITransforms() const
Return a reference to the AMI transforms.
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
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
const word & neighbPatchName() const
Neighbour patch name.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
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.
label singlePatchProc() const
Index of processor that holds all of both sides, or -1 if.
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.
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
static const zero Zero
Definition: zero.H:97
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
virtual const scalarField & weightsSum() const
Return the weights sum for this patch.
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).
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
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 vectorTensorTransform I
static const Tensor I
Definition: Tensor.H:80
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:265
const AMIPatchToPatchInterpolation::interpolationMethod AMIMethod_
AMI Method.
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 R(A, B, C, D, E, F, K, M)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
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 class for managing temporary objects.
Definition: PtrList.H:53
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.
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