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 
70 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
71 
73 {
74  const cyclicAMIPolyPatch& half0 = *this;
75  vectorField half0Areas(half0.size());
76  forAll(half0, facei)
77  {
78  half0Areas[facei] = half0[facei].normal(half0.points());
79  }
80 
81  const cyclicAMIPolyPatch& half1 = neighbPatch();
82  vectorField half1Areas(half1.size());
83  forAll(half1, facei)
84  {
85  half1Areas[facei] = half1[facei].normal(half1.points());
86  }
87 
88  calcTransforms
89  (
90  half0,
91  half0.faceCentres(),
92  half0Areas,
93  half1.faceCentres(),
94  half1Areas
95  );
96 
97  if (debug)
98  {
99  Pout<< "calcTransforms() : patch: " << name() << nl
100  << " forwardT = " << forwardT() << nl
101  << " reverseT = " << reverseT() << nl
102  << " separation = " << separation() << nl
103  << " collocated = " << collocated() << nl << endl;
104  }
105 }
106 
107 
109 (
110  const primitivePatch& half0,
111  const pointField& half0Ctrs,
112  const vectorField& half0Areas,
113  const pointField& half1Ctrs,
114  const vectorField& half1Areas
115 )
116 {
117  if (transform() != neighbPatch().transform())
118  {
120  << "Patch " << name()
121  << " has transform type " << transformTypeNames[transform()]
122  << ", neighbour patch " << neighbPatchName()
123  << " has transform type "
124  << neighbPatch().transformTypeNames[neighbPatch().transform()]
125  << exit(FatalError);
126  }
127 
128 
129  // Calculate transformation tensors
130 
131  switch (transform())
132  {
133  case ROTATIONAL:
134  {
135  tensor revT = Zero;
136 
137  if (rotationAngleDefined_)
138  {
139  const tensor T(rotationAxis_*rotationAxis_);
140 
141  const tensor S
142  (
143  0, -rotationAxis_.z(), rotationAxis_.y(),
144  rotationAxis_.z(), 0, -rotationAxis_.x(),
145  -rotationAxis_.y(), rotationAxis_.x(), 0
146  );
147 
148  const tensor revTPos
149  (
150  T
151  + cos(rotationAngle_)*(tensor::I - T)
152  + sin(rotationAngle_)*S
153  );
154 
155  const tensor revTNeg
156  (
157  T
158  + cos(-rotationAngle_)*(tensor::I - T)
159  + sin(-rotationAngle_)*S
160  );
161 
162  // Check - assume correct angle when difference in face areas
163  // is the smallest
164  const vector transformedAreaPos = gSum(half1Areas & revTPos);
165  const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
166  const vector area0 = gSum(half0Areas);
167  const scalar magArea0 = mag(area0) + ROOTVSMALL;
168 
169  // Areas have opposite sign, so sum should be zero when correct
170  // rotation applied
171  const scalar errorPos = mag(transformedAreaPos + area0);
172  const scalar errorNeg = mag(transformedAreaNeg + area0);
173 
174  const scalar normErrorPos = errorPos/magArea0;
175  const scalar normErrorNeg = errorNeg/magArea0;
176 
177  if (errorPos > errorNeg && normErrorNeg < matchTolerance())
178  {
179  revT = revTNeg;
180  rotationAngle_ *= -1;
181  }
182  else
183  {
184  revT = revTPos;
185  }
186 
187  const scalar areaError = min(normErrorPos, normErrorNeg);
188 
189  if (areaError > matchTolerance())
190  {
192  << "Patch areas are not consistent within "
193  << 100*matchTolerance()
194  << " % indicating a possible error in the specified "
195  << "angle of rotation" << nl
196  << " owner patch : " << name() << nl
197  << " neighbour patch : " << neighbPatch().name()
198  << nl
199  << " angle : "
200  << radToDeg(rotationAngle_) << " deg" << nl
201  << " area error : " << 100*areaError << " %"
202  << " match tolerance : " << matchTolerance()
203  << endl;
204  }
205 
206  if (debug)
207  {
208  scalar theta = radToDeg(rotationAngle_);
209 
210  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
211  << name()
212  << " Specified rotation:"
213  << " swept angle: " << theta << " [deg]"
214  << " reverse transform: " << revT
215  << endl;
216  }
217  }
218  else
219  {
220  point n0 = Zero;
221  point n1 = Zero;
222  if (half0Ctrs.size())
223  {
224  n0 = findFaceNormalMaxRadius(half0Ctrs);
225  }
226  if (half1Ctrs.size())
227  {
228  n1 = -findFaceNormalMaxRadius(half1Ctrs);
229  }
230 
231  reduce(n0, maxMagSqrOp<point>());
232  reduce(n1, maxMagSqrOp<point>());
233 
234  n0 /= mag(n0) + VSMALL;
235  n1 /= mag(n1) + VSMALL;
236 
237  // Extended tensor from two local coordinate systems calculated
238  // using normal and rotation axis
239  const tensor E0
240  (
241  rotationAxis_,
242  (n0 ^ rotationAxis_),
243  n0
244  );
245  const tensor E1
246  (
247  rotationAxis_,
248  (-n1 ^ rotationAxis_),
249  -n1
250  );
251  revT = E1.T() & E0;
252 
253  if (debug)
254  {
255  scalar theta = radToDeg(acos(-(n0 & n1)));
256 
257  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
258  << name()
259  << " Specified rotation:"
260  << " n0:" << n0 << " n1:" << n1
261  << " swept angle: " << theta << " [deg]"
262  << " reverse transform: " << revT
263  << endl;
264  }
265  }
266 
267  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
268  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
269  const_cast<vectorField&>(separation()).setSize(0);
270  const_cast<boolList&>(collocated()) = boolList(1, false);
271 
272  break;
273  }
274  case TRANSLATIONAL:
275  {
276  if (debug)
277  {
278  Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
279  << " Specified translation : " << separationVector_
280  << endl;
281  }
282 
283  const_cast<tensorField&>(forwardT()).clear();
284  const_cast<tensorField&>(reverseT()).clear();
285  const_cast<vectorField&>(separation()) = vectorField
286  (
287  1,
288  separationVector_
289  );
290  const_cast<boolList&>(collocated()) = boolList(1, false);
291 
292  break;
293  }
294  default:
295  {
296  if (debug)
297  {
298  Pout<< "patch:" << name()
299  << " Assuming cyclic AMI pairs are colocated" << endl;
300  }
301 
302  const_cast<tensorField&>(forwardT()).clear();
303  const_cast<tensorField&>(reverseT()).clear();
304  const_cast<vectorField&>(separation()).setSize(0);
305  const_cast<boolList&>(collocated()) = boolList(1, true);
306 
307  break;
308  }
309  }
310 
311  if (debug)
312  {
313  Pout<< "patch: " << name() << nl
314  << " forwardT = " << forwardT() << nl
315  << " reverseT = " << reverseT() << nl
316  << " separation = " << separation() << nl
317  << " collocated = " << collocated() << nl << endl;
318  }
319 }
320 
321 
322 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
323 
325 (
327 ) const
328 {
329  if (owner())
330  {
331  AMIPtr_.clear();
332 
333  const polyPatch& nbr = neighbPatch();
334  pointField nbrPoints
335  (
336  neighbPatch().boundaryMesh().mesh().points(),
337  neighbPatch().meshPoints()
338  );
339 
340  if (debug)
341  {
342  const Time& t = boundaryMesh().mesh().time();
343  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
344  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
345  }
346 
347  // Transform neighbour patch to local system
348  transformPosition(nbrPoints);
349  primitivePatch nbrPatch0
350  (
352  (
353  nbr.localFaces(),
354  nbr.size()
355  ),
356  nbrPoints
357  );
358 
359  if (debug)
360  {
361  const Time& t = boundaryMesh().mesh().time();
362  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
363  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
364 
365  OFstream osO(t.path()/name() + "_ownerPatch.obj");
366  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
367  }
368 
369  // Construct/apply AMI interpolation to determine addressing and weights
370  AMIPtr_.reset
371  (
373  (
374  *this,
375  nbrPatch0,
376  surfPtr(),
378  AMIRequireMatch_,
379  AMIMethod,
380  AMILowWeightCorrection_,
381  AMIReverse_
382  )
383  );
384 
385  if (debug)
386  {
387  Pout<< "cyclicAMIPolyPatch : " << name()
388  << " constructed AMI with " << nl
389  << " " << "srcAddress:" << AMIPtr_().srcAddress().size()
390  << nl
391  << " " << "tgAddress :" << AMIPtr_().tgtAddress().size()
392  << nl << endl;
393  }
394  }
395 }
396 
397 
399 {
400  // Clear the invalid AMI
401  AMIPtr_.clear();
402 
404 }
405 
406 
408 {
409  calcGeometry
410  (
411  *this,
412  faceCentres(),
413  faceAreas(),
414  faceCellCentres(),
415  neighbPatch().faceCentres(),
416  neighbPatch().faceAreas(),
417  neighbPatch().faceCellCentres()
418  );
419 }
420 
421 
423 (
424  PstreamBuffers& pBufs,
425  const pointField& p
426 )
427 {
428  // Clear the invalid AMI
429  AMIPtr_.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 AMI
453  AMIPtr_.clear();
454 
456 }
457 
458 
460 {
461  polyPatch::updateMesh(pBufs);
462 }
463 
464 
466 {
467  AMIPtr_.clear();
469 }
470 
471 
472 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
473 
475 (
476  const word& name,
477  const label size,
478  const label start,
479  const label index,
480  const polyBoundaryMesh& bm,
481  const word& patchType,
482  const transformType transform
483 )
484 :
485  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
486  nbrPatchName_(word::null),
487  nbrPatchID_(-1),
488  rotationAxis_(Zero),
489  rotationCentre_(point::zero),
490  rotationAngleDefined_(false),
491  rotationAngle_(0.0),
492  separationVector_(Zero),
493  AMIPtr_(NULL),
494  AMIReverse_(false),
495  AMIRequireMatch_(true),
496  AMILowWeightCorrection_(-1.0),
497  surfPtr_(NULL),
498  surfDict_(fileName("surface"))
499 {
500  // Neighbour patch might not be valid yet so no transformation
501  // calculation possible
502 }
503 
504 
506 (
507  const word& name,
508  const dictionary& dict,
509  const label index,
510  const polyBoundaryMesh& bm,
511  const word& patchType
512 )
513 :
514  coupledPolyPatch(name, dict, index, bm, patchType),
515  nbrPatchName_(dict.lookupOrDefault<word>("neighbourPatch", "")),
516  coupleGroup_(dict),
517  nbrPatchID_(-1),
518  rotationAxis_(Zero),
519  rotationCentre_(point::zero),
520  rotationAngleDefined_(false),
521  rotationAngle_(0.0),
522  separationVector_(Zero),
523  AMIPtr_(NULL),
524  AMIReverse_(dict.lookupOrDefault<bool>("flipNormals", false)),
525  AMIRequireMatch_(true),
526  AMILowWeightCorrection_(dict.lookupOrDefault("lowWeightCorrection", -1.0)),
527  surfPtr_(NULL),
528  surfDict_(dict.subOrEmptyDict("surface"))
529 {
530  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
531  {
533  (
534  dict
535  ) << "No \"neighbourPatch\" or \"coupleGroup\" provided."
536  << exit(FatalIOError);
537  }
538 
539  if (nbrPatchName_ == name)
540  {
542  (
543  dict
544  ) << "Neighbour patch name " << nbrPatchName_
545  << " cannot be the same as this patch " << name
546  << exit(FatalIOError);
547  }
548 
549  switch (transform())
550  {
551  case ROTATIONAL:
552  {
553  dict.lookup("rotationAxis") >> rotationAxis_;
554  dict.lookup("rotationCentre") >> rotationCentre_;
555  if (dict.readIfPresent("rotationAngle", rotationAngle_))
556  {
557  rotationAngleDefined_ = true;
558  rotationAngle_ = degToRad(rotationAngle_);
559 
560  if (debug)
561  {
562  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
563  << endl;
564  }
565  }
566 
567  scalar magRot = mag(rotationAxis_);
568  if (magRot < SMALL)
569  {
571  (
572  dict
573  ) << "Illegal rotationAxis " << rotationAxis_ << endl
574  << "Please supply a non-zero vector."
575  << exit(FatalIOError);
576  }
577  rotationAxis_ /= magRot;
578 
579  break;
580  }
581  case TRANSLATIONAL:
582  {
583  dict.lookup("separationVector") >> separationVector_;
584  break;
585  }
586  default:
587  {
588  // No additional info required
589  }
590  }
591 
592  // Neighbour patch might not be valid yet so no transformation
593  // calculation possible
594 }
595 
596 
598 (
599  const cyclicAMIPolyPatch& pp,
600  const polyBoundaryMesh& bm
601 )
602 :
603  coupledPolyPatch(pp, bm),
604  nbrPatchName_(pp.nbrPatchName_),
605  coupleGroup_(pp.coupleGroup_),
606  nbrPatchID_(-1),
607  rotationAxis_(pp.rotationAxis_),
608  rotationCentre_(pp.rotationCentre_),
609  rotationAngleDefined_(pp.rotationAngleDefined_),
610  rotationAngle_(pp.rotationAngle_),
611  separationVector_(pp.separationVector_),
612  AMIPtr_(NULL),
613  AMIReverse_(pp.AMIReverse_),
614  AMIRequireMatch_(pp.AMIRequireMatch_),
615  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
616  surfPtr_(NULL),
617  surfDict_(pp.surfDict_)
618 {
619  // Neighbour patch might not be valid yet so no transformation
620  // calculation possible
621 }
622 
623 
625 (
626  const cyclicAMIPolyPatch& pp,
627  const polyBoundaryMesh& bm,
628  const label index,
629  const label newSize,
630  const label newStart,
631  const word& nbrPatchName
632 )
633 :
634  coupledPolyPatch(pp, bm, index, newSize, newStart),
635  nbrPatchName_(nbrPatchName),
636  coupleGroup_(pp.coupleGroup_),
637  nbrPatchID_(-1),
638  rotationAxis_(pp.rotationAxis_),
639  rotationCentre_(pp.rotationCentre_),
640  rotationAngleDefined_(pp.rotationAngleDefined_),
641  rotationAngle_(pp.rotationAngle_),
642  separationVector_(pp.separationVector_),
643  AMIPtr_(NULL),
644  AMIReverse_(pp.AMIReverse_),
645  AMIRequireMatch_(pp.AMIRequireMatch_),
646  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
647  surfPtr_(NULL),
648  surfDict_(pp.surfDict_)
649 {
650  if (nbrPatchName_ == name())
651  {
653  << "Neighbour patch name " << nbrPatchName_
654  << " cannot be the same as this patch " << name()
655  << exit(FatalError);
656  }
657 
658  // Neighbour patch might not be valid yet so no transformation
659  // calculation possible
660 }
661 
662 
664 (
665  const cyclicAMIPolyPatch& pp,
666  const polyBoundaryMesh& bm,
667  const label index,
668  const labelUList& mapAddressing,
669  const label newStart
670 )
671 :
672  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
673  nbrPatchName_(pp.nbrPatchName_),
674  coupleGroup_(pp.coupleGroup_),
675  nbrPatchID_(-1),
676  rotationAxis_(pp.rotationAxis_),
677  rotationCentre_(pp.rotationCentre_),
678  rotationAngleDefined_(pp.rotationAngleDefined_),
679  rotationAngle_(pp.rotationAngle_),
680  separationVector_(pp.separationVector_),
681  AMIPtr_(NULL),
682  AMIReverse_(pp.AMIReverse_),
683  AMIRequireMatch_(pp.AMIRequireMatch_),
684  AMILowWeightCorrection_(pp.AMILowWeightCorrection_),
685  surfPtr_(NULL),
686  surfDict_(pp.surfDict_)
687 {}
688 
689 
690 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
691 
693 {}
694 
695 
696 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
697 
699 {
700  if (nbrPatchID_ == -1)
701  {
702  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
703 
704  if (nbrPatchID_ == -1)
705  {
707  << "Illegal neighbourPatch name " << neighbPatchName()
708  << nl << "Valid patch names are "
709  << this->boundaryMesh().names()
710  << exit(FatalError);
711  }
712 
713  // Check that it is a cyclic AMI patch
714  const cyclicAMIPolyPatch& nbrPatch =
715  refCast<const cyclicAMIPolyPatch>
716  (
717  this->boundaryMesh()[nbrPatchID_]
718  );
719 
720  if (nbrPatch.neighbPatchName() != name())
721  {
723  << "Patch " << name()
724  << " specifies neighbour patch " << neighbPatchName()
725  << nl << " but that in return specifies "
726  << nbrPatch.neighbPatchName() << endl;
727  }
728  }
729 
730  return nbrPatchID_;
731 }
732 
733 
735 {
736  return index() < neighbPatchID();
737 }
738 
739 
741 {
742  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
743  return refCast<const cyclicAMIPolyPatch>(pp);
744 }
745 
746 
749 {
750  const word surfType(surfDict_.lookupOrDefault<word>("type", "none"));
751 
752  if (!surfPtr_.valid() && owner() && surfType != "none")
753  {
754  word surfName(surfDict_.lookupOrDefault("name", name()));
755 
756  const polyMesh& mesh = boundaryMesh().mesh();
757 
758  surfPtr_ =
760  (
761  surfType,
762  IOobject
763  (
764  surfName,
765  mesh.time().constant(),
766  "triSurface",
767  mesh,
770  ),
771  surfDict_
772  );
773  }
774 
775  return surfPtr_;
776 }
777 
778 
780 {
781  if (!owner())
782  {
784  << "AMI interpolator only available to owner patch"
785  << abort(FatalError);
786  }
787 
788  if (!AMIPtr_.valid())
789  {
790  resetAMI();
791  }
792 
793  return AMIPtr_();
794 }
795 
796 
798 {
799  if (owner())
800  {
801  return AMI().applyLowWeightCorrection();
802  }
803  else
804  {
805  return neighbPatch().AMI().applyLowWeightCorrection();
806  }
807 }
808 
809 
811 {
812  if (!parallel())
813  {
814  if (transform() == ROTATIONAL)
815  {
816  l = Foam::transform(forwardT(), l - rotationCentre_)
817  + rotationCentre_;
818  }
819  else
820  {
821  l = Foam::transform(forwardT(), l);
822  }
823  }
824  else if (separated())
825  {
826  // transformPosition gets called on the receiving side,
827  // separation gets calculated on the sending side so subtract
828 
829  const vectorField& s = separation();
830  if (s.size() == 1)
831  {
832  forAll(l, i)
833  {
834  l[i] -= s[0];
835  }
836  }
837  else
838  {
839  l -= s;
840  }
841  }
842 }
843 
844 
846 (
847  point& l,
848  const label facei
849 ) const
850 {
851  if (!parallel())
852  {
853  const tensor& T =
854  (
855  forwardT().size() == 1
856  ? forwardT()[0]
857  : forwardT()[facei]
858  );
859 
860  if (transform() == ROTATIONAL)
861  {
862  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
863  }
864  else
865  {
866  l = Foam::transform(T, l);
867  }
868  }
869  else if (separated())
870  {
871  const vector& s =
872  (
873  separation().size() == 1
874  ? separation()[0]
875  : separation()[facei]
876  );
877 
878  l -= s;
879  }
880 }
881 
882 
884 (
885  point& l,
886  const label facei
887 ) const
888 {
889  if (!parallel())
890  {
891  const tensor& T =
892  (
893  reverseT().size() == 1
894  ? reverseT()[0]
895  : reverseT()[facei]
896  );
897 
898  if (transform() == ROTATIONAL)
899  {
900  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
901  }
902  else
903  {
904  l = Foam::transform(T, l);
905  }
906  }
907  else if (separated())
908  {
909  const vector& s =
910  (
911  separation().size() == 1
912  ? separation()[0]
913  : separation()[facei]
914  );
915 
916  l += s;
917  }
918 }
919 
920 
922 (
923  vector& d,
924  const label facei
925 ) const
926 {
927  if (!parallel())
928  {
929  const tensor& T =
930  (
931  reverseT().size() == 1
932  ? reverseT()[0]
933  : reverseT()[facei]
934  );
935 
936  d = Foam::transform(T, d);
937  }
938 }
939 
940 
942 (
943  const primitivePatch& referPatch,
944  const pointField& thisCtrs,
945  const vectorField& thisAreas,
946  const pointField& thisCc,
947  const pointField& nbrCtrs,
948  const vectorField& nbrAreas,
949  const pointField& nbrCc
950 )
951 {
952  calcTransforms
953  (
954  referPatch,
955  thisCtrs,
956  thisAreas,
957  nbrCtrs,
958  nbrAreas
959  );
960 }
961 
962 
964 (
965  PstreamBuffers& pBufs,
966  const primitivePatch& pp
967 ) const
968 {}
969 
970 
972 (
973  PstreamBuffers& pBufs,
974  const primitivePatch& pp,
975  labelList& faceMap,
976  labelList& rotation
977 ) const
978 {
979  faceMap.setSize(pp.size());
980  faceMap = -1;
981 
982  rotation.setSize(pp.size());
983  rotation = 0;
984 
985  return false;
986 }
987 
988 
990 (
991  const label facei,
992  const vector& n,
993  point& p
994 ) const
995 {
996  point prt(p);
997  reverseTransformPosition(prt, facei);
998 
999  vector nrt(n);
1000  reverseTransformDirection(nrt, facei);
1001 
1002  label nbrFacei = -1;
1003 
1004  if (owner())
1005  {
1006  nbrFacei = AMI().tgtPointFace
1007  (
1008  *this,
1009  neighbPatch(),
1010  nrt,
1011  facei,
1012  prt
1013  );
1014  }
1015  else
1016  {
1017  nbrFacei = neighbPatch().AMI().srcPointFace
1018  (
1019  neighbPatch(),
1020  *this,
1021  nrt,
1022  facei,
1023  prt
1024  );
1025  }
1026 
1027  if (nbrFacei >= 0)
1028  {
1029  p = prt;
1030  }
1031 
1032  return nbrFacei;
1033 }
1034 
1035 
1037 {
1039  if (!nbrPatchName_.empty())
1040  {
1041  os.writeKeyword("neighbourPatch") << nbrPatchName_
1042  << token::END_STATEMENT << nl;
1043  }
1044  coupleGroup_.write(os);
1045 
1046  switch (transform())
1047  {
1048  case ROTATIONAL:
1049  {
1050  os.writeKeyword("rotationAxis") << rotationAxis_
1051  << token::END_STATEMENT << nl;
1052  os.writeKeyword("rotationCentre") << rotationCentre_
1053  << token::END_STATEMENT << nl;
1054 
1055  if (rotationAngleDefined_)
1056  {
1057  os.writeKeyword("rotationAngle") << radToDeg(rotationAngle_)
1058  << token::END_STATEMENT << nl;
1059  }
1060 
1061  break;
1062  }
1063  case TRANSLATIONAL:
1064  {
1065  os.writeKeyword("separationVector") << separationVector_
1066  << token::END_STATEMENT << nl;
1067  break;
1068  }
1069  case NOORDERING:
1070  {
1071  break;
1072  }
1073  default:
1074  {
1075  // No additional info to write
1076  }
1077  }
1078 
1079  if (AMIReverse_)
1080  {
1081  os.writeKeyword("flipNormals") << AMIReverse_
1082  << token::END_STATEMENT << nl;
1083  }
1084 
1085  if (AMILowWeightCorrection_ > 0)
1086  {
1087  os.writeKeyword("lowWeightCorrection") << AMILowWeightCorrection_
1088  << token::END_STATEMENT << nl;
1089  }
1090 
1091  if (!surfDict_.empty())
1092  {
1093  os.writeKeyword(surfDict_.dictName());
1094  os << surfDict_;
1095  }
1096 }
1097 
1098 
1099 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:117
vector separationVector_
Translation vector.
dimensionedScalar acos(const dimensionedScalar &ds)
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
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()
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:668
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Output to file stream.
Definition: OFstream.H:81
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
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
const Field< PointType > & points() const
Return reference to global points.
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 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
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.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
const bool AMIReverse_
Flag to indicate that slave patch should be reversed for AMI.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar cos(const dimensionedScalar &ds)
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const scalar AMILowWeightCorrection_
Low weight correction threshold for AMI.
const autoPtr< searchableSurface > & surfPtr() const
Return a reference to the projection surface.
virtual label neighbPatchID() const
Neighbour patch ID.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
static const word null
An empty word.
Definition: word.H:77
virtual bool owner() const
Does this side own the patch?
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)
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
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: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
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.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
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:321
defineTypeNameAndDebug(combustionModel, 0)
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
const bMesh & mesh() const
Definition: boundaryMesh.H:202
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 > &)
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.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
void setSize(const label)
Reset size of List.
Definition: List.C:295
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p.
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.
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
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 vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
fileName path() const
Return path.
Definition: Time.H:269
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const word & neighbPatchName() const
Neighbour patch name.
word nbrPatchName_
Name of other half.
scalar rotationAngle_
Rotation angle.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
IOerror FatalIOError