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-2021 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 #include "cyclicPolyPatch.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(cyclicAMIPolyPatch, 0);
43 
44  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, word);
45  addToRunTimeSelectionTable(polyPatch, cyclicAMIPolyPatch, dictionary);
46 }
47 
48 const Foam::scalar Foam::cyclicAMIPolyPatch::tolerance_ = 1e-10;
49 
50 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
51 
52 Foam::vector Foam::cyclicAMIPolyPatch::findFaceNormalMaxRadius
53 (
54  const pointField& faceCentres
55 ) const
56 {
57  // Determine a face furthest away from the axis
58 
60 
61  const scalarField magRadSqr(magSqr(n));
62 
63  label facei = findMax(magRadSqr);
64 
66  << "Patch: " << name() << nl
67  << " rotFace = " << facei << nl
68  << " point = " << faceCentres[facei] << nl
69  << " distance = " << Foam::sqrt(magRadSqr[facei])
70  << endl;
71 
72  return n[facei];
73 }
74 
75 
77 (
78  const primitivePatch& half0,
79  const pointField& half0Ctrs,
80  const vectorField& half0Areas,
81  const pointField& half1Ctrs,
82  const vectorField& half1Areas
83 )
84 {
85  if (transform() != neighbPatch().transform())
86  {
88  << "Patch " << name()
89  << " has transform type " << transformTypeNames[transform()]
90  << ", neighbour patch " << neighbPatchName()
91  << " has transform type "
92  << neighbPatch().transformTypeNames[neighbPatch().transform()]
93  << exit(FatalError);
94  }
95 
96 
97  // Calculate transformation tensors
98 
99  switch (transform())
100  {
101  case ROTATIONAL:
102  {
103  tensor revT = Zero;
104 
105  if (rotationAngleDefined_)
106  {
107  const tensor T(rotationAxis_*rotationAxis_);
108 
109  const tensor S
110  (
111  0, -rotationAxis_.z(), rotationAxis_.y(),
112  rotationAxis_.z(), 0, -rotationAxis_.x(),
113  -rotationAxis_.y(), rotationAxis_.x(), 0
114  );
115 
116  const tensor revTPos
117  (
118  T
119  + cos(rotationAngle_)*(tensor::I - T)
120  + sin(rotationAngle_)*S
121  );
122 
123  const tensor revTNeg
124  (
125  T
126  + cos(-rotationAngle_)*(tensor::I - T)
127  + sin(-rotationAngle_)*S
128  );
129 
130  // Check - assume correct angle when difference in face areas
131  // is the smallest
132  const vector transformedAreaPos = gSum(half1Areas & revTPos);
133  const vector transformedAreaNeg = gSum(half1Areas & revTNeg);
134  const vector area0 = gSum(half0Areas);
135  const scalar magArea0 = mag(area0) + ROOTVSMALL;
136 
137  // Areas have opposite sign, so sum should be zero when correct
138  // rotation applied
139  const scalar errorPos = mag(transformedAreaPos + area0);
140  const scalar errorNeg = mag(transformedAreaNeg + area0);
141 
142  const scalar normErrorPos = errorPos/magArea0;
143  const scalar normErrorNeg = errorNeg/magArea0;
144 
145  if (errorPos > errorNeg && normErrorNeg < matchTolerance())
146  {
147  revT = revTNeg;
148  rotationAngle_ *= -1;
149  }
150  else
151  {
152  revT = revTPos;
153  }
154 
155  const scalar areaError = min(normErrorPos, normErrorNeg);
156 
157  if (areaError > matchTolerance())
158  {
160  << "Patch areas are not consistent within "
161  << 100*matchTolerance()
162  << " % indicating a possible error in the specified "
163  << "angle of rotation" << nl
164  << " owner patch : " << name() << nl
165  << " neighbour patch : " << neighbPatch().name()
166  << nl
167  << " angle : "
168  << radToDeg(rotationAngle_) << " deg" << nl
169  << " area error : " << 100*areaError << " %"
170  << " match tolerance : " << matchTolerance()
171  << endl;
172  }
173 
174  if (debug)
175  {
176  scalar theta = radToDeg(rotationAngle_);
177 
178  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
179  << name()
180  << " Specified rotation:"
181  << " swept angle: " << theta << " [deg]"
182  << " reverse transform: " << revT
183  << endl;
184  }
185  }
186  else
187  {
188  point n0 = Zero;
189  point n1 = Zero;
190  if (half0Ctrs.size())
191  {
192  n0 = findFaceNormalMaxRadius(half0Ctrs);
193  }
194  if (half1Ctrs.size())
195  {
196  n1 = -findFaceNormalMaxRadius(half1Ctrs);
197  }
198 
199  reduce(n0, maxMagSqrOp<point>());
200  reduce(n1, maxMagSqrOp<point>());
201 
202  n0.normalise();
203  n1.normalise();
204 
205  // Extended tensor from two local coordinate systems calculated
206  // using normal and rotation axis
207  const tensor E0
208  (
209  rotationAxis_,
210  (n0 ^ rotationAxis_),
211  n0
212  );
213  const tensor E1
214  (
215  rotationAxis_,
216  (-n1 ^ rotationAxis_),
217  -n1
218  );
219  revT = E1.T() & E0;
220 
221  if (debug)
222  {
223  scalar theta = radToDeg(acos(-(n0 & n1)));
224 
225  Pout<< "cyclicAMIPolyPatch::calcTransforms: patch:"
226  << name()
227  << " Specified rotation:"
228  << " n0:" << n0 << " n1:" << n1
229  << " swept angle: " << theta << " [deg]"
230  << " reverse transform: " << revT
231  << endl;
232  }
233  }
234 
235  const_cast<tensorField&>(forwardT()) = tensorField(1, revT.T());
236  const_cast<tensorField&>(reverseT()) = tensorField(1, revT);
237  const_cast<vectorField&>(separation()).setSize(0);
238  const_cast<boolList&>(collocated()) = boolList(1, false);
239 
240  break;
241  }
242  case TRANSLATIONAL:
243  {
244  if (debug)
245  {
246  Pout<< "cyclicAMIPolyPatch::calcTransforms : patch:" << name()
247  << " Specified translation : " << separationVector_
248  << endl;
249  }
250 
251  const_cast<tensorField&>(forwardT()).clear();
252  const_cast<tensorField&>(reverseT()).clear();
253  const_cast<vectorField&>(separation()) = vectorField
254  (
255  1,
256  separationVector_
257  );
258  const_cast<boolList&>(collocated()) = boolList(1, false);
259 
260  break;
261  }
262  default:
263  {
264  if (debug)
265  {
266  Pout<< "patch:" << name()
267  << " Assuming cyclic AMI pairs are colocated" << endl;
268  }
269 
270  const_cast<tensorField&>(forwardT()).clear();
271  const_cast<tensorField&>(reverseT()).clear();
272  const_cast<vectorField&>(separation()).setSize(0);
273  const_cast<boolList&>(collocated()) = boolList(1, true);
274 
275  break;
276  }
277  }
278 
279  if (debug)
280  {
281  Pout<< "patch: " << name() << nl
282  << " forwardT = " << forwardT() << nl
283  << " reverseT = " << reverseT() << nl
284  << " separation = " << separation() << nl
285  << " collocated = " << collocated() << nl << endl;
286  }
287 }
288 
289 
292 {
294 
295  const label periodicID = periodicPatchID();
296  if (periodicID != -1)
297  {
298  // Get the periodic patch
299  const coupledPolyPatch& perPp
300  (
301  refCast<const coupledPolyPatch>
302  (
303  boundaryMesh()[periodicID]
304  )
305  );
306  if (!perPp.parallel())
307  {
308  vector axis(Zero);
309  point axisPoint(Zero);
310  if (isA<cyclicPolyPatch>(perPp))
311  {
312  axis =
313  refCast<const cyclicPolyPatch>(perPp).rotationAxis();
314  axisPoint =
315  refCast<const cyclicPolyPatch>(perPp).rotationCentre();
316  }
317  else if (isA<cyclicAMIPolyPatch>(perPp))
318  {
319  axis =
320  refCast<const cyclicAMIPolyPatch>(perPp).rotationAxis();
321  axisPoint =
322  refCast<const cyclicAMIPolyPatch>(perPp).rotationCentre();
323  }
324  else
325  {
326  FatalErrorInFunction << "On patch " << name()
327  << " have unsupported periodicPatch " << perPp.name()
328  << exit(FatalError);
329  }
330 
331  csPtr.set(new coordSystem::cylindrical(axisPoint, axis));
332  }
333  }
334  return csPtr;
335 }
336 
337 
338 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
339 
342 {
343  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
344 
345  if (!surfPtr_ && owner() && surfType != "none")
346  {
347  word surfName(surfDict_.getOrDefault("name", name()));
348 
349  const polyMesh& mesh = boundaryMesh().mesh();
350 
351  surfPtr_ =
353  (
354  surfType,
355  IOobject
356  (
357  surfName,
358  mesh.time().constant(),
359  "triSurface",
360  mesh,
363  ),
364  surfDict_
365  );
366  }
367 
368  return surfPtr_;
369 }
370 
371 
373 {
374  resetAMI(boundaryMesh().mesh().points());
375 }
376 
377 
379 {
381 
382  if (!owner())
383  {
384  return;
385  }
386 
387  const cyclicAMIPolyPatch& nbr = neighbPatch();
388  pointField srcPoints(localPoints());
389  pointField nbrPoints(nbr.localPoints());
390 
391  if (debug)
392  {
393  const Time& t = boundaryMesh().mesh().time();
394  OFstream os(t.path()/name() + "_neighbourPatch-org.obj");
395  meshTools::writeOBJ(os, neighbPatch().localFaces(), nbrPoints);
396  }
397 
398  label patchSize0 = size();
399  label nbrPatchSize0 = nbr.size();
400 
401  if (createAMIFaces_)
402  {
403  // AMI is created based on the original patch faces (non-extended patch)
404  if (srcFaceIDs_.size())
405  {
406  patchSize0 = srcFaceIDs_.size();
407  }
408  if (tgtFaceIDs_.size())
409  {
410  nbrPatchSize0 = tgtFaceIDs_.size();
411  }
412  }
413 
414  // Transform neighbour patch to local system
415  transformPosition(nbrPoints);
416  primitivePatch nbrPatch0
417  (
418  SubList<face>(nbr.localFaces(), nbrPatchSize0),
419  nbrPoints
420  );
421  primitivePatch patch0
422  (
423  SubList<face>(localFaces(), patchSize0),
424  srcPoints
425  );
426 
427 
428  if (debug)
429  {
430  const Time& t = boundaryMesh().mesh().time();
431  OFstream osN(t.path()/name() + "_neighbourPatch-trans.obj");
432  meshTools::writeOBJ(osN, nbrPatch0.localFaces(), nbrPoints);
433 
434  OFstream osO(t.path()/name() + "_ownerPatch.obj");
435  meshTools::writeOBJ(osO, this->localFaces(), localPoints());
436  }
437 
438  // Construct/apply AMI interpolation to determine addressing and weights
439  AMIPtr_->upToDate() = false;
440  AMIPtr_->calculate(patch0, nbrPatch0, surfPtr());
441 
442  if (debug)
443  {
444  AMIPtr_->checkSymmetricWeights(true);
445  }
446 }
447 
448 
450 {
452 
453  const cyclicAMIPolyPatch& half0 = *this;
454  vectorField half0Areas(half0.size());
455  forAll(half0, facei)
456  {
457  half0Areas[facei] = half0[facei].areaNormal(half0.points());
458  }
459 
460  const cyclicAMIPolyPatch& half1 = neighbPatch();
461  vectorField half1Areas(half1.size());
462  forAll(half1, facei)
463  {
464  half1Areas[facei] = half1[facei].areaNormal(half1.points());
465  }
466 
467  calcTransforms
468  (
469  half0,
470  half0.faceCentres(),
471  half0Areas,
472  half1.faceCentres(),
473  half1Areas
474  );
475 
476  DebugPout
477  << "calcTransforms() : patch: " << name() << nl
478  << " forwardT = " << forwardT() << nl
479  << " reverseT = " << reverseT() << nl
480  << " separation = " << separation() << nl
481  << " collocated = " << collocated() << nl << endl;
482 }
483 
484 
486 {
488 
489  // Flag AMI as needing update
490  AMIPtr_->upToDate() = false;
491 
493 
494  // Early calculation of transforms so e.g. cyclicACMI can use them.
495  // Note: also triggers primitiveMesh face centre. Note that cell
496  // centres should -not- be calculated
497  // since e.g. cyclicACMI overrides face areas
498  calcTransforms();
499 }
500 
501 
503 {
505 }
506 
507 
509 (
510  PstreamBuffers& pBufs,
511  const pointField& p
512 )
513 {
515 
516  // See below. Clear out any local geometry
518 
519  // Note: processorPolyPatch::initMovePoints calls
520  // processorPolyPatch::initGeometry which will trigger calculation of
521  // patch faceCentres() and cell volumes...
522 
523  if (createAMIFaces_)
524  {
525  // Note: AMI should have been updated in setTopology
526 
527  // faceAreas() and faceCentres() have been reset and will be
528  // recalculated on-demand using the mesh points and no longer
529  // correspond to the scaled areas!
530  restoreScaledGeometry();
531 
532  // deltas need to be recalculated to use new face centres!
533  }
534  else
535  {
536  AMIPtr_->upToDate() = false;
537  }
538 
539  // Early calculation of transforms. See above.
540  calcTransforms();
541 }
542 
543 
545 (
546  PstreamBuffers& pBufs,
547  const pointField& p
548 )
549 {
551 
552 // Note: not calling movePoints since this will undo our manipulations!
553 // polyPatch::movePoints(pBufs, p);
554 /*
555  polyPatch::movePoints
556  -> primitivePatch::movePoints
557  -> primitivePatch::clearGeom:
558  deleteDemandDrivenData(localPointsPtr_);
559  deleteDemandDrivenData(faceCentresPtr_);
560  deleteDemandDrivenData(faceAreasPtr_);
561  deleteDemandDrivenData(magFaceAreasPtr_);
562  deleteDemandDrivenData(faceNormalsPtr_);
563  deleteDemandDrivenData(pointNormalsPtr_);
564 */
565 }
566 
567 
569 {
571 
573 
574  if (createAMIFaces_ && boundaryMesh().mesh().topoChanging() && owner())
575  {
576  setAMIFaces();
577  }
578 }
579 
580 
582 {
584 
585  // Note: this clears out cellCentres(), faceCentres() and faceAreas()
586  polyPatch::updateMesh(pBufs);
587 }
588 
589 
591 {
593 
594  if (!updatingAMI_)
595  {
596  AMIPtr_->upToDate() = false;
597  }
598 
600 }
601 
602 
603 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
604 
606 (
607  const word& name,
608  const label size,
609  const label start,
610  const label index,
611  const polyBoundaryMesh& bm,
612  const word& patchType,
613  const transformType transform,
614  const word& defaultAMIMethod
615 )
616 :
617  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
618  nbrPatchName_(word::null),
619  nbrPatchID_(-1),
620  fraction_(Zero),
621  rotationAxis_(Zero),
622  rotationCentre_(Zero),
623  rotationAngleDefined_(false),
624  rotationAngle_(0.0),
625  separationVector_(Zero),
626  periodicPatchName_(word::null),
627  periodicPatchID_(-1),
628  AMIPtr_(AMIInterpolation::New(defaultAMIMethod)),
629  surfDict_(fileName("surface")),
630  surfPtr_(nullptr),
631  createAMIFaces_(false),
632  moveFaceCentres_(false),
633  updatingAMI_(true),
634  srcFaceIDs_(),
635  tgtFaceIDs_(),
636  faceAreas0_(),
637  faceCentres0_()
638 {
639  // Neighbour patch might not be valid yet so no transformation
640  // calculation possible
641 }
642 
643 
645 (
646  const word& name,
647  const dictionary& dict,
648  const label index,
649  const polyBoundaryMesh& bm,
650  const word& patchType,
651  const word& defaultAMIMethod
652 )
653 :
654  coupledPolyPatch(name, dict, index, bm, patchType),
655  nbrPatchName_(dict.getOrDefault<word>("neighbourPatch", word::null)),
656  coupleGroup_(dict),
657  nbrPatchID_(-1),
658  fraction_(dict.getOrDefault<scalar>("fraction", Zero)),
659  rotationAxis_(Zero),
660  rotationCentre_(Zero),
661  rotationAngleDefined_(false),
662  rotationAngle_(0.0),
663  separationVector_(Zero),
664  periodicPatchName_(dict.getOrDefault<word>("periodicPatch", word::null)),
665  periodicPatchID_(-1),
666  AMIPtr_
667  (
669  (
670  dict.getOrDefault<word>("AMIMethod", defaultAMIMethod),
671  dict,
672  dict.getOrDefault("flipNormals", false)
673  )
674  ),
675  surfDict_(dict.subOrEmptyDict("surface")),
676  surfPtr_(nullptr),
677  createAMIFaces_(dict.getOrDefault("createAMIFaces", false)),
678  moveFaceCentres_(false),
679  updatingAMI_(true),
680  srcFaceIDs_(),
681  tgtFaceIDs_(),
682  faceAreas0_(),
683  faceCentres0_()
684 {
685  if (nbrPatchName_ == word::null && !coupleGroup_.valid())
686  {
688  << "No \"neighbourPatch\" or \"coupleGroup\" provided."
689  << exit(FatalIOError);
690  }
691 
692  if (nbrPatchName_ == name)
693  {
695  << "Neighbour patch name " << nbrPatchName_
696  << " cannot be the same as this patch " << name
697  << exit(FatalIOError);
698  }
699 
700  switch (transform())
701  {
702  case ROTATIONAL:
703  {
704  dict.readEntry("rotationAxis", rotationAxis_);
705  dict.readEntry("rotationCentre", rotationCentre_);
706  if (dict.readIfPresent("rotationAngle", rotationAngle_))
707  {
708  rotationAngleDefined_ = true;
709  rotationAngle_ = degToRad(rotationAngle_);
710 
711  if (debug)
712  {
713  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
714  << endl;
715  }
716  }
717 
718  scalar magRot = mag(rotationAxis_);
719  if (magRot < SMALL)
720  {
722  << "Illegal rotationAxis " << rotationAxis_ << endl
723  << "Please supply a non-zero vector."
724  << exit(FatalIOError);
725  }
726  rotationAxis_ /= magRot;
727 
728  break;
729  }
730  case TRANSLATIONAL:
731  {
732  dict.readEntry("separationVector", separationVector_);
733  break;
734  }
735  default:
736  {
737  // No additional info required
738  }
739  }
740 
741  // Neighbour patch might not be valid yet so no transformation
742  // calculation possible
743 
744  // If topology change, recover the sizes of the original patches and
745  // read additional controls
746  if (createAMIFaces_)
747  {
748  srcFaceIDs_.setSize(dict.get<label>("srcSize"));
749  tgtFaceIDs_.setSize(dict.get<label>("tgtSize"));
750  moveFaceCentres_ = dict.getOrDefault("moveFaceCentres", true);
751  }
752 }
753 
754 
756 (
757  const cyclicAMIPolyPatch& pp,
758  const polyBoundaryMesh& bm
759 )
760 :
761  coupledPolyPatch(pp, bm),
762  nbrPatchName_(pp.nbrPatchName_),
763  coupleGroup_(pp.coupleGroup_),
764  nbrPatchID_(-1),
765  fraction_(pp.fraction_),
766  rotationAxis_(pp.rotationAxis_),
767  rotationCentre_(pp.rotationCentre_),
768  rotationAngleDefined_(pp.rotationAngleDefined_),
769  rotationAngle_(pp.rotationAngle_),
770  separationVector_(pp.separationVector_),
771  periodicPatchName_(pp.periodicPatchName_),
772  periodicPatchID_(-1),
773  AMIPtr_(pp.AMIPtr_->clone()),
774  surfDict_(pp.surfDict_),
775  surfPtr_(nullptr),
776  createAMIFaces_(pp.createAMIFaces_),
777  moveFaceCentres_(pp.moveFaceCentres_),
778  updatingAMI_(true),
779  srcFaceIDs_(),
780  tgtFaceIDs_(),
781  faceAreas0_(),
782  faceCentres0_()
783 {
784  // Neighbour patch might not be valid yet so no transformation
785  // calculation possible
786 }
787 
788 
790 (
791  const cyclicAMIPolyPatch& pp,
792  const polyBoundaryMesh& bm,
793  const label index,
794  const label newSize,
795  const label newStart,
796  const word& nbrPatchName
797 )
798 :
799  coupledPolyPatch(pp, bm, index, newSize, newStart),
800  nbrPatchName_(nbrPatchName),
801  coupleGroup_(pp.coupleGroup_),
802  nbrPatchID_(-1),
803  fraction_(pp.fraction_),
804  rotationAxis_(pp.rotationAxis_),
805  rotationCentre_(pp.rotationCentre_),
806  rotationAngleDefined_(pp.rotationAngleDefined_),
807  rotationAngle_(pp.rotationAngle_),
808  separationVector_(pp.separationVector_),
809  periodicPatchName_(pp.periodicPatchName_),
810  periodicPatchID_(-1),
811  AMIPtr_(pp.AMIPtr_->clone()),
812  surfDict_(pp.surfDict_),
813  surfPtr_(nullptr),
814  createAMIFaces_(pp.createAMIFaces_),
815  moveFaceCentres_(pp.moveFaceCentres_),
816  updatingAMI_(true),
817  srcFaceIDs_(),
818  tgtFaceIDs_(),
819  faceAreas0_(),
820  faceCentres0_()
821 {
822  if (nbrPatchName_ == name())
823  {
825  << "Neighbour patch name " << nbrPatchName_
826  << " cannot be the same as this patch " << name()
827  << exit(FatalError);
828  }
829 
830  // Neighbour patch might not be valid yet so no transformation
831  // calculation possible
832 }
833 
834 
836 (
837  const cyclicAMIPolyPatch& pp,
838  const polyBoundaryMesh& bm,
839  const label index,
840  const labelUList& mapAddressing,
841  const label newStart
842 )
843 :
844  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
845  nbrPatchName_(pp.nbrPatchName_),
846  coupleGroup_(pp.coupleGroup_),
847  nbrPatchID_(-1),
848  fraction_(pp.fraction_),
849  rotationAxis_(pp.rotationAxis_),
850  rotationCentre_(pp.rotationCentre_),
851  rotationAngleDefined_(pp.rotationAngleDefined_),
852  rotationAngle_(pp.rotationAngle_),
853  separationVector_(pp.separationVector_),
854  periodicPatchName_(pp.periodicPatchName_),
855  periodicPatchID_(-1),
856  AMIPtr_(pp.AMIPtr_->clone()),
857  surfDict_(pp.surfDict_),
858  surfPtr_(nullptr),
859  createAMIFaces_(pp.createAMIFaces_),
860  moveFaceCentres_(pp.moveFaceCentres_),
861  updatingAMI_(true),
862  srcFaceIDs_(),
863  tgtFaceIDs_(),
864  faceAreas0_(),
865  faceCentres0_()
866 {}
867 
868 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
869 
871 (
872  label& newFaces,
873  label& newProcFaces
874 ) const
875 {
876  const labelListList& addSourceFaces = AMI().srcAddress();
877 
878  // Add new faces as many weights for AMI
879  forAll (addSourceFaces, faceI)
880  {
881  const labelList& nbrFaceIs = addSourceFaces[faceI];
882 
883  forAll (nbrFaceIs, j)
884  {
885  label nbrFaceI = nbrFaceIs[j];
886 
887  if (nbrFaceI < neighbPatch().size())
888  {
889  // local faces
890  newFaces++;
891  }
892  else
893  {
894  // Proc faces
895  newProcFaces++;
896  }
897  }
898  }
899 }
900 
901 
903 {
904  if (nbrPatchID_ == -1)
905  {
906  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
907 
908  if (nbrPatchID_ == -1)
909  {
911  << "Illegal neighbourPatch name " << neighbPatchName()
912  << nl << "Valid patch names are "
913  << this->boundaryMesh().names()
914  << exit(FatalError);
915  }
916 
917  // Check that it is a cyclic AMI patch
918  const cyclicAMIPolyPatch& nbrPatch =
919  refCast<const cyclicAMIPolyPatch>
920  (
921  this->boundaryMesh()[nbrPatchID_]
922  );
923 
924  if (nbrPatch.neighbPatchName() != name())
925  {
927  << "Patch " << name()
928  << " specifies neighbour patch " << neighbPatchName()
929  << nl << " but that in return specifies "
930  << nbrPatch.neighbPatchName() << endl;
931  }
932  }
933 
934  return nbrPatchID_;
935 }
936 
937 
939 {
940  if (periodicPatchName_ == word::null)
941  {
942  return -1;
943  }
944  else
945  {
946  if (periodicPatchID_ == -1)
947  {
948  periodicPatchID_ = boundaryMesh().findPatchID(periodicPatchName_);
949 
950  if (periodicPatchID_ == -1)
951  {
953  << "Illegal neighbourPatch name " << periodicPatchName_
954  << nl << "Valid patch names are "
955  << this->boundaryMesh().names()
956  << exit(FatalError);
957  }
958  }
959  return periodicPatchID_;
960  }
961 }
962 
963 
965 {
966  return index() < neighbPatchID();
967 }
968 
969 
971 {
972  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
973  return refCast<const cyclicAMIPolyPatch>(pp);
974 }
975 
976 
978 {
979  if (!owner())
980  {
982  << "AMI interpolator only available to owner patch"
983  << abort(FatalError);
984  }
985 
986  if (!AMIPtr_->upToDate())
987  {
988  resetAMI();
989  }
990 
991  return *AMIPtr_;
992 }
993 
994 
996 {
997  if (owner())
998  {
999  return AMI().applyLowWeightCorrection();
1000  }
1001  else
1002  {
1003  return neighbPatch().AMI().applyLowWeightCorrection();
1004  }
1005 }
1006 
1007 
1009 {
1010  if (!parallel())
1011  {
1012  if (transform() == ROTATIONAL)
1013  {
1014  l = Foam::transform(forwardT(), l - rotationCentre_)
1015  + rotationCentre_;
1016  }
1017  else
1018  {
1019  l = Foam::transform(forwardT(), l);
1020  }
1021  }
1022  else if (separated())
1023  {
1024  // transformPosition gets called on the receiving side,
1025  // separation gets calculated on the sending side so subtract
1026 
1027  const vectorField& s = separation();
1028  if (s.size() == 1)
1029  {
1030  forAll(l, i)
1031  {
1032  l[i] -= s[0];
1033  }
1034  }
1035  else
1036  {
1037  l -= s;
1038  }
1039  }
1040 }
1041 
1042 
1045  point& l,
1046  const label facei
1047 ) const
1048 {
1049  if (!parallel())
1050  {
1051  const tensor& T =
1052  (
1053  forwardT().size() == 1
1054  ? forwardT()[0]
1055  : forwardT()[facei]
1056  );
1057 
1058  if (transform() == ROTATIONAL)
1059  {
1060  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1061  }
1062  else
1063  {
1064  l = Foam::transform(T, l);
1065  }
1066  }
1067  else if (separated())
1068  {
1069  const vector& s =
1070  (
1071  separation().size() == 1
1072  ? separation()[0]
1073  : separation()[facei]
1074  );
1075 
1076  l -= s;
1077  }
1078 }
1079 
1080 
1083  point& l,
1084  const label facei
1085 ) const
1086 {
1087  if (!parallel())
1088  {
1089  const tensor& T =
1090  (
1091  reverseT().size() == 1
1092  ? reverseT()[0]
1093  : reverseT()[facei]
1094  );
1095 
1096  if (transform() == ROTATIONAL)
1097  {
1098  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1099  }
1100  else
1101  {
1102  l = Foam::transform(T, l);
1103  }
1104  }
1105  else if (separated())
1106  {
1107  const vector& s =
1108  (
1109  separation().size() == 1
1110  ? separation()[0]
1111  : separation()[facei]
1112  );
1113 
1114  l += s;
1115  }
1116 }
1117 
1118 
1121  vector& d,
1122  const label facei
1123 ) const
1124 {
1125  if (!parallel())
1126  {
1127  const tensor& T =
1128  (
1129  reverseT().size() == 1
1130  ? reverseT()[0]
1131  : reverseT()[facei]
1132  );
1133 
1134  d = Foam::transform(T, d);
1135  }
1136 }
1137 
1138 
1141  const primitivePatch& referPatch,
1142  const pointField& thisCtrs,
1143  const vectorField& thisAreas,
1144  const pointField& thisCc,
1145  const pointField& nbrCtrs,
1146  const vectorField& nbrAreas,
1147  const pointField& nbrCc
1148 )
1149 {}
1150 
1151 
1154  PstreamBuffers& pBufs,
1155  const primitivePatch& pp
1156 ) const
1157 {}
1158 
1159 
1162  PstreamBuffers& pBufs,
1163  const primitivePatch& pp,
1164  labelList& faceMap,
1165  labelList& rotation
1166 ) const
1167 {
1168  faceMap.setSize(pp.size());
1169  faceMap = -1;
1170 
1171  rotation.setSize(pp.size());
1172  rotation = 0;
1173 
1174  return false;
1175 }
1176 
1177 
1180  const label facei,
1181  const vector& n,
1182  point& p
1183 ) const
1184 {
1185  point prt(p);
1186  reverseTransformPosition(prt, facei);
1187 
1188  vector nrt(n);
1189  reverseTransformDirection(nrt, facei);
1190 
1191  label nbrFacei = -1;
1192 
1193  if (owner())
1194  {
1195  nbrFacei = AMI().tgtPointFace
1196  (
1197  *this,
1198  neighbPatch(),
1199  nrt,
1200  facei,
1201  prt
1202  );
1203  }
1204  else
1205  {
1206  nbrFacei = neighbPatch().AMI().srcPointFace
1207  (
1208  neighbPatch(),
1209  *this,
1210  nrt,
1211  facei,
1212  prt
1213  );
1214  }
1215 
1216  if (nbrFacei >= 0)
1217  {
1218  p = prt;
1219  }
1220 
1221  return nbrFacei;
1222 }
1223 
1224 
1226 {
1228  if (!nbrPatchName_.empty())
1229  {
1230  os.writeEntry("neighbourPatch", nbrPatchName_);
1231  }
1232  coupleGroup_.write(os);
1233 
1234  switch (transform())
1235  {
1236  case ROTATIONAL:
1237  {
1238  os.writeEntry("rotationAxis", rotationAxis_);
1239  os.writeEntry("rotationCentre", rotationCentre_);
1240 
1241  if (rotationAngleDefined_)
1242  {
1243  os.writeEntry("rotationAngle", radToDeg(rotationAngle_));
1244  }
1245 
1246  break;
1247  }
1248  case TRANSLATIONAL:
1249  {
1250  os.writeEntry("separationVector", separationVector_);
1251  break;
1252  }
1253  case NOORDERING:
1254  {
1255  break;
1256  }
1257  default:
1258  {
1259  // No additional info to write
1260  }
1261  }
1262 
1263  if (periodicPatchName_ != word::null)
1264  {
1265  os.writeEntry("periodicPatch", periodicPatchName_);
1266  }
1267 
1268  AMIPtr_->write(os);
1269 
1270  if (!surfDict_.empty())
1271  {
1272  surfDict_.writeEntry(surfDict_.dictName(), os);
1273  }
1274 
1275  if (createAMIFaces_)
1276  {
1277  os.writeEntry("createAMIFaces", createAMIFaces_);
1278  os.writeEntry("srcSize", srcFaceIDs_.size());
1279  os.writeEntry("tgtSize", tgtFaceIDs_.size());
1280  os.writeEntry("moveFaceCentres", moveFaceCentres_);
1281  }
1282 
1283  os.writeEntryIfDifferent<scalar>("fraction", Zero, fraction_);
1284 }
1285 
1286 
1287 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::cyclicAMIPolyPatch::order
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
Definition: cyclicAMIPolyPatch.C:1161
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:248
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:169
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:372
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:964
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
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:63
cyclicPolyPatch.H
Foam::cyclicAMIPolyPatch::nbrPatchName_
word nbrPatchName_
Name of other half.
Definition: cyclicAMIPolyPatch.H:92
Foam::cyclicAMIPolyPatch::initOrder
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Definition: cyclicAMIPolyPatch.C:1153
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::periodicPatchName_
word periodicPatchName_
Periodic patch name.
Definition: cyclicAMIPolyPatch.H:130
Foam::cyclicAMIPolyPatch::rotationCentre_
point rotationCentre_
Point on axis of rotation for rotational cyclics.
Definition: cyclicAMIPolyPatch.H:112
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:54
Foam::cyclicAMIPolyPatch::initGeometry
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: cyclicAMIPolyPatch.C:485
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
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:1082
Foam::cyclicAMIPolyPatch::surfDict_
const dictionary surfDict_
Dictionary used during projection surface construction.
Definition: cyclicAMIPolyPatch.H:140
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:341
unitConversion.H
Unit conversion functions.
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::cyclicAMIPolyPatch::coupleGroup_
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
Definition: cyclicAMIPolyPatch.H:95
Foam::cyclicAMIPolyPatch::tolerance_
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
Definition: cyclicAMIPolyPatch.H:376
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:54
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:107
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:521
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:509
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:581
Foam::coupledPolyPatch::write
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
Definition: coupledPolyPatch.C:577
Foam::coupledPolyPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: coupledPolyPatch.H:295
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:1008
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::cyclicAMIPolyPatch::newInternalProcFaces
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
Definition: cyclicAMIPolyPatch.C:871
Foam::polyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:114
Foam::autoPtr::set
void set(T *p) noexcept
Deprecated(2018-02) Identical to reset().
Definition: autoPtr.H:288
Foam::cyclicAMIPolyPatch::separationVector_
vector separationVector_
Translation vector.
Definition: cyclicAMIPolyPatch.H:124
Foam::cyclicAMIPolyPatch::calcGeometry
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
Definition: cyclicAMIPolyPatch.C:502
Foam::Field< vector >
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1225
Foam::cyclicAMIPolyPatch::AMIPtr_
autoPtr< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
Definition: cyclicAMIPolyPatch.H:137
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::cyclicAMIPolyPatch::pointFace
label pointFace(const label facei, const vector &n, point &p) const
Definition: cyclicAMIPolyPatch.C:1179
Foam::cyclicAMIPolyPatch::movePoints
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
Definition: cyclicAMIPolyPatch.C:545
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:302
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:99
Foam::cyclicAMIPolyPatch::calcTransforms
virtual void calcTransforms()
Recalculate the transformation tensors.
Definition: cyclicAMIPolyPatch.C:449
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
Foam::cyclicAMIPolyPatch::rotationAxis_
vector rotationAxis_
Axis of rotation for rotational cyclics.
Definition: cyclicAMIPolyPatch.H:109
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicAMIPolyPatch::rotationAngle_
scalar rotationAngle_
Rotation angle.
Definition: cyclicAMIPolyPatch.H:118
Foam::cyclicAMIPolyPatch::neighbPatchID
virtual label neighbPatchID() const
Neighbour patch ID.
Definition: cyclicAMIPolyPatch.C:902
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:123
Foam::cyclicAMIPolyPatch::neighbPatch
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Definition: cyclicAMIPolyPatch.C:970
Foam::word::valid
static bool valid(char c)
Is this character valid for a word?
Definition: wordI.H:59
os
OBJstream os(runTime.globalPath()/outputName)
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:317
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:606
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:115
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:393
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
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:359
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:153
Foam::cyclicAMIPolyPatch::fraction_
const scalar fraction_
Particle displacement fraction accross AMI.
Definition: cyclicAMIPolyPatch.H:101
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::cyclicAMIPolyPatch::AMI
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
Definition: cyclicAMIPolyPatch.C:977
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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:321
Foam::cyclicAMIPolyPatch::applyLowWeightCorrection
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
Definition: cyclicAMIPolyPatch.C:995
Foam::cyclicAMIPolyPatch::initUpdateMesh
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: cyclicAMIPolyPatch.C:568
Foam::cyclicAMIPolyPatch::clearGeom
virtual void clearGeom()
Clear geometry.
Definition: cyclicAMIPolyPatch.C:590
Foam::Vector< scalar >
Foam::findMax
label findMax(const ListType &input, label start=0)
Foam::List< labelList >
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.
Definition: cyclicAMIPolyPatch.C:1120
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:80
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
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:473
Foam::PrimitivePatch::movePoints
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
Definition: PrimitivePatch.C:171
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::cyclicAMIPolyPatch::createAMIFaces_
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
Definition: cyclicAMIPolyPatch.H:150
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:60
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::cyclicAMIPolyPatch::cylindricalCS
autoPtr< coordSystem::cylindrical > cylindricalCS() const
Create a coordinate system from the periodic patch (or nullptr)
Definition: cyclicAMIPolyPatch.C:291
Foam::cyclicAMIPolyPatch::periodicPatchID
label periodicPatchID() const
Periodic patch ID (or -1)
Definition: cyclicAMIPolyPatch.C:938
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
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:405
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:68