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