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-2019 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"
29 #include "unitConversion.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(cyclicAMIPolyPatch, 0);
37 
38  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, word);
39  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, dictionary);
40 }
41 
42 
43 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
44 
45 Foam::vector Foam::cyclicAMIPolyPatch::findFaceNormalMaxRadius
46 (
47  const pointField& faceCentres
48 ) const
49 {
50  // Determine a face furthest away from the axis
51 
52  const vectorField n((faceCentres - rotationCentre_) ^ rotationAxis_);
53 
54  const scalarField magRadSqr(magSqr(n));
55 
56  label facei = findMax(magRadSqr);
57 
58  if (debug)
59  {
60  Info<< "findFaceMaxRadius(const pointField&) : patch: " << name() << nl
61  << " rotFace = " << facei << nl
62  << " point = " << faceCentres[facei] << nl
63  << " distance = " << Foam::sqrt(magRadSqr[facei])
64  << endl;
65  }
66 
67  return n[facei];
68 }
69 
70 
72 (
73  const primitivePatch& half0,
74  const pointField& half0Ctrs,
75  const vectorField& half0Areas,
76  const pointField& half1Ctrs,
77  const vectorField& half1Areas
78 )
79 {
80  if (transform() != neighbPatch().transform())
81  {
83  << "Patch " << name()
84  << " has transform type " << transformTypeNames[transform()]
85  << ", neighbour patch " << neighbPatchName()
86  << " has transform type "
87  << neighbPatch().transformTypeNames[neighbPatch().transform()]
88  << exit(FatalError);
89  }
90 
91 
92  // Calculate transformation tensors
93 
94  switch (transform())
95  {
96  case ROTATIONAL:
97  {
98  tensor revT = Zero;
99 
100  if (rotationAngleDefined_)
101  {
102  const tensor T(rotationAxis_*rotationAxis_);
103 
104  const tensor S
105  (
106  0, -rotationAxis_.z(), rotationAxis_.y(),
107  rotationAxis_.z(), 0, -rotationAxis_.x(),
108  -rotationAxis_.y(), rotationAxis_.x(), 0
109  );
110 
111  const tensor revTPos
112  (
113  T
114  + cos(rotationAngle_)*(tensor::I - T)
115  + sin(rotationAngle_)*S
116  );
117 
118  const tensor revTNeg
119  (
120  T
121  + cos(-rotationAngle_)*(tensor::I - T)
122  + sin(-rotationAngle_)*S
123  );
124 
125  // Check - assume correct angle when difference in face areas
126  // is the smallest
127  const vector transformedAreaPos = gSum(half1Areas & revTPos);
128  const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
129  const vector area0 = gSum(half0Areas);
130  const scalar magArea0 = mag(area0) + rootVSmall;
131 
132  // Areas have opposite sign, so sum should be zero when correct
133  // rotation applied
134  const scalar errorPos = mag(transformedAreaPos + area0);
135  const scalar errorNeg = mag(transformedAreaNeg + area0);
136 
137  const scalar normErrorPos = errorPos/magArea0;
138  const scalar normErrorNeg = errorNeg/magArea0;
139 
140  if (errorPos > errorNeg && normErrorNeg < matchTolerance())
141  {
142  revT = revTNeg;
143  rotationAngle_ *= -1;
144  }
145  else
146  {
147  revT = revTPos;
148  }
149 
150  const scalar areaError = min(normErrorPos, normErrorNeg);
151 
152  if (areaError > matchTolerance())
153  {
155  << "Patch areas are not consistent within "
156  << 100*matchTolerance()
157  << " % indicating a possible error in the specified "
158  << "angle of rotation" << nl
159  << " owner patch : " << name() << nl
160  << " neighbour patch : " << neighbPatch().name()
161  << nl
162  << " angle : "
163  << radToDeg(rotationAngle_) << " deg" << nl
164  << " area error : " << 100*areaError << " %"
165  << " match tolerance : " << matchTolerance()
166  << endl;
167  }
168 
169  if (debug)
170  {
171  scalar theta = radToDeg(rotationAngle_);
172 
173  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
174  << name()
175  << " Specified rotation:"
176  << " swept angle: " << theta << " [deg]"
177  << " reverse transform: " << revT
178  << endl;
179  }
180  }
181  else
182  {
183  point n0 = Zero;
184  point n1 = Zero;
185  if (half0Ctrs.size())
186  {
187  n0 = findFaceNormalMaxRadius(half0Ctrs);
188  }
189  if (half1Ctrs.size())
190  {
191  n1 = -findFaceNormalMaxRadius(half1Ctrs);
192  }
193 
194  reduce(n0, maxMagSqrOp<point>());
195  reduce(n1, maxMagSqrOp<point>());
196 
197  n0 /= mag(n0) + vSmall;
198  n1 /= mag(n1) + vSmall;
199 
200  // Extended tensor from two local coordinate systems calculated
201  // using normal and rotation axis
202  const tensor E0
203  (
204  rotationAxis_,
205  (n0 ^ rotationAxis_),
206  n0
207  );
208  const tensor E1
209  (
210  rotationAxis_,
211  (-n1 ^ rotationAxis_),
212  -n1
213  );
214  revT = E1.T() & E0;
215 
216  if (debug)
217  {
218  scalar theta = radToDeg(acos(-(n0 & n1)));
219 
220  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
221  << name()
222  << " Specified rotation:"
223  << " n0:" << n0 << " n1:" << n1
224  << " swept angle: " << theta << " [deg]"
225  << " reverse transform: " << revT
226  << endl;
227  }
228  }
229 
230  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
231  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
232  const_cast<vectorField&>(separation()).setSize(0);
233  const_cast<boolList&>(collocated()) = boolList(1, false);
234 
235  break;
236  }
237  case TRANSLATIONAL:
238  {
239  if (debug)
240  {
241  Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
242  << " Specified translation : " << separationVector_
243  << endl;
244  }
245 
246  const_cast<tensorField&>(forwardT()).clear();
247  const_cast<tensorField&>(reverseT()).clear();
248  const_cast<vectorField&>(separation()) = vectorField
249  (
250  1,
251  separationVector_
252  );
253  const_cast<boolList&>(collocated()) = boolList(1, false);
254 
255  break;
256  }
257  default:
258  {
259  if (debug)
260  {
261  Pout<< "patch:" << name()
262  << " Assuming cyclic AMI pairs are colocated" << endl;
263  }
264 
265  const_cast<tensorField&>(forwardT()).clear();
266  const_cast<tensorField&>(reverseT()).clear();
267  const_cast<vectorField&>(separation()).setSize(0);
268  const_cast<boolList&>(collocated()) = boolList(1, true);
269 
270  break;
271  }
272  }
273 
274  if (debug)
275  {
276  Pout<< "patch: " << name() << nl
277  << " forwardT = " << forwardT() << nl
278  << " reverseT = " << reverseT() << nl
279  << " separation = " << separation() << nl
280  << " collocated = " << collocated() << nl << endl;
281  }
282 }
283 
284 
285 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
286 
288 {
289  if (owner())
290  {
291  const polyPatch& nbr = neighbPatch();
292  pointField nbrPoints
293  (
294  neighbPatch().boundaryMesh().mesh().points(),
295  neighbPatch().meshPoints()
296  );
297 
298  if (debug)
299  {
300  const Time& t = boundaryMesh().mesh().time();
301  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
302  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
303  }
304 
305  // Transform neighbour patch to local system
306  transformPosition(nbrPoints);
307  primitivePatch nbrPatch0
308  (
310  (
311  nbr.localFaces(),
312  nbr.size()
313  ),
314  nbrPoints
315  );
316 
317  if (debug)
318  {
319  const Time& t = boundaryMesh().mesh().time();
320  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
321  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
322 
323  OFstream osO(t.path()/name() + "_ownerPatch.obj");
324  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
325  }
326 
327  // Construct/apply AMI interpolation to determine addressing and weights
328  AMIs_.resize(1);
329  AMIs_.set
330  (
331  0,
332  new AMIInterpolation
333  (
334  *this,
335  nbrPatch0,
336  surfPtr(),
338  AMIRequireMatch_,
339  AMIMethod_,
340  AMILowWeightCorrection_,
341  AMIReverse_
342  )
343  );
344 
345  AMITransforms_.resize(1, vectorTensorTransform::I);
346 
347  if (debug)
348  {
349  Pout<< "cyclicAMIPolyPatch : " << name()
350  << " constructed AMI with " << nl
351  << " " << "srcAddress:" << AMIs_[0].srcAddress().size()
352  << nl
353  << " " << "tgAddress :" << AMIs_[0].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].area(half0.points());
367  }
368 
369  const cyclicAMIPolyPatch& half1 = neighbPatch();
370  vectorField half1Areas(half1.size());
371  forAll(half1, facei)
372  {
373  half1Areas[facei] = half1[facei].area(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 AMIs and transforms
399  AMIs_.clear();
400  AMITransforms_.clear();
401 
403 }
404 
405 
407 {
408  calcGeometry
409  (
410  *this,
411  faceCentres(),
412  faceAreas(),
413  faceCellCentres(),
414  neighbPatch().faceCentres(),
415  neighbPatch().faceAreas(),
416  neighbPatch().faceCellCentres()
417  );
418 }
419 
420 
422 (
423  PstreamBuffers& pBufs,
424  const pointField& p
425 )
426 {
427  // Clear the invalid AMIs and transforms
428  AMIs_.clear();
429  AMITransforms_.clear();
430 
431  polyPatch::initMovePoints(pBufs, p);
432 
433  // See below. Clear out any local geometry
435 }
436 
437 
439 (
440  PstreamBuffers& pBufs,
441  const pointField& p
442 )
443 {
444  polyPatch::movePoints(pBufs, p);
445 
446  calcTransforms();
447 }
448 
449 
451 {
452  // Clear the invalid AMIs and transforms
453  AMIs_.clear();
454  AMITransforms_.clear();
455 
457 }
458 
459 
461 {
462  polyPatch::updateMesh(pBufs);
463 }
464 
465 
467 {
468  // Clear the invalid AMIs and transforms
469  AMIs_.clear();
470  AMITransforms_.clear();
471 
473 }
474 
475 
476 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
477 
479 (
480  const word& name,
481  const label size,
482  const label start,
483  const label index,
484  const polyBoundaryMesh& bm,
485  const word& patchType,
486  const transformType transform,
487  const bool AMIRequireMatch,
489 )
490 :
491  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
492  nbrPatchName_(word::null),
493  nbrPatchID_(-1),
494  rotationAxis_(Zero),
495  rotationCentre_(point::zero),
496  rotationAngleDefined_(false),
497  rotationAngle_(0.0),
498  separationVector_(Zero),
499  AMIs_(),
500  AMITransforms_(),
501  AMIReverse_(false),
502  AMIRequireMatch_(AMIRequireMatch),
503  AMILowWeightCorrection_(-1.0),
504  AMIMethod_(AMIMethod),
505  surfPtr_(nullptr),
506  surfDict_(fileName("surface"))
507 {
508  // Neighbour patch might not be valid yet so no transformation
509  // calculation possible
510 }
511 
512 
514 (
515  const word& name,
516  const dictionary& dict,
517  const label index,
518  const polyBoundaryMesh& bm,
519  const word& patchType,
520  const bool AMIRequireMatch,
522 )
523 :
524  coupledPolyPatch(name, dict, index, bm, patchType),
525  nbrPatchName_(dict.lookupOrDefault<word>("neighbourPatch", "")),
526  coupleGroup_(dict),
527  nbrPatchID_(-1),
528  rotationAxis_(Zero),
529  rotationCentre_(point::zero),
530  rotationAngleDefined_(false),
531  rotationAngle_(0.0),
532  separationVector_(Zero),
533  AMIs_(),
534  AMITransforms_(),
535  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
536  AMIRequireMatch_(AMIRequireMatch),
537  AMILowWeightCorrection_(dict.lookupOrDefault("lowWeightCorrection", -1.0)),
538  AMIMethod_
539  (
540  dict.found("method")
542  (
543  dict.lookup("method")
544  )
545  : AMIMethod
546  ),
547  surfPtr_(nullptr),
548  surfDict_(dict.subOrEmptyDict("surface"))
549 {
550  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
551  {
553  (
554  dict
555  ) << "No \"neighbourPatch\" or \"coupleGroup\" provided."
556  << exit(FatalIOError);
557  }
558 
559  if (nbrPatchName_ == name)
560  {
562  (
563  dict
564  ) << "Neighbour patch name " << nbrPatchName_
565  << " cannot be the same as this patch " << name
566  << exit(FatalIOError);
567  }
568 
569  switch (transform())
570  {
571  case ROTATIONAL:
572  {
573  dict.lookup("rotationAxis") >> rotationAxis_;
574  dict.lookup("rotationCentre") >> rotationCentre_;
575  if (dict.readIfPresent("rotationAngle", rotationAngle_))
576  {
577  rotationAngleDefined_ = true;
578  rotationAngle_ = degToRad(rotationAngle_);
579 
580  if (debug)
581  {
582  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
583  << endl;
584  }
585  }
586 
587  scalar magRot = mag(rotationAxis_);
588  if (magRot < small)
589  {
591  (
592  dict
593  ) << "Illegal rotationAxis " << rotationAxis_ << endl
594  << "Please supply a non-zero vector."
595  << exit(FatalIOError);
596  }
597  rotationAxis_ /= magRot;
598 
599  break;
600  }
601  case TRANSLATIONAL:
602  {
603  dict.lookup("separationVector") >> separationVector_;
604  break;
605  }
606  default:
607  {
608  // No additional info required
609  }
610  }
611 
612  // Neighbour patch might not be valid yet so no transformation
613  // calculation possible
614 }
615 
616 
618 (
619  const cyclicAMIPolyPatch& pp,
620  const polyBoundaryMesh& bm
621 )
622 :
623  coupledPolyPatch(pp, bm),
624  nbrPatchName_(pp.nbrPatchName_),
625  coupleGroup_(pp.coupleGroup_),
626  nbrPatchID_(-1),
627  rotationAxis_(pp.rotationAxis_),
628  rotationCentre_(pp.rotationCentre_),
629  rotationAngleDefined_(pp.rotationAngleDefined_),
630  rotationAngle_(pp.rotationAngle_),
631  separationVector_(pp.separationVector_),
632  AMIs_(),
633  AMITransforms_(),
634  AMIReverse_(pp.AMIReverse_),
635  AMIRequireMatch_(pp.AMIRequireMatch_),
636  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
637  AMIMethod_(pp.AMIMethod_),
638  surfPtr_(nullptr),
639  surfDict_(pp.surfDict_)
640 {
641  // Neighbour patch might not be valid yet so no transformation
642  // calculation possible
643 }
644 
645 
647 (
648  const cyclicAMIPolyPatch& pp,
649  const polyBoundaryMesh& bm,
650  const label index,
651  const label newSize,
652  const label newStart,
653  const word& nbrPatchName
654 )
655 :
656  coupledPolyPatch(pp, bm, index, newSize, newStart),
657  nbrPatchName_(nbrPatchName),
658  coupleGroup_(pp.coupleGroup_),
659  nbrPatchID_(-1),
660  rotationAxis_(pp.rotationAxis_),
661  rotationCentre_(pp.rotationCentre_),
662  rotationAngleDefined_(pp.rotationAngleDefined_),
663  rotationAngle_(pp.rotationAngle_),
664  separationVector_(pp.separationVector_),
665  AMIs_(),
666  AMITransforms_(),
667  AMIReverse_(pp.AMIReverse_),
668  AMIRequireMatch_(pp.AMIRequireMatch_),
669  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
670  AMIMethod_(pp.AMIMethod_),
671  surfPtr_(nullptr),
672  surfDict_(pp.surfDict_)
673 {
674  if (nbrPatchName_ == name())
675  {
677  << "Neighbour patch name " << nbrPatchName_
678  << " cannot be the same as this patch " << name()
679  << exit(FatalError);
680  }
681 
682  // Neighbour patch might not be valid yet so no transformation
683  // calculation possible
684 }
685 
686 
688 (
689  const cyclicAMIPolyPatch& pp,
690  const polyBoundaryMesh& bm,
691  const label index,
692  const labelUList& mapAddressing,
693  const label newStart
694 )
695 :
696  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
697  nbrPatchName_(pp.nbrPatchName_),
698  coupleGroup_(pp.coupleGroup_),
699  nbrPatchID_(-1),
700  rotationAxis_(pp.rotationAxis_),
701  rotationCentre_(pp.rotationCentre_),
702  rotationAngleDefined_(pp.rotationAngleDefined_),
703  rotationAngle_(pp.rotationAngle_),
704  separationVector_(pp.separationVector_),
705  AMIs_(),
706  AMITransforms_(),
707  AMIReverse_(pp.AMIReverse_),
708  AMIRequireMatch_(pp.AMIRequireMatch_),
709  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
710  AMIMethod_(pp.AMIMethod_),
711  surfPtr_(nullptr),
712  surfDict_(pp.surfDict_)
713 {}
714 
715 
716 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
717 
719 {}
720 
721 
722 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
723 
725 {
726  if (nbrPatchID_ == -1)
727  {
728  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
729 
730  if (nbrPatchID_ == -1)
731  {
733  << "Illegal neighbourPatch name " << neighbPatchName()
734  << nl << "Valid patch names are "
735  << this->boundaryMesh().names()
736  << exit(FatalError);
737  }
738 
739  // Check that it is a cyclic AMI patch
740  const cyclicAMIPolyPatch& nbrPatch =
741  refCast<const cyclicAMIPolyPatch>
742  (
743  this->boundaryMesh()[nbrPatchID_]
744  );
745 
746  if (nbrPatch.neighbPatchName() != name())
747  {
749  << "Patch " << name()
750  << " specifies neighbour patch " << neighbPatchName()
751  << nl << " but that in return specifies "
752  << nbrPatch.neighbPatchName() << endl;
753  }
754  }
755 
756  return nbrPatchID_;
757 }
758 
759 
761 {
762  return index() < neighbPatchID();
763 }
764 
765 
767 {
768  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
769  return refCast<const cyclicAMIPolyPatch>(pp);
770 }
771 
772 
775 {
776  const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
777 
778  if (!surfPtr_.valid() && owner() && surfType != "none")
779  {
780  word surfName(surfDict_.lookupOrDefault("name", name()));
781 
782  const polyMesh& mesh = boundaryMesh().mesh();
783 
784  surfPtr_ =
786  (
787  surfType,
788  IOobject
789  (
790  surfName,
791  mesh.time().constant(),
792  "triSurface",
793  mesh,
796  ),
797  surfDict_
798  );
799  }
800 
801  return surfPtr_;
802 }
803 
804 
807 {
808  if (!owner())
809  {
811  << "AMI interpolators only available to owner patch"
812  << abort(FatalError);
813  }
814 
815  if (AMIs_.empty())
816  {
817  resetAMI();
818  }
819 
820  return AMIs_;
821 }
822 
823 
826 {
827  if (!owner())
828  {
830  << "AMI transforms only available to owner patch"
831  << abort(FatalError);
832  }
833 
834  if (AMIs_.empty())
835  {
836  resetAMI();
837  }
838 
839  return AMITransforms_;
840 }
841 
842 
844 {
845  if (owner())
846  {
847  return AMILowWeightCorrection_ > 0;
848  }
849  else
850  {
851  return neighbPatch().AMILowWeightCorrection_ > 0;
852  }
853 }
854 
855 
857 {
858  if (owner())
859  {
860  return AMIs()[0].srcWeightsSum();
861  }
862  else
863  {
864  return neighbPatch().AMIs()[0].tgtWeightsSum();
865  }
866 }
867 
868 
870 {
871  if (owner())
872  {
873  return AMIs()[0].tgtWeightsSum();
874  }
875  else
876  {
877  return neighbPatch().AMIs()[0].srcWeightsSum();
878  }
879 }
880 
881 
883 {
884  if (!parallel())
885  {
886  if (transform() == ROTATIONAL)
887  {
888  l = Foam::transform(forwardT(), l - rotationCentre_)
889  + rotationCentre_;
890  }
891  else
892  {
893  l = Foam::transform(forwardT(), l);
894  }
895  }
896  else if (separated())
897  {
898  // transformPosition gets called on the receiving side,
899  // separation gets calculated on the sending side so subtract
900 
901  const vectorField& s = separation();
902  if (s.size() == 1)
903  {
904  forAll(l, i)
905  {
906  l[i] -= s[0];
907  }
908  }
909  else
910  {
911  l -= s;
912  }
913  }
914 }
915 
916 
918 (
919  point& l,
920  const label facei
921 ) const
922 {
923  if (!parallel())
924  {
925  const tensor& T =
926  (
927  forwardT().size() == 1
928  ? forwardT()[0]
929  : forwardT()[facei]
930  );
931 
932  if (transform() == ROTATIONAL)
933  {
934  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
935  }
936  else
937  {
938  l = Foam::transform(T, l);
939  }
940  }
941  else if (separated())
942  {
943  const vector& s =
944  (
945  separation().size() == 1
946  ? separation()[0]
947  : separation()[facei]
948  );
949 
950  l -= s;
951  }
952 }
953 
954 
956 (
957  vector& d,
958  const label facei
959 ) const
960 {
961  if (!parallel())
962  {
963  const tensor& T =
964  (
965  forwardT().size() == 1
966  ? forwardT()[0]
967  : forwardT()[facei]
968  );
969 
970  d = Foam::transform(T, d);
971  }
972 }
973 
974 
976 (
977  point& l,
978  const label facei
979 ) const
980 {
981  if (!parallel())
982  {
983  const tensor& T =
984  (
985  reverseT().size() == 1
986  ? reverseT()[0]
987  : reverseT()[facei]
988  );
989 
990  if (transform() == ROTATIONAL)
991  {
992  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
993  }
994  else
995  {
996  l = Foam::transform(T, l);
997  }
998  }
999  else if (separated())
1000  {
1001  const vector& s =
1002  (
1003  separation().size() == 1
1004  ? separation()[0]
1005  : separation()[facei]
1006  );
1007 
1008  l += s;
1009  }
1010 }
1011 
1012 
1015  vector& d,
1016  const label facei
1017 ) const
1018 {
1019  if (!parallel())
1020  {
1021  const tensor& T =
1022  (
1023  reverseT().size() == 1
1024  ? reverseT()[0]
1025  : reverseT()[facei]
1026  );
1027 
1028  d = Foam::transform(T, d);
1029  }
1030 }
1031 
1032 
1035  const scalarField& fld,
1036  const direction cmpt,
1037  const direction rank,
1038  const scalarUList& defaultValues
1039 ) const
1040 {
1041  const cyclicAMIPolyPatch& nei = neighbPatch();
1042 
1043  tmp<scalarField> result(new scalarField(size(), Zero));
1044 
1045  if (owner())
1046  {
1047  forAll(AMIs(), i)
1048  {
1049  const scalar r =
1050  pow(inv(AMITransforms()[i]).R()(cmpt, cmpt), rank);
1051 
1052  result.ref() +=
1053  AMIs()[i].interpolateToSource(r*fld, defaultValues);
1054  }
1055  }
1056  else
1057  {
1058  forAll(nei.AMIs(), i)
1059  {
1060  const scalar r =
1061  pow(nei.AMITransforms()[i].R()(cmpt, cmpt), rank);
1062 
1063  result.ref() +=
1064  nei.AMIs()[i].interpolateToTarget(r*fld, defaultValues);
1065  }
1066  }
1067 
1068  return result;
1069 }
1070 
1071 
1074  const primitivePatch& referPatch,
1075  const pointField& thisCtrs,
1076  const vectorField& thisAreas,
1077  const pointField& thisCc,
1078  const pointField& nbrCtrs,
1079  const vectorField& nbrAreas,
1080  const pointField& nbrCc
1081 )
1082 {
1083  calcTransforms
1084  (
1085  referPatch,
1086  thisCtrs,
1087  thisAreas,
1088  nbrCtrs,
1089  nbrAreas
1090  );
1091 }
1092 
1093 
1096  PstreamBuffers& pBufs,
1097  const primitivePatch& pp
1098 ) const
1099 {}
1100 
1101 
1104  PstreamBuffers& pBufs,
1105  const primitivePatch& pp,
1106  labelList& faceMap,
1107  labelList& rotation
1108 ) const
1109 {
1110  faceMap.setSize(pp.size());
1111  faceMap = -1;
1112 
1113  rotation.setSize(pp.size());
1114  rotation = 0;
1115 
1116  return false;
1117 }
1118 
1119 
1122  const label facei,
1123  const vector& n,
1124  point& p
1125 ) const
1126 {
1127  point pt(p);
1128  reverseTransformPosition(pt, facei);
1129 
1130  vector nt(n);
1131  reverseTransformDirection(nt, facei);
1132 
1133  if (owner())
1134  {
1135  forAll(AMIs(), i)
1136  {
1137  point ptt = AMITransforms()[i].transformPosition(pt);
1138  const vector ntt = AMITransforms()[i].transform(nt);
1139 
1140  const label nbrFacei =
1141  AMIs()[i].tgtPointFace(*this, neighbPatch(), ntt, facei, ptt);
1142 
1143  if (nbrFacei >= 0)
1144  {
1145  p = ptt;
1146  return labelPair(i, nbrFacei);
1147  }
1148  }
1149  }
1150  else
1151  {
1152  forAll(neighbPatch().AMIs(), i)
1153  {
1154  point ptt =
1155  neighbPatch().AMITransforms()[i].invTransformPosition(pt);
1156  const vector ntt =
1157  neighbPatch().AMITransforms()[i].invTransform(nt);
1158 
1159  const label nbrFacei =
1160  neighbPatch().AMIs()[i].srcPointFace
1161  (
1162  neighbPatch(),
1163  *this,
1164  ntt,
1165  facei,
1166  ptt
1167  );
1168 
1169  if (nbrFacei >= 0)
1170  {
1171  p = ptt;
1172  return labelPair(i, nbrFacei);
1173  }
1174  }
1175  }
1176 
1177  return labelPair(-1, -1);
1178 }
1179 
1180 
1182 {
1183  const cyclicAMIPolyPatch& patch = owner() ? *this : neighbPatch();
1184 
1185  const label proc = patch.AMIs()[0].singlePatchProc();
1186 
1187  for (label i = 1; i < patch.AMIs().size(); ++ i)
1188  {
1189  if (patch.AMIs()[i].singlePatchProc() != proc)
1190  {
1191  return -1;
1192  }
1193  }
1194 
1195  return proc;
1196 }
1197 
1198 
1200 {
1202  if (!nbrPatchName_.empty())
1203  {
1204  writeEntry(os, "neighbourPatch", nbrPatchName_);
1205  }
1206  coupleGroup_.write(os);
1207 
1208  switch (transform())
1209  {
1210  case ROTATIONAL:
1211  {
1212  writeEntry(os, "rotationAxis", rotationAxis_);
1213  writeEntry(os, "rotationCentre", rotationCentre_);
1214 
1215  if (rotationAngleDefined_)
1216  {
1217  writeEntry(os, "rotationAngle", radToDeg(rotationAngle_));
1218  }
1219 
1220  break;
1221  }
1222  case TRANSLATIONAL:
1223  {
1224  writeEntry(os, "separationVector", separationVector_);
1225  break;
1226  }
1227  case NOORDERING:
1228  {
1229  break;
1230  }
1231  default:
1232  {
1233  // No additional info to write
1234  }
1235  }
1236 
1237  if (AMIReverse_)
1238  {
1239  writeEntry(os, "flipNormals", AMIReverse_);
1240  }
1241 
1242  if (AMILowWeightCorrection_ > 0)
1243  {
1244  writeEntry(os, "lowWeightCorrection", AMILowWeightCorrection_);
1245  }
1246 
1247  writeEntry
1248  (
1249  os,
1250  "method",
1252  );
1253 
1254  if (!surfDict_.empty())
1255  {
1256  os.writeKeyword(surfDict_.dictName());
1257  os << surfDict_;
1258  }
1259 }
1260 
1261 
1262 // ************************************************************************* //
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:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:79
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:158
#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
Unit conversion functions.
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
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.
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 AMIInterpolation::interpolationMethod AMIMethod=AMIInterpolation::imFaceAreaWeight)
Construct from (base couped patch) components.
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
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< 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
const AMIInterpolation::interpolationMethod AMIMethod_
AMI Method.
const PtrList< AMIInterpolation > & AMIs() const
Return a reference to the AMI interpolators.
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:60
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
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)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PrimitivePatch< SubList< face >, 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
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
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.
static interpolationMethod wordTointerpolationMethod(const word &method)
Convert word to interpolationMethod.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const bMesh & mesh() const
Definition: boundaryMesh.H:199
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:734
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
virtual void transformDirection(vector &d, const label facei) const
Transform a patch-based direction from nbr side to this side.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
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