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-2023 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()).clear();
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()).clear();
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 {
293  const label periodicID = periodicPatchID();
294  if (periodicID != -1)
295  {
296  // Get the periodic patch
297  const coupledPolyPatch& perPp =
298  refCast<const coupledPolyPatch>(boundaryMesh()[periodicID]);
299 
300  if (!perPp.parallel())
301  {
302  vector axis(Zero);
303  point axisPoint(Zero);
304 
305  if
306  (
307  const cyclicPolyPatch* cpp = isA<cyclicPolyPatch>(perPp)
308  )
309  {
310  axis = cpp->rotationAxis();
311  axisPoint = cpp->rotationCentre();
312  }
313  else if
314  (
315  const cyclicAMIPolyPatch* cpp = isA<cyclicAMIPolyPatch>(perPp)
316  )
317  {
318  axis = cpp->rotationAxis();
319  axisPoint = cpp->rotationCentre();
320  }
321  else
322  {
323  FatalErrorInFunction << "On patch " << name()
324  << " have unsupported periodicPatch " << perPp.name()
325  << exit(FatalError);
326  }
327 
328  return autoPtr<coordSystem::cylindrical>::New(axisPoint, axis);
329  }
330  }
331  return nullptr;
332 }
333 
334 
335 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
336 
339 {
340  const word surfType(surfDict_.getOrDefault<word>("type", "none"));
341 
342  if (!surfPtr_ && owner() && surfType != "none")
343  {
344  word surfName(surfDict_.getOrDefault("name", name()));
345 
346  const polyMesh& mesh = boundaryMesh().mesh();
347 
348  surfPtr_ =
350  (
351  surfType,
352  IOobject
353  (
354  surfName,
355  mesh.time().constant(),
356  "triSurface",
357  mesh,
360  ),
361  surfDict_
362  );
363  }
364 
365  return surfPtr_;
366 }
367 
370 {
371  resetAMI(boundaryMesh().mesh().points());
372 }
373 
374 
376 {
378 
379  if (!owner())
380  {
381  return;
382  }
383 
384  const cyclicAMIPolyPatch& nbr = neighbPatch();
385  const pointField srcPoints(points, meshPoints());
386  pointField nbrPoints(points, nbr.meshPoints());
387 
388  Info<< "AMI: Creating AMI for source:" << name()
389  << " and target:" << nbr.name() << endl;
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  // Note: e.g. redistributePar might construct the AMI with a different
441  // worldComm so reset it to the mesh.comm.
442  AMIPtr_->comm(boundaryMesh().mesh().comm());
443  AMIPtr_->calculate(patch0, nbrPatch0, surfPtr());
444 
445  if (debug)
446  {
447  AMIPtr_->checkSymmetricWeights(true);
448  }
449 }
450 
451 
453 {
455 
456  const cyclicAMIPolyPatch& half0 = *this;
457  vectorField half0Areas(half0.size());
458  forAll(half0, facei)
459  {
460  half0Areas[facei] = half0[facei].areaNormal(half0.points());
461  }
462 
463  const cyclicAMIPolyPatch& half1 = neighbPatch();
464  vectorField half1Areas(half1.size());
465  forAll(half1, facei)
466  {
467  half1Areas[facei] = half1[facei].areaNormal(half1.points());
468  }
469 
470  calcTransforms
471  (
472  half0,
473  half0.faceCentres(),
474  half0Areas,
475  half1.faceCentres(),
476  half1Areas
477  );
478 
479  DebugPout
480  << "calcTransforms() : patch: " << name() << nl
481  << " forwardT = " << forwardT() << nl
482  << " reverseT = " << reverseT() << nl
483  << " separation = " << separation() << nl
484  << " collocated = " << collocated() << nl << endl;
485 }
486 
487 
489 {
491 
492  // Flag AMI as needing update
493  AMIPtr_->upToDate(false);
494 
496 
497  // Early calculation of transforms so e.g. cyclicACMI can use them.
498  // Note: also triggers primitiveMesh face centre. Note that cell
499  // centres should -not- be calculated
500  // since e.g. cyclicACMI overrides face areas
501  calcTransforms();
502 }
503 
504 
506 {
508 }
509 
510 
512 (
513  PstreamBuffers& pBufs,
514  const pointField& p
515 )
516 {
518 
519  // See below. Clear out any local geometry
521 
522  // Note: processorPolyPatch::initMovePoints calls
523  // processorPolyPatch::initGeometry which will trigger calculation of
524  // patch faceCentres() and cell volumes...
525 
526  if (createAMIFaces_)
527  {
528  // Note: AMI should have been updated in setTopology
529 
530  // faceAreas() and faceCentres() have been reset and will be
531  // recalculated on-demand using the mesh points and no longer
532  // correspond to the scaled areas!
533  restoreScaledGeometry();
534 
535  // deltas need to be recalculated to use new face centres!
536  }
537  else
538  {
539  AMIPtr_->upToDate(false);
540  }
542  // Early calculation of transforms. See above.
543  calcTransforms();
544 }
545 
546 
548 (
549  PstreamBuffers& pBufs,
550  const pointField& p
551 )
552 {
554 
555 // Note: not calling movePoints since this will undo our manipulations!
556 // polyPatch::movePoints(pBufs, p);
557 /*
558  polyPatch::movePoints
559  -> primitivePatch::movePoints
560  -> primitivePatch::clearGeom:
561  deleteDemandDrivenData(localPointsPtr_);
562  deleteDemandDrivenData(faceCentresPtr_);
563  deleteDemandDrivenData(faceAreasPtr_);
564  deleteDemandDrivenData(magFaceAreasPtr_);
565  deleteDemandDrivenData(faceNormalsPtr_);
566  deleteDemandDrivenData(pointNormalsPtr_);
567 */
568 }
569 
570 
572 {
574 
576 
577  if (createAMIFaces_ && boundaryMesh().mesh().topoChanging() && owner())
578  {
579  setAMIFaces();
580  }
581 }
582 
583 
584 void Foam::cyclicAMIPolyPatch::updateMesh(PstreamBuffers& pBufs)
585 {
587 
588  // Note: this clears out cellCentres(), faceCentres() and faceAreas()
589  polyPatch::updateMesh(pBufs);
590 }
591 
592 
594 {
596 
597  if (!updatingAMI_)
598  {
599  AMIPtr_->upToDate(false);
600  }
601 
603 }
604 
605 
606 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
607 
609 (
610  const word& name,
611  const label size,
612  const label start,
613  const label index,
614  const polyBoundaryMesh& bm,
615  const word& patchType,
616  const transformType transform,
617  const word& defaultAMIMethod
618 )
619 :
620  coupledPolyPatch(name, size, start, index, bm, patchType, transform),
621  nbrPatchName_(),
622  nbrPatchID_(-1),
623  fraction_(Zero),
624  rotationAxis_(Zero),
625  rotationCentre_(Zero),
626  rotationAngleDefined_(false),
627  rotationAngle_(0.0),
628  separationVector_(Zero),
629  periodicPatchName_(),
630  periodicPatchID_(-1),
631  AMIPtr_(AMIInterpolation::New(defaultAMIMethod)),
632  surfDict_(fileName("surface")),
633  surfPtr_(nullptr),
634  createAMIFaces_(false),
635  moveFaceCentres_(false),
636  updatingAMI_(true),
637  srcFaceIDs_(),
638  tgtFaceIDs_(),
639  faceAreas0_(),
640  faceCentres0_()
641 {
642  // Neighbour patch might not be valid yet so no transformation
643  // calculation possible
644 }
645 
646 
648 (
649  const word& name,
650  const dictionary& dict,
651  const label index,
652  const polyBoundaryMesh& bm,
653  const word& patchType,
654  const word& defaultAMIMethod
655 )
656 :
657  coupledPolyPatch(name, dict, index, bm, patchType),
658  nbrPatchName_(dict.getOrDefault<word>("neighbourPatch", word::null)),
659  coupleGroup_(dict),
660  nbrPatchID_(-1),
661  fraction_(dict.getOrDefault<scalar>("fraction", Zero)),
662  rotationAxis_(Zero),
663  rotationCentre_(Zero),
664  rotationAngleDefined_(false),
665  rotationAngle_(0.0),
666  separationVector_(Zero),
667  periodicPatchName_(dict.getOrDefault<word>("periodicPatch", word::null)),
668  periodicPatchID_(-1),
669  AMIPtr_
670  (
672  (
673  dict.getOrDefault<word>("AMIMethod", defaultAMIMethod),
674  dict,
675  dict.getOrDefault("flipNormals", false)
676  )
677  ),
678  surfDict_(dict.subOrEmptyDict("surface")),
679  surfPtr_(nullptr),
680  createAMIFaces_(dict.getOrDefault("createAMIFaces", false)),
681  moveFaceCentres_(false),
682  updatingAMI_(true),
683  srcFaceIDs_(),
684  tgtFaceIDs_(),
685  faceAreas0_(),
686  faceCentres0_()
687 {
688  if (nbrPatchName_.empty() && !coupleGroup_.good())
689  {
691  << "No \"neighbourPatch\" or \"coupleGroup\" provided."
692  << exit(FatalIOError);
693  }
694 
695  if (nbrPatchName_ == name)
696  {
698  << "Neighbour patch name " << nbrPatchName_
699  << " cannot be the same as this patch " << name
700  << exit(FatalIOError);
701  }
702 
703  switch (transform())
704  {
705  case ROTATIONAL:
706  {
707  dict.readEntry("rotationAxis", rotationAxis_);
708  dict.readEntry("rotationCentre", rotationCentre_);
709  if (dict.readIfPresent("rotationAngle", rotationAngle_))
710  {
711  rotationAngleDefined_ = true;
713 
714  if (debug)
715  {
716  Info<< "rotationAngle: " << rotationAngle_ << " [rad]"
717  << endl;
718  }
719  }
720 
721  scalar magRot = mag(rotationAxis_);
722  if (magRot < SMALL)
723  {
725  << "Illegal rotationAxis " << rotationAxis_ << endl
726  << "Please supply a non-zero vector."
727  << exit(FatalIOError);
728  }
729  rotationAxis_ /= magRot;
730 
731  break;
732  }
733  case TRANSLATIONAL:
734  {
735  dict.readEntry("separationVector", separationVector_);
736  break;
737  }
738  default:
739  {
740  // No additional info required
741  }
742  }
743 
744  // Neighbour patch might not be valid yet so no transformation
745  // calculation possible
746 
747  // If topology change, recover the sizes of the original patches and
748  // read additional controls
749  if (createAMIFaces_)
750  {
751  srcFaceIDs_.setSize(dict.get<label>("srcSize"));
752  tgtFaceIDs_.setSize(dict.get<label>("tgtSize"));
753  moveFaceCentres_ = dict.getOrDefault("moveFaceCentres", true);
754  }
755 }
756 
757 
759 (
760  const cyclicAMIPolyPatch& pp,
761  const polyBoundaryMesh& bm
762 )
763 :
764  coupledPolyPatch(pp, bm),
765  nbrPatchName_(pp.nbrPatchName_),
766  coupleGroup_(pp.coupleGroup_),
767  nbrPatchID_(-1),
768  fraction_(pp.fraction_),
769  rotationAxis_(pp.rotationAxis_),
770  rotationCentre_(pp.rotationCentre_),
771  rotationAngleDefined_(pp.rotationAngleDefined_),
772  rotationAngle_(pp.rotationAngle_),
773  separationVector_(pp.separationVector_),
774  periodicPatchName_(pp.periodicPatchName_),
775  periodicPatchID_(-1),
776  AMIPtr_(pp.AMIPtr_->clone()),
777  surfDict_(pp.surfDict_),
778  surfPtr_(nullptr),
779  createAMIFaces_(pp.createAMIFaces_),
780  moveFaceCentres_(pp.moveFaceCentres_),
781  updatingAMI_(true),
782  srcFaceIDs_(),
783  tgtFaceIDs_(),
784  faceAreas0_(),
785  faceCentres0_()
786 {
787  // Neighbour patch might not be valid yet so no transformation
788  // calculation possible
789 }
790 
791 
793 (
794  const cyclicAMIPolyPatch& pp,
795  const polyBoundaryMesh& bm,
796  const label index,
797  const label newSize,
798  const label newStart,
799  const word& nbrPatchName
800 )
801 :
802  coupledPolyPatch(pp, bm, index, newSize, newStart),
803  nbrPatchName_(nbrPatchName),
804  coupleGroup_(pp.coupleGroup_),
805  nbrPatchID_(-1),
806  fraction_(pp.fraction_),
807  rotationAxis_(pp.rotationAxis_),
808  rotationCentre_(pp.rotationCentre_),
809  rotationAngleDefined_(pp.rotationAngleDefined_),
810  rotationAngle_(pp.rotationAngle_),
811  separationVector_(pp.separationVector_),
812  periodicPatchName_(pp.periodicPatchName_),
813  periodicPatchID_(-1),
814  AMIPtr_(pp.AMIPtr_->clone()),
815  surfDict_(pp.surfDict_),
816  surfPtr_(nullptr),
817  createAMIFaces_(pp.createAMIFaces_),
818  moveFaceCentres_(pp.moveFaceCentres_),
819  updatingAMI_(true),
820  srcFaceIDs_(),
821  tgtFaceIDs_(),
822  faceAreas0_(),
823  faceCentres0_()
824 {
825  if (nbrPatchName_ == name())
826  {
828  << "Neighbour patch name " << nbrPatchName_
829  << " cannot be the same as this patch " << name()
830  << exit(FatalError);
831  }
833  // Neighbour patch might not be valid yet so no transformation
834  // calculation possible
835 }
836 
837 
839 (
840  const cyclicAMIPolyPatch& pp,
841  const polyBoundaryMesh& bm,
842  const label index,
843  const labelUList& mapAddressing,
844  const label newStart
845 )
846 :
847  coupledPolyPatch(pp, bm, index, mapAddressing, newStart),
848  nbrPatchName_(pp.nbrPatchName_),
849  coupleGroup_(pp.coupleGroup_),
850  nbrPatchID_(-1),
851  fraction_(pp.fraction_),
852  rotationAxis_(pp.rotationAxis_),
853  rotationCentre_(pp.rotationCentre_),
854  rotationAngleDefined_(pp.rotationAngleDefined_),
855  rotationAngle_(pp.rotationAngle_),
856  separationVector_(pp.separationVector_),
857  periodicPatchName_(pp.periodicPatchName_),
858  periodicPatchID_(-1),
859  AMIPtr_(pp.AMIPtr_->clone()),
860  surfDict_(pp.surfDict_),
861  surfPtr_(nullptr),
862  createAMIFaces_(pp.createAMIFaces_),
863  moveFaceCentres_(pp.moveFaceCentres_),
864  updatingAMI_(true),
865  srcFaceIDs_(),
866  tgtFaceIDs_(),
867  faceAreas0_(),
868  faceCentres0_()
869 {}
870 
871 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
872 
874 (
875  label& newFaces,
876  label& newProcFaces
877 ) const
878 {
879  const labelListList& addSourceFaces = AMI().srcAddress();
880 
881  // Add new faces as many weights for AMI
882  forAll (addSourceFaces, faceI)
883  {
884  const labelList& nbrFaceIs = addSourceFaces[faceI];
885 
886  forAll (nbrFaceIs, j)
887  {
888  label nbrFaceI = nbrFaceIs[j];
889 
890  if (nbrFaceI < neighbPatch().size())
891  {
892  // local faces
893  newFaces++;
894  }
895  else
896  {
897  // Proc faces
898  newProcFaces++;
899  }
900  }
901  }
902 }
903 
904 
905 Foam::label Foam::cyclicAMIPolyPatch::neighbPatchID() const
906 {
907  if (nbrPatchID_ == -1)
908  {
909  nbrPatchID_ = this->boundaryMesh().findPatchID(neighbPatchName());
910 
911  if (nbrPatchID_ == -1)
912  {
914  << "Illegal neighbourPatch name " << neighbPatchName()
915  << nl << "Valid patch names are "
916  << this->boundaryMesh().names()
917  << exit(FatalError);
918  }
919 
920  // Check that it is a cyclic AMI patch
921  const cyclicAMIPolyPatch& nbrPatch =
922  refCast<const cyclicAMIPolyPatch>
923  (
924  this->boundaryMesh()[nbrPatchID_]
925  );
926 
927  if (nbrPatch.neighbPatchName() != name())
928  {
930  << "Patch " << name()
931  << " specifies neighbour patch " << neighbPatchName()
932  << nl << " but that in return specifies "
933  << nbrPatch.neighbPatchName() << endl;
934  }
935  }
936 
937  return nbrPatchID_;
938 }
939 
940 
942 {
943  if (periodicPatchName_.empty())
944  {
945  return -1;
946  }
947  else
948  {
949  if (periodicPatchID_ == -1)
950  {
951  periodicPatchID_ = boundaryMesh().findPatchID(periodicPatchName_);
952 
953  if (periodicPatchID_ == -1)
954  {
956  << "Illegal neighbourPatch name " << periodicPatchName_
957  << nl << "Valid patch names are "
958  << this->boundaryMesh().names()
959  << exit(FatalError);
960  }
961  }
962  return periodicPatchID_;
963  }
964 }
965 
968 {
969  return index() < neighbPatchID();
970 }
971 
972 
974 {
975  const polyPatch& pp = this->boundaryMesh()[neighbPatchID()];
976  return refCast<const cyclicAMIPolyPatch>(pp);
977 }
978 
979 
981 {
982  if (!owner())
983  {
985  << "AMI interpolator only available to owner patch"
986  << abort(FatalError);
987  }
988 
989  if (!AMIPtr_->upToDate())
990  {
991  resetAMI();
992  }
993 
994  return *AMIPtr_;
995 }
996 
997 
999 {
1000  if (owner())
1001  {
1002  return AMI().applyLowWeightCorrection();
1003  }
1004  else
1005  {
1006  return neighbPatch().AMI().applyLowWeightCorrection();
1007  }
1008 }
1009 
1010 
1012 {
1013  if (!parallel())
1014  {
1015  if (transform() == ROTATIONAL)
1016  {
1017  l = Foam::transform(forwardT(), l - rotationCentre_)
1018  + rotationCentre_;
1019  }
1020  else
1021  {
1022  l = Foam::transform(forwardT(), l);
1023  }
1024  }
1025  else if (separated())
1026  {
1027  // transformPosition gets called on the receiving side,
1028  // separation gets calculated on the sending side so subtract
1029 
1030  const vectorField& s = separation();
1031  if (s.size() == 1)
1032  {
1033  forAll(l, i)
1034  {
1035  l[i] -= s[0];
1036  }
1037  }
1038  else
1039  {
1040  l -= s;
1041  }
1042  }
1043 }
1044 
1045 
1047 (
1048  point& l,
1049  const label facei
1050 ) const
1051 {
1052  if (!parallel())
1053  {
1054  const tensor& T =
1055  (
1056  forwardT().size() == 1
1057  ? forwardT()[0]
1058  : forwardT()[facei]
1059  );
1060 
1061  if (transform() == ROTATIONAL)
1062  {
1063  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1064  }
1065  else
1066  {
1067  l = Foam::transform(T, l);
1068  }
1069  }
1070  else if (separated())
1071  {
1072  const vector& s =
1073  (
1074  separation().size() == 1
1075  ? separation()[0]
1076  : separation()[facei]
1077  );
1079  l -= s;
1080  }
1081 }
1082 
1083 
1085 (
1086  point& l,
1087  const label facei
1088 ) const
1089 {
1090  if (!parallel())
1091  {
1092  const tensor& T =
1093  (
1094  reverseT().size() == 1
1095  ? reverseT()[0]
1096  : reverseT()[facei]
1097  );
1098 
1099  if (transform() == ROTATIONAL)
1100  {
1101  l = Foam::transform(T, l - rotationCentre_) + rotationCentre_;
1102  }
1103  else
1104  {
1105  l = Foam::transform(T, l);
1106  }
1107  }
1108  else if (separated())
1109  {
1110  const vector& s =
1111  (
1112  separation().size() == 1
1113  ? separation()[0]
1114  : separation()[facei]
1115  );
1117  l += s;
1118  }
1119 }
1120 
1121 
1123 (
1124  vector& d,
1125  const label facei
1126 ) const
1127 {
1128  if (!parallel())
1129  {
1130  const tensor& T =
1131  (
1132  reverseT().size() == 1
1133  ? reverseT()[0]
1134  : reverseT()[facei]
1135  );
1137  d = Foam::transform(T, d);
1138  }
1139 }
1140 
1141 
1143 (
1144  const primitivePatch& referPatch,
1145  const pointField& thisCtrs,
1146  const vectorField& thisAreas,
1147  const pointField& thisCc,
1148  const pointField& nbrCtrs,
1149  const vectorField& nbrAreas,
1150  const pointField& nbrCc
1151 )
1152 {}
1153 
1154 
1156 (
1158  const primitivePatch& pp
1159 ) const
1160 {}
1161 
1162 
1164 (
1165  PstreamBuffers& pBufs,
1166  const primitivePatch& pp,
1167  labelList& faceMap,
1168  labelList& rotation
1169 ) const
1170 {
1171  faceMap.setSize(pp.size());
1172  faceMap = -1;
1173 
1174  rotation.setSize(pp.size());
1175  rotation = 0;
1176 
1177  return false;
1178 }
1179 
1180 
1182 (
1183  const label facei,
1184  const vector& n,
1185  point& p
1186 ) const
1187 {
1188  point prt(p);
1189  reverseTransformPosition(prt, facei);
1190 
1191  vector nrt(n);
1192  reverseTransformDirection(nrt, facei);
1193 
1194  label nbrFacei = -1;
1195 
1196  if (owner())
1197  {
1198  nbrFacei = AMI().tgtPointFace
1199  (
1200  *this,
1201  neighbPatch(),
1202  nrt,
1203  facei,
1204  prt
1205  );
1206  }
1207  else
1208  {
1209  nbrFacei = neighbPatch().AMI().srcPointFace
1210  (
1211  neighbPatch(),
1212  *this,
1213  nrt,
1214  facei,
1215  prt
1216  );
1217  }
1218 
1219  if (nbrFacei >= 0)
1220  {
1221  p = prt;
1222  }
1223 
1224  return nbrFacei;
1225 }
1226 
1227 
1228 void Foam::cyclicAMIPolyPatch::write(Ostream& os) const
1229 {
1231  if (!nbrPatchName_.empty())
1232  {
1233  os.writeEntry("neighbourPatch", nbrPatchName_);
1234  }
1235  coupleGroup_.write(os);
1236 
1237  switch (transform())
1238  {
1239  case ROTATIONAL:
1240  {
1241  os.writeEntry("rotationAxis", rotationAxis_);
1242  os.writeEntry("rotationCentre", rotationCentre_);
1243 
1244  if (rotationAngleDefined_)
1245  {
1246  os.writeEntry("rotationAngle", radToDeg(rotationAngle_));
1247  }
1248 
1249  break;
1250  }
1251  case TRANSLATIONAL:
1252  {
1253  os.writeEntry("separationVector", separationVector_);
1254  break;
1255  }
1256  case NOORDERING:
1257  {
1258  break;
1259  }
1260  default:
1261  {
1262  // No additional info to write
1263  }
1264  }
1265 
1266  if (!periodicPatchName_.empty())
1267  {
1268  os.writeEntry("periodicPatch", periodicPatchName_);
1269  }
1270 
1271  AMIPtr_->write(os);
1272 
1273  if (!surfDict_.empty())
1274  {
1275  surfDict_.writeEntry(surfDict_.dictName(), os);
1276  }
1277 
1278  if (createAMIFaces_)
1279  {
1280  os.writeEntry("createAMIFaces", createAMIFaces_);
1281  os.writeEntry("srcSize", srcFaceIDs_.size());
1282  os.writeEntry("tgtSize", tgtFaceIDs_.size());
1283  os.writeEntry("moveFaceCentres", moveFaceCentres_);
1284  }
1285 
1286  os.writeEntryIfDifferent<scalar>("fraction", Zero, fraction_);
1287 }
1288 
1289 
1290 // ************************************************************************* //
virtual void clearGeom()
Clear geometry.
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
vector separationVector_
Translation vector.
dictionary dict
dimensionedScalar acos(const dimensionedScalar &ds)
const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
virtual label neighbPatchID() const
Neighbour patch ID.
virtual void write(Ostream &os) const
Write the polyPatch data as a dictionary.
virtual bool order(PstreamBuffers &, const primitivePatch &, labelList &faceMap, labelList &rotation) const
Return new ordering for primitivePatch.
A class for handling file names.
Definition: fileName.H:72
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.
autoPtr< coordSystem::cylindrical > cylindricalCS() const
Create a coordinate system from the periodic patch (or nullptr)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual void initUpdateMesh(PstreamBuffers &)
Initialise the update of the patch topology.
Definition: polyPatch.H:135
bool moveFaceCentres_
Move face centres (default = no)
Unit conversion functions.
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static const scalar tolerance_
Tolerance used e.g. for area calculations/limits.
virtual void calcTransforms()
Recalculate the transformation tensors.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
dimensionedScalar sqrt(const dimensionedScalar &ds)
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void transformPosition(pointField &) const
Transform patch-based positions from nbr side to this side.
bool applyLowWeightCorrection() const
Return true if applying the low weight correction.
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
vector rotationAxis_
Axis of rotation for rotational cyclics.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
bool rotationAngleDefined_
Flag to show whether the rotation angle is defined.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to nbr side.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
virtual void newInternalProcFaces(label &, label &) const
Return number of new internal of this polyPatch faces.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual void clearGeom()
Clear geometry.
Definition: polyPatch.C:66
bool createAMIFaces_
Flag to indicate that new AMI faces will created.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void calcGeometry(PstreamBuffers &)
Calculate the patch geometry.
virtual void movePoints(const Field< point_type > &)
Correct patch after moving points.
const coupleGroupIdentifier coupleGroup_
Optional patchGroup to find neighbPatch.
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
bool good() const noexcept
The patchGroup has a non-empty name.
const pointField & points
static const Identity< scalar > I
Definition: Identity.H:100
A class for handling words, derived from Foam::string.
Definition: word.H:63
const autoPtr< searchableSurface > & surfPtr() const
Create and return pointer to the projection surface.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
#define DebugInFunction
Report an information message using Foam::Info.
point rotationCentre_
Point on axis of rotation for rotational cyclics.
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
patchWriters clear()
virtual void updateMesh(PstreamBuffers &)
Update of the patch topology.
Definition: polyPatch.C:59
Vector< scalar > vector
Definition: vector.H:57
Cyclic patch for Arbitrary Mesh Interface (AMI)
label periodicPatchID() const
Periodic patch ID (or -1)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
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:336
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:105
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const word & name() const noexcept
The patch name.
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
virtual void initGeometry(PstreamBuffers &)
Initialise the calculation of the patch geometry.
Definition: polyPatch.H:112
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
virtual void initOrder(PstreamBuffers &, const primitivePatch &) const
Initialize ordering for primitivePatch. Does not refer to *this (except for name() and type() etc...
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
constexpr scalar radToDeg(const scalar rad) noexcept
Conversion from radians to degrees.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
virtual bool owner() const
Does this side own the patch?
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
#define DebugPout
Report an information message using Foam::Pout.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
virtual transformType transform() const
Type of transform.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const bMesh & mesh() const
Definition: boundaryMesh.H:259
volScalarField & p
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
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))
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
List< bool > boolList
A List of bools.
Definition: List.H:60
word nbrPatchName_
Name of other half.
virtual void reverseTransformPosition(point &l, const label facei) const
Transform a patch-based position from this side to nbr side.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
scalar rotationAngle_
Rotation angle.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
label pointFace(const label facei, const vector &n, point &p) const
Return face index on neighbour patch which shares point p following trajectory vector n...
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual void movePoints(PstreamBuffers &pBufs, const pointField &)
Correct patches after moving points.
virtual void initMovePoints(PstreamBuffers &pBufs, const pointField &)
Initialise the patches for moving points.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:315