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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "cyclicAMIPolyPatch.H"
30 #include "SubField.H"
31 #include "Time.H"
32 #include "unitConversion.H"
33 #include "OFstream.H"
34 #include "meshTools.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(cyclicAMIPolyPatch, 0);
42 
43  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, word);
44  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, dictionary);
45 }
46 
47 const Foam::scalar Foam::cyclicAMIPolyPatch::tolerance_ = 1e-10;
48 
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 
51 Foam::vector Foam::cyclicAMIPolyPatch::findFaceNormalMaxRadius
52 (
53  const pointField& faceCentres
54 ) const
55 {
56  // Determine a face furthest away from the axis
57 
59 
60  const scalarField magRadSqr(magSqr(n));
61 
62  label facei = findMax(magRadSqr);
63 
65  << "Patch: " << name() << nl
66  << " rotFace = " << facei << nl
67  << " point = " << faceCentres[facei] << nl
68  << " distance = " << Foam::sqrt(magRadSqr[facei])
69  << endl;
70 
71  return n[facei];
72 }
73 
74 
76 (
77  const primitivePatch& half0,
78  const pointField& half0Ctrs,
79  const vectorField& half0Areas,
80  const pointField& half1Ctrs,
81  const vectorField& half1Areas
82 )
83 {
84  if (transform() != neighbPatch().transform())
85  {
87  << "Patch " << name()
88  << " has transform type " << transformTypeNames[transform()]
89  << ", neighbour patch " << neighbPatchName()
90  << " has transform type "
91  << neighbPatch().transformTypeNames[neighbPatch().transform()]
92  << exit(FatalError);
93  }
94 
95 
96  // Calculate transformation tensors
97 
98  switch (transform())
99  {
100  case ROTATIONAL:
101  {
102  tensor revT = Zero;
103 
104  if (rotationAngleDefined_)
105  {
106  const tensor T(rotationAxis_*rotationAxis_);
107 
108  const tensor S
109  (
110  0, -rotationAxis_.z(), rotationAxis_.y(),
111  rotationAxis_.z(), 0, -rotationAxis_.x(),
112  -rotationAxis_.y(), rotationAxis_.x(), 0
113  );
114 
115  const tensor revTPos
116  (
117  T
118  + cos(rotationAngle_)*(tensor::I - T)
119  + sin(rotationAngle_)*S
120  );
121 
122  const tensor revTNeg
123  (
124  T
125  + cos(-rotationAngle_)*(tensor::I - T)
126  + sin(-rotationAngle_)*S
127  );
128 
129  // Check - assume correct angle when difference in face areas
130  // is the smallest
131  const vector transformedAreaPos = gSum(half1Areas & revTPos);
132  const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
133  const vector area0 = gSum(half0Areas);
134  const scalar magArea0 = mag(area0) + ROOTVSMALL;
135 
136  // Areas have opposite sign, so sum should be zero when correct
137  // rotation applied
138  const scalar errorPos = mag(transformedAreaPos + area0);
139  const scalar errorNeg = mag(transformedAreaNeg + area0);
140 
141  const scalar normErrorPos = errorPos/magArea0;
142  const scalar normErrorNeg = errorNeg/magArea0;
143 
144  if (errorPos > errorNeg && normErrorNeg < matchTolerance())
145  {
146  revT = revTNeg;
147  rotationAngle_ *= -1;
148  }
149  else
150  {
151  revT = revTPos;
152  }
153 
154  const scalar areaError = min(normErrorPos, normErrorNeg);
155 
156  if (areaError > matchTolerance())
157  {
159  << "Patch areas are not consistent within "
160  << 100*matchTolerance()
161  << " % indicating a possible error in the specified "
162  << "angle of rotation" << nl
163  << " owner patch : " << name() << nl
164  << " neighbour patch : " << neighbPatch().name()
165  << nl
166  << " angle : "
167  << radToDeg(rotationAngle_) << " deg" << nl
168  << " area error : " << 100*areaError << " %"
169  << " match tolerance : " << matchTolerance()
170  << endl;
171  }
172 
173  if (debug)
174  {
175  scalar theta = radToDeg(rotationAngle_);
176 
177  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
178  << name()
179  << " Specified rotation:"
180  << " swept angle: " << theta << " [deg]"
181  << " reverse transform: " << revT
182  << endl;
183  }
184  }
185  else
186  {
187  point n0 = Zero;
188  point n1 = Zero;
189  if (half0Ctrs.size())
190  {
191  n0 = findFaceNormalMaxRadius(half0Ctrs);
192  }
193  if (half1Ctrs.size())
194  {
195  n1 = -findFaceNormalMaxRadius(half1Ctrs);
196  }
197 
198  reduce(n0, maxMagSqrOp<point>());
199  reduce(n1, maxMagSqrOp<point>());
200 
201  n0.normalise();
202  n1.normalise();
203 
204  // Extended tensor from two local coordinate systems calculated
205  // using normal and rotation axis
206  const tensor E0
207  (
208  rotationAxis_,
209  (n0 ^ rotationAxis_),
210  n0
211  );
212  const tensor E1
213  (
214  rotationAxis_,
215  (-n1 ^ rotationAxis_),
216  -n1
217  );
218  revT = E1.T() & E0;
219 
220  if (debug)
221  {
222  scalar theta = radToDeg(acos(-(n0 & n1)));
223 
224  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
225  << name()
226  << " Specified rotation:"
227  << " n0:" << n0 << " n1:" << n1
228  << " swept angle: " << theta << " [deg]"
229  << " reverse transform: " << revT
230  << endl;
231  }
232  }
233 
234  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
235  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
236  const_cast<vectorField&>(separation()).setSize(0);
237  const_cast<boolList&>(collocated()) = boolList(1, false);
238 
239  break;
240  }
241  case TRANSLATIONAL:
242  {
243  if (debug)
244  {
245  Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
246  << " Specified translation : " << separationVector_
247  << endl;
248  }
249 
250  const_cast<tensorField&>(forwardT()).clear();
251  const_cast<tensorField&>(reverseT()).clear();
252  const_cast<vectorField&>(separation()) = vectorField
253  (
254  1,
255  separationVector_
256  );
257  const_cast<boolList&>(collocated()) = boolList(1, false);
258 
259  break;
260  }
261  default:
262  {
263  if (debug)
264  {
265  Pout<< "patch:" << name()
266  << " Assuming cyclic AMI pairs are colocated" << endl;
267  }
268 
269  const_cast<tensorField&>(forwardT()).clear();
270  const_cast<tensorField&>(reverseT()).clear();
271  const_cast<vectorField&>(separation()).setSize(0);
272  const_cast<boolList&>(collocated()) = boolList(1, true);
273 
274  break;
275  }
276  }
277 
278  if (debug)
279  {
280  Pout<< "patch: " << name() << nl
281  << " forwardT = " << forwardT() << nl
282  << " reverseT = " << reverseT() << nl
283  << " separation = " << separation() << nl
284  << " collocated = " << collocated() << nl << endl;
285  }
286 }
287 
288 
289 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
290 
293 {
294  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
295 
296  if (!surfPtr_ && owner() && surfType != "none")
297  {
298  word surfName(surfDict_.getOrDefault("name", name()));
299 
300  const polyMesh& mesh = boundaryMesh().mesh();
301 
302  surfPtr_ =
304  (
305  surfType,
306  IOobject
307  (
308  surfName,
309  mesh.time().constant(),
310  "triSurface",
311  mesh,
314  ),
315  surfDict_
316  );
317  }
318 
319  return surfPtr_;
320 }
321 
322 
324 {
325  resetAMI(boundaryMesh().mesh().points());
326 }
327 
328 
330 {
332 
333  if (!owner())
334  {
335  return;
336  }
337 
338  const cyclicAMIPolyPatch& nbr = neighbPatch();
339  pointField srcPoints(localPoints());
340  pointField nbrPoints(nbr.localPoints());
341 
342  if (debug)
343  {
344  const Time& t = boundaryMesh().mesh().time();
345  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
346  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
347  }
348 
349  label patchSize0 = size();
350  label nbrPatchSize0 = nbr.size();
351 
352  if (createAMIFaces_)
353  {
354  // AMI is created based on the original patch faces (non-extended patch)
355  if (srcFaceIDs_.size())
356  {
357  patchSize0 = srcFaceIDs_.size();
358  }
359  if (tgtFaceIDs_.size())
360  {
361  nbrPatchSize0 = tgtFaceIDs_.size();
362  }
363  }
364 
365  // Transform neighbour patch to local system
366  transformPosition(nbrPoints);
367  primitivePatch nbrPatch0
368  (
369  SubList<face>(nbr.localFaces(), nbrPatchSize0),
370  nbrPoints
371  );
372  primitivePatch patch0
373  (
374  SubList<face>(localFaces(), patchSize0),
375  srcPoints
376  );
377 
378 
379  if (debug)
380  {
381  const Time& t = boundaryMesh().mesh().time();
382  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
383  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
384 
385  OFstream osO(t.path()/name() + "_ownerPatch.obj");
386  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
387  }
388 
389  // Construct/apply AMI interpolation to determine addressing and weights
390  AMIPtr_->upToDate() = false;
391  AMIPtr_->calculate(patch0, nbrPatch0, surfPtr());
392 
393  if (debug)
394  {
395  AMIPtr_->checkSymmetricWeights(true);
396  }
397 }
398 
399 
401 {
403 
404  const cyclicAMIPolyPatch& half0 = *this;
405  vectorField half0Areas(half0.size());
406  forAll(half0, facei)
407  {
408  half0Areas[facei] = half0[facei].areaNormal(half0.points());
409  }
410 
411  const cyclicAMIPolyPatch& half1 = neighbPatch();
412  vectorField half1Areas(half1.size());
413  forAll(half1, facei)
414  {
415  half1Areas[facei] = half1[facei].areaNormal(half1.points());
416  }
417 
418  calcTransforms
419  (
420  half0,
421  half0.faceCentres(),
422  half0Areas,
423  half1.faceCentres(),
424  half1Areas
425  );
426 
427  DebugPout
428  << "calcTransforms() : patch: " << name() << nl
429  << " forwardT = " << forwardT() << nl
430  << " reverseT = " << reverseT() << nl
431  << " separation = " << separation() << nl
432  << " collocated = " << collocated() << nl << endl;
433 }
434 
435 
437 {
439 
440  // Flag AMI as needing update
441  AMIPtr_->upToDate() = false;
442 
444 
445  // Early calculation of transforms so e.g. cyclicACMI can use them.
446  // Note: also triggers primitiveMesh face centre. Note that cell
447  // centres should -not- be calculated
448  // since e.g. cyclicACMI overrides face areas
449  calcTransforms();
450 }
451 
452 
454 {
456 }
457 
458 
460 (
461  PstreamBuffers& pBufs,
462  const pointField& p
463 )
464 {
466 
467  // See below. Clear out any local geometry
469 
470  // Note: processorPolyPatch::initMovePoints calls
471  // processorPolyPatch::initGeometry which will trigger calculation of
472  // patch faceCentres() and cell volumes...
473 
474  if (createAMIFaces_)
475  {
476  // Note: AMI should have been updated in setTopology
477 
478  // faceAreas() and faceCentres() have been reset and will be
479  // recalculated on-demand using the mesh points and no longer
480  // correspond to the scaled areas!
481  restoreScaledGeometry();
482 
483  // deltas need to be recalculated to use new face centres!
484  }
485  else
486  {
487  AMIPtr_->upToDate() = false;
488  }
489 
490  // Early calculation of transforms. See above.
491  calcTransforms();
492 }
493 
494 
496 (
497  PstreamBuffers& pBufs,
498  const pointField& p
499 )
500 {
502 
503 // Note: not calling movePoints since this will undo our manipulations!
504 // polyPatch::movePoints(pBufs, p);
505 /*
506  polyPatch::movePoints
507  -> primitivePatch::movePoints
508  -> primitivePatch::clearGeom:
509  deleteDemandDrivenData(localPointsPtr_);
510  deleteDemandDrivenData(faceCentresPtr_);
511  deleteDemandDrivenData(faceAreasPtr_);
512  deleteDemandDrivenData(magFaceAreasPtr_);
513  deleteDemandDrivenData(faceNormalsPtr_);
514  deleteDemandDrivenData(pointNormalsPtr_);
515 */
516 }
517 
518 
520 {
522 
524 
525  if (createAMIFaces_ && boundaryMesh().mesh().topoChanging() && owner())
526  {
527  setAMIFaces();
528  }
529 }
530 
531 
533 {
535 
536  // Note: this clears out cellCentres(), faceCentres() and faceAreas()
537  polyPatch::updateMesh(pBufs);
538 }
539 
540 
542 {
544 
545  if (!updatingAMI_)
546  {
547  AMIPtr_->upToDate() = false;
548  }
549 
551 }
552 
553 
554 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
555 
557 (
558  const word& name,
559  const label size,
560  const label start,
561  const label index,
562  const polyBoundaryMesh& bm,
563  const word& patchType,
564  const transformType transform,
565  const word& defaultAMIMethod
566 )
567 :
568  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
569  nbrPatchName_(word::null),
570  nbrPatchID_(-1),
571  fraction_(Zero),
572  rotationAxis_(Zero),
573  rotationCentre_(Zero),
574  rotationAngleDefined_(false),
575  rotationAngle_(0.0),
576  separationVector_(Zero),
577  AMIPtr_(AMIInterpolation::New(defaultAMIMethod)),
578  surfDict_(fileName("surface")),
579  surfPtr_(nullptr),
580  createAMIFaces_(false),
581  moveFaceCentres_(false),
582  updatingAMI_(true),
583  srcFaceIDs_(),
584  tgtFaceIDs_(),
585  faceAreas0_(),
586  faceCentres0_()
587 {
588  // Neighbour patch might not be valid yet so no transformation
589  // calculation possible
590 }
591 
592 
594 (
595  const word& name,
596  const dictionary& dict,
597  const label index,
598  const polyBoundaryMesh& bm,
599  const word& patchType,
600  const word& defaultAMIMethod
601 )
602 :
603  coupledPolyPatch(name, dict, index, bm, patchType),
604  nbrPatchName_(dict.getOrDefault<word>("neighbourPatch", word::null)),
605  coupleGroup_(dict),
606  nbrPatchID_(-1),
607  fraction_(dict.getOrDefault<scalar>("fraction", Zero)),
608  rotationAxis_(Zero),
609  rotationCentre_(Zero),
610  rotationAngleDefined_(false),
611  rotationAngle_(0.0),
612  separationVector_(Zero),
613  AMIPtr_
614  (
616  (
617  dict.getOrDefault<word>("AMIMethod", defaultAMIMethod),
618  dict,
619  dict.getOrDefault("flipNormals", false)
620  )
621  ),
622  surfDict_(dict.subOrEmptyDict("surface")),
623  surfPtr_(nullptr),
624  createAMIFaces_(dict.getOrDefault("createAMIFaces", false)),
625  moveFaceCentres_(false),
626  updatingAMI_(true),
627  srcFaceIDs_(),
628  tgtFaceIDs_(),
629  faceAreas0_(),
630  faceCentres0_()
631 {
632  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
633  {
635  << "No \"neighbourPatch\" or \"coupleGroup\" provided."
636  << exit(FatalIOError);
637  }
638 
639  if (nbrPatchName_ == name)
640  {
642  << "Neighbour patch name " << nbrPatchName_
643  << " cannot be the same as this patch " << name
644  << exit(FatalIOError);
645  }
646 
647  switch (transform())
648  {
649  case ROTATIONAL:
650  {
651  dict.readEntry("rotationAxis", rotationAxis_);
652  dict.readEntry("rotationCentre", rotationCentre_);
653  if (dict.readIfPresent("rotationAngle", rotationAngle_))
654  {
655  rotationAngleDefined_ = true;
656  rotationAngle_ = degToRad(rotationAngle_);
657 
658  if (debug)
659  {
660  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
661  << endl;
662  }
663  }
664 
665  scalar magRot = mag(rotationAxis_);
666  if (magRot < SMALL)
667  {
669  << "Illegal rotationAxis " << rotationAxis_ << endl
670  << "Please supply a non-zero vector."
671  << exit(FatalIOError);
672  }
673  rotationAxis_ /= magRot;
674 
675  break;
676  }
677  case TRANSLATIONAL:
678  {
679  dict.readEntry("separationVector", separationVector_);
680  break;
681  }
682  default:
683  {
684  // No additional info required
685  }
686  }
687 
688  // Neighbour patch might not be valid yet so no transformation
689  // calculation possible
690 
691  // If topology change, recover the sizes of the original patches and
692  // read additional controls
693  if (createAMIFaces_)
694  {
695  srcFaceIDs_.setSize(dict.get<label>("srcSize"));
696  tgtFaceIDs_.setSize(dict.get<label>("tgtSize"));
697  moveFaceCentres_ = dict.getOrDefault("moveFaceCentres", true);
698  }
699 }
700 
701 
703 (
704  const cyclicAMIPolyPatch& pp,
705  const polyBoundaryMesh& bm
706 )
707 :
708  coupledPolyPatch(pp, bm),
709  nbrPatchName_(pp.nbrPatchName_),
710  coupleGroup_(pp.coupleGroup_),
711  nbrPatchID_(-1),
712  fraction_(pp.fraction_),
713  rotationAxis_(pp.rotationAxis_),
714  rotationCentre_(pp.rotationCentre_),
715  rotationAngleDefined_(pp.rotationAngleDefined_),
716  rotationAngle_(pp.rotationAngle_),
717  separationVector_(pp.separationVector_),
718  AMIPtr_(pp.AMIPtr_->clone()),
719  surfDict_(pp.surfDict_),
720  surfPtr_(nullptr),
721  createAMIFaces_(pp.createAMIFaces_),
722  moveFaceCentres_(pp.moveFaceCentres_),
723  updatingAMI_(true),
724  srcFaceIDs_(),
725  tgtFaceIDs_(),
726  faceAreas0_(),
727  faceCentres0_()
728 {
729  // Neighbour patch might not be valid yet so no transformation
730  // calculation possible
731 }
732 
733 
735 (
736  const cyclicAMIPolyPatch& pp,
737  const polyBoundaryMesh& bm,
738  const label index,
739  const label newSize,
740  const label newStart,
741  const word& nbrPatchName
742 )
743 :
744  coupledPolyPatch(pp, bm, index, newSize, newStart),
745  nbrPatchName_(nbrPatchName),
746  coupleGroup_(pp.coupleGroup_),
747  nbrPatchID_(-1),
748  fraction_(pp.fraction_),
749  rotationAxis_(pp.rotationAxis_),
750  rotationCentre_(pp.rotationCentre_),
751  rotationAngleDefined_(pp.rotationAngleDefined_),
752  rotationAngle_(pp.rotationAngle_),
753  separationVector_(pp.separationVector_),
754  AMIPtr_(pp.AMIPtr_->clone()),
755  surfDict_(pp.surfDict_),
756  surfPtr_(nullptr),
757  createAMIFaces_(pp.createAMIFaces_),
758  moveFaceCentres_(pp.moveFaceCentres_),
759  updatingAMI_(true),
760  srcFaceIDs_(),
761  tgtFaceIDs_(),
762  faceAreas0_(),
763  faceCentres0_()
764 {
765  if (nbrPatchName_ == name())
766  {
768  << "Neighbour patch name " << nbrPatchName_
769  << " cannot be the same as this patch " << name()
770  << exit(FatalError);
771  }
772 
773  // Neighbour patch might not be valid yet so no transformation
774  // calculation possible
775 }
776 
777 
779 (
780  const cyclicAMIPolyPatch& pp,
781  const polyBoundaryMesh& bm,
782  const label index,
783  const labelUList& mapAddressing,
784  const label newStart
785 )
786 :
787  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
788  nbrPatchName_(pp.nbrPatchName_),
789  coupleGroup_(pp.coupleGroup_),
790  nbrPatchID_(-1),
791  fraction_(pp.fraction_),
792  rotationAxis_(pp.rotationAxis_),
793  rotationCentre_(pp.rotationCentre_),
794  rotationAngleDefined_(pp.rotationAngleDefined_),
795  rotationAngle_(pp.rotationAngle_),
796  separationVector_(pp.separationVector_),
797  AMIPtr_(pp.AMIPtr_->clone()),
798  surfDict_(pp.surfDict_),
799  surfPtr_(nullptr),
800  createAMIFaces_(pp.createAMIFaces_),
801  moveFaceCentres_(pp.moveFaceCentres_),
802  updatingAMI_(true),
803  srcFaceIDs_(),
804  tgtFaceIDs_(),
805  faceAreas0_(),
806  faceCentres0_()
807 {}
808 
809 
810 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
811 
813 {
814  if (nbrPatchID_ == -1)
815  {
816  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
817 
818  if (nbrPatchID_ == -1)
819  {
821  << "Illegal neighbourPatch name " << neighbPatchName()
822  << nl << "Valid patch names are "
823  << this->boundaryMesh().names()
824  << exit(FatalError);
825  }
826 
827  // Check that it is a cyclic AMI patch
828  const cyclicAMIPolyPatch& nbrPatch =
829  refCast<const cyclicAMIPolyPatch>
830  (
831  this->boundaryMesh()[nbrPatchID_]
832  );
833 
834  if (nbrPatch.neighbPatchName() != name())
835  {
837  << "Patch " << name()
838  << " specifies neighbour patch " << neighbPatchName()
839  << nl << " but that in return specifies "
840  << nbrPatch.neighbPatchName() << endl;
841  }
842  }
843 
844  return nbrPatchID_;
845 }
846 
847 
849 {
850  return index() < neighbPatchID();
851 }
852 
853 
855 {
856  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
857  return refCast<const cyclicAMIPolyPatch>(pp);
858 }
859 
860 
862 {
863  if (!owner())
864  {
866  << "AMI interpolator only available to owner patch"
867  << abort(FatalError);
868  }
869 
870  if (!AMIPtr_->upToDate())
871  {
872  resetAMI();
873  }
874 
875  return *AMIPtr_;
876 }
877 
878 
880 {
881  if (owner())
882  {
883  return AMI().applyLowWeightCorrection();
884  }
885  else
886  {
887  return neighbPatch().AMI().applyLowWeightCorrection();
888  }
889 }
890 
891 
893 {
894  if (!parallel())
895  {
896  if (transform() == ROTATIONAL)
897  {
898  l = Foam::transform(forwardT(), l - rotationCentre_)
899  + rotationCentre_;
900  }
901  else
902  {
903  l = Foam::transform(forwardT(), l);
904  }
905  }
906  else if (separated())
907  {
908  // transformPosition gets called on the receiving side,
909  // separation gets calculated on the sending side so subtract
910 
911  const vectorField& s = separation();
912  if (s.size() == 1)
913  {
914  forAll(l, i)
915  {
916  l[i] -= s[0];
917  }
918  }
919  else
920  {
921  l -= s;
922  }
923  }
924 }
925 
926 
928 (
929  point& l,
930  const label facei
931 ) const
932 {
933  if (!parallel())
934  {
935  const tensor& T =
936  (
937  forwardT().size() == 1
938  ? forwardT()[0]
939  : forwardT()[facei]
940  );
941 
942  if (transform() == ROTATIONAL)
943  {
944  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
945  }
946  else
947  {
948  l = Foam::transform(T, l);
949  }
950  }
951  else if (separated())
952  {
953  const vector& s =
954  (
955  separation().size() == 1
956  ? separation()[0]
957  : separation()[facei]
958  );
959 
960  l -= s;
961  }
962 }
963 
964 
966 (
967  point& l,
968  const label facei
969 ) const
970 {
971  if (!parallel())
972  {
973  const tensor& T =
974  (
975  reverseT().size() == 1
976  ? reverseT()[0]
977  : reverseT()[facei]
978  );
979 
980  if (transform() == ROTATIONAL)
981  {
982  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
983  }
984  else
985  {
986  l = Foam::transform(T, l);
987  }
988  }
989  else if (separated())
990  {
991  const vector& s =
992  (
993  separation().size() == 1
994  ? separation()[0]
995  : separation()[facei]
996  );
997 
998  l += s;
999  }
1000 }
1001 
1002 
1005  vector& d,
1006  const label facei
1007 ) const
1008 {
1009  if (!parallel())
1010  {
1011  const tensor& T =
1012  (
1013  reverseT().size() == 1
1014  ? reverseT()[0]
1015  : reverseT()[facei]
1016  );
1017 
1018  d = Foam::transform(T, d);
1019  }
1020 }
1021 
1022 
1025  const primitivePatch& referPatch,
1026  const pointField& thisCtrs,
1027  const vectorField& thisAreas,
1028  const pointField& thisCc,
1029  const pointField& nbrCtrs,
1030  const vectorField& nbrAreas,
1031  const pointField& nbrCc
1032 )
1033 {}
1034 
1035 
1038  PstreamBuffers& pBufs,
1039  const primitivePatch& pp
1040 ) const
1041 {}
1042 
1043 
1046  PstreamBuffers& pBufs,
1047  const primitivePatch& pp,
1048  labelList& faceMap,
1049  labelList& rotation
1050 ) const
1051 {
1052  faceMap.setSize(pp.size());
1053  faceMap = -1;
1054 
1055  rotation.setSize(pp.size());
1056  rotation = 0;
1057 
1058  return false;
1059 }
1060 
1061 
1064  const label facei,
1065  const vector& n,
1066  point& p
1067 ) const
1068 {
1069  point prt(p);
1070  reverseTransformPosition(prt, facei);
1071 
1072  vector nrt(n);
1073  reverseTransformDirection(nrt, facei);
1074 
1075  label nbrFacei = -1;
1076 
1077  if (owner())
1078  {
1079  nbrFacei = AMI().tgtPointFace
1080  (
1081  *this,
1082  neighbPatch(),
1083  nrt,
1084  facei,
1085  prt
1086  );
1087  }
1088  else
1089  {
1090  nbrFacei = neighbPatch().AMI().srcPointFace
1091  (
1092  neighbPatch(),
1093  *this,
1094  nrt,
1095  facei,
1096  prt
1097  );
1098  }
1099 
1100  if (nbrFacei >= 0)
1101  {
1102  p = prt;
1103  }
1104 
1105  return nbrFacei;
1106 }
1107 
1108 
1110 {
1112  if (!nbrPatchName_.empty())
1113  {
1114  os.writeEntry("neighbourPatch", nbrPatchName_);
1115  }
1116  coupleGroup_.write(os);
1117 
1118  switch (transform())
1119  {
1120  case ROTATIONAL:
1121  {
1122  os.writeEntry("rotationAxis", rotationAxis_);
1123  os.writeEntry("rotationCentre", rotationCentre_);
1124 
1125  if (rotationAngleDefined_)
1126  {
1127  os.writeEntry("rotationAngle", radToDeg(rotationAngle_));
1128  }
1129 
1130  break;
1131  }
1132  case TRANSLATIONAL:
1133  {
1134  os.writeEntry("separationVector", separationVector_);
1135  break;
1136  }
1137  case NOORDERING:
1138  {
1139  break;
1140  }
1141  default:
1142  {
1143  // No additional info to write
1144  }
1145  }
1146 
1147  AMIPtr_->write(os);
1148 
1149  if (!surfDict_.empty())
1150  {
1151  surfDict_.writeEntry(surfDict_.dictName(), os);
1152  }
1153 
1154  if (createAMIFaces_)
1155  {
1156  os.writeEntry("createAMIFaces", createAMIFaces_);
1157  os.writeEntry("srcSize", srcFaceIDs_.size());
1158  os.writeEntry("tgtSize", tgtFaceIDs_.size());
1159  os.writeEntry("moveFaceCentres", moveFaceCentres_);
1160  }
1161 
1162  os.writeEntryIfDifferent<scalar>("fraction", Zero, fraction_);
1163 }
1164 
1165 
1166 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::cyclicAMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicAMIPolyPatch.C:1045
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:206
Foam::Tensor< scalar >
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
setSize
points setSize(newPointi)
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
meshTools.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::cyclicAMIPolyPatch::resetAMI
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
Definition: cyclicAMIPolyPatch.C:323
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:848
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
SubField.H
Foam::Tensor::I
static const Tensor I
Definition: Tensor.H:85
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::cyclicAMIPolyPatch::nbrPatchName_
word nbrPatchName_
Name of other half.
Definition: cyclicAMIPolyPatch.H:91
Foam::cyclicAMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Definition: cyclicAMIPolyPatch.C:1037
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::cyclicAMIPolyPatch::rotationCentre_
point rotationCentre_
Point on axis of rotation for rotational cyclics.
Definition: cyclicAMIPolyPatch.H:111
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::cyclicAMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicAMIPolyPatch.C:436
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:87
Foam::polyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:74
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:966
Foam::cyclicAMIPolyPatch::surfDict_
const dictionary surfDict_
Dictionary used during projection surface construction.
Definition: cyclicAMIPolyPatch.H:130
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::cyclicAMIPolyPatch::surfPtr
const autoPtr< searchableSurface > & surfPtr() const
Create and return pointer to the projection surface.
Definition: cyclicAMIPolyPatch.C:292
unitConversion.H
Unit conversion functions.
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::cyclicAMIPolyPatch::coupleGroup_
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
Definition: cyclicAMIPolyPatch.H:94
Foam::cyclicAMIPolyPatch::tolerance_
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
Definition: cyclicAMIPolyPatch.H:331
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:53
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::AMIInterpolation::New
static autoPtr< AMIInterpolation > New(const word &modelName, const dictionary &dict, const bool reverseTarget=false)
Selector for dictionary.
Definition: AMIInterpolationNew.C:33
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::polyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:67
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::cyclicAMIPolyPatch::initMovePoints
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
Definition: cyclicAMIPolyPatch.C:460
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::cyclicAMIPolyPatch::updateMesh
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: cyclicAMIPolyPatch.C:532
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:892
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::polyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:117
Foam::cyclicAMIPolyPatch::separationVector_
vector separationVector_
Translation vector.
Definition: cyclicAMIPolyPatch.H:123
Foam::cyclicAMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicAMIPolyPatch.C:453
Foam::Field< vector >
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1109
Foam::cyclicAMIPolyPatch::AMIPtr_
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
Definition: cyclicAMIPolyPatch.H:127
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cyclicAMIPolyPatch::pointFace
label pointFace(const label facei, const vector &n, point &p) const
Definition: cyclicAMIPolyPatch.C:1063
Foam::cyclicAMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicAMIPolyPatch.C:496
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::radToDeg
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
Definition: unitConversion.H:54
Foam::polyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:102
Foam::coupledPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: coupledPolyPatch.C:567
Foam::cyclicAMIPolyPatch::calcTransforms
virtual void calcTransforms()
Recalculate the transformation tensors.
Definition: cyclicAMIPolyPatch.C:400
Foam::cyclicAMIPolyPatch::rotationAxis_
vector rotationAxis_
Axis of rotation for rotational cyclics.
Definition: cyclicAMIPolyPatch.H:108
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicAMIPolyPatch::rotationAngle_
scalar rotationAngle_
Rotation angle.
Definition: cyclicAMIPolyPatch.H:117
Foam::cyclicAMIPolyPatch::neighbPatchID
virtual label neighbPatchID() const
Neighbour patch ID.
Definition: cyclicAMIPolyPatch.C:812
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::cyclicAMIPolyPatch::neighbPatch
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Definition: cyclicAMIPolyPatch.C:854
Foam::word::valid
static bool valid(char c)
Is this character valid for a word?
Definition: wordI.H:130
Foam::PrimitivePatch::points
const Field< point_type > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:305
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:297
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, const word &defaultAMIMethod=faceAreaWeightAMI::typeName)
Construct from (base coupled patch) components.
Definition: cyclicAMIPolyPatch.C:557
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cyclicAMIPolyPatch::rotationAngleDefined_
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
Definition: cyclicAMIPolyPatch.H:114
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
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:43
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:370
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:608
cyclicAMIPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:339
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::cyclicAMIPolyPatch::moveFaceCentres_
bool moveFaceCentres_
Move face centres (default = no)
Definition: cyclicAMIPolyPatch.H:143
Foam::cyclicAMIPolyPatch::fraction_
const scalar fraction_
Particle displacement fraction accross AMI.
Definition: cyclicAMIPolyPatch.H:100
Time.H
Foam::autoPtr< Foam::searchableSurface >
Foam::cyclicAMIPolyPatch::AMI
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
Definition: cyclicAMIPolyPatch.C:861
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:358
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:313
Foam::cyclicAMIPolyPatch::applyLowWeightCorrection
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
Definition: cyclicAMIPolyPatch.C:879
Foam::cyclicAMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicAMIPolyPatch.C:519
Foam::cyclicAMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicAMIPolyPatch.C:541
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::findMax
label findMax(const ListType &input, label start=0)
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
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:1004
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
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:103
Foam::cyclicAMIPolyPatch::neighbPatchName
const word & neighbPatchName() const
Neighbour patch name.
Definition: cyclicAMIPolyPatchI.H:49
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:62
Foam::primitivePatch
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Definition: primitivePatch.H:51
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
Foam::PrimitivePatch::movePoints
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
Definition: PrimitivePatch.C:172
Foam::cyclicAMIPolyPatch::createAMIFaces_
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
Definition: cyclicAMIPolyPatch.H:140
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:59
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::patchIdentifier::name
const word & name() const
The patch name.
Definition: patchIdentifier.H:134
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:61
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:67