cyclicPeriodicAMIPolyPatch.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 OpenCFD Ltd.
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 
29 // For debugging
30 #include "OBJstream.H"
31 #include "PatchTools.H"
32 #include "Time.H"
33 //Note: cannot use vtkSurfaceWriter here - circular linkage
34 //#include "vtkSurfaceWriter.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(cyclicPeriodicAMIPolyPatch, 0);
41 
42  addToRunTimeSelectionTable(polyPatch, cyclicPeriodicAMIPolyPatch, word);
44  (
45  polyPatch,
46  cyclicPeriodicAMIPolyPatch,
47  dictionary
48  );
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
55 {
56  if (owner())
57  {
58  // At this point we guarantee that the transformations have been
59  // updated. There is one particular case, where if the periodic patch
60  // locally has zero faces then its transformation will not be set. We
61  // need to synchronise the transforms over the zero-sized patches as
62  // well.
63  //
64  // We can't put the logic into cyclicPolyPatch as
65  // processorCyclicPolyPatch uses cyclicPolyPatch functionality.
66  // Synchronisation will fail because processor-type patches do not exist
67  // on all processors.
68  //
69  // The use in cyclicPeriodicAMI is special; we use the patch even
70  // though we have no faces in it. Ideally, in future, the transformation
71  // logic will be abstracted out, and we won't need a periodic patch
72  // here. Until then, we can just fix the behaviour of the zero-sized
73  // coupled patches here
74 
75  // Get the periodic patch
76  const coupledPolyPatch& periodicPatch
77  (
78  refCast<const coupledPolyPatch>
79  (
81  )
82  );
83 
84  // If there are any zero-sized periodic patches
85  if (returnReduce((size() && !periodicPatch.size()), orOp<bool>()))
86  {
87  if (periodicPatch.separation().size() > 1)
88  {
90  << "Periodic patch " << periodicPatchName_
91  << " has non-uniform separation vector "
92  << periodicPatch.separation()
93  << "This is not allowed inside " << type()
94  << " patch " << name()
95  << exit(FatalError);
96  }
97 
98  if (periodicPatch.forwardT().size() > 1)
99  {
101  << "Periodic patch " << periodicPatchName_
102  << " has non-uniform transformation tensor "
103  << periodicPatch.forwardT()
104  << "This is not allowed inside " << type()
105  << " patch " << name()
106  << exit(FatalError);
107  }
108 
109  // Note that zero-sized patches will have zero-sized fields for the
110  // separation vector, forward and reverse transforms. These need
111  // replacing with the transformations from other processors.
112 
113  // Parallel in this context refers to a parallel transformation,
114  // rather than a rotational transformation.
115 
116  // Note that a cyclic with zero faces is considered parallel so
117  // explicitly check for that.
118  bool isParallel =
119  (
120  periodicPatch.size()
121  && periodicPatch.parallel()
122  );
123  reduce(isParallel, orOp<bool>());
124 
125  if (isParallel)
126  {
127  // Sync a list of separation vectors
129  sep[Pstream::myProcNo()] = periodicPatch.separation();
130  Pstream::gatherList(sep);
132 
134  coll[Pstream::myProcNo()] = periodicPatch.collocated();
135  Pstream::gatherList(coll);
136  Pstream::scatterList(coll);
137 
138  // If locally we have zero faces pick the first one that has a
139  // separation vector
140  if (!periodicPatch.size())
141  {
142  forAll(sep, procI)
143  {
144  if (sep[procI].size())
145  {
146  const_cast<vectorField&>
147  (
148  periodicPatch.separation()
149  ) = sep[procI];
150  const_cast<boolList&>
151  (
152  periodicPatch.collocated()
153  ) = coll[procI];
154 
155  break;
156  }
157  }
158  }
159  }
160  else
161  {
162  // Sync a list of forward and reverse transforms
164  forwardT[Pstream::myProcNo()] = periodicPatch.forwardT();
167 
169  reverseT[Pstream::myProcNo()] = periodicPatch.reverseT();
172 
173  // If locally we have zero faces pick the first one that has a
174  // transformation vector
175  if (!periodicPatch.size())
176  {
177  forAll(forwardT, procI)
178  {
179  if (forwardT[procI].size())
180  {
181  const_cast<tensorField&>
182  (
183  periodicPatch.forwardT()
184  ) = forwardT[procI];
185  const_cast<tensorField&>
186  (
187  periodicPatch.reverseT()
188  ) = reverseT[procI];
189 
190  break;
191  }
192  }
193  }
194  }
195  }
196  }
197 }
198 
199 
201 (
202  const primitivePatch& p,
203  OBJstream& str
204 ) const
205 {
206  // Collect faces and points
208  faceList allFaces;
209  labelList pointMergeMap;
211  (
212  -1.0, // do not merge points
213  p,
214  allPoints,
215  allFaces,
216  pointMergeMap
217  );
218 
219  if (Pstream::master())
220  {
221  // Write base geometry
222  str.write(allFaces, allPoints);
223  }
224 }
225 
226 
228 (
230 ) const
231 {
232  if (owner())
233  {
234  // Get the periodic patch
235  const coupledPolyPatch& periodicPatch
236  (
237  refCast<const coupledPolyPatch>
238  (
239  boundaryMesh()[periodicPatchID()]
240  )
241  );
242 
243  // Synchronise the transforms
244  syncTransforms();
245 
246  // Create copies of both patches' points, transformed to the owner
247  pointField thisPoints0(localPoints());
248  pointField nbrPoints0(neighbPatch().localPoints());
249  transformPosition(nbrPoints0);
250 
251  // Reset the stored number of periodic transformations to a lower
252  // absolute value if possible
253  if (nSectors_ > 0)
254  {
255  if (nTransforms_ > nSectors_/2)
256  {
257  nTransforms_ -= nSectors_;
258  }
259  else if (nTransforms_ < - nSectors_/2)
260  {
261  nTransforms_ += nSectors_;
262  }
263  }
264 
265  // Apply the stored number of periodic transforms
266  for (label i = 0; i < nTransforms_; ++ i)
267  {
268  periodicPatch.transformPosition(thisPoints0);
269  }
270  for (label i = 0; i > nTransforms_; -- i)
271  {
272  periodicPatch.transformPosition(nbrPoints0);
273  }
274 
275  autoPtr<OBJstream> ownStr;
276  autoPtr<OBJstream> neiStr;
277  if (debug)
278  {
279  const Time& runTime = boundaryMesh().mesh().time();
280 
281  fileName dir(runTime.rootPath()/runTime.globalCaseName());
282  fileName postfix("_" + runTime.timeName()+"_expanded.obj");
283 
284  ownStr.reset(new OBJstream(dir/name() + postfix));
285  neiStr.reset(new OBJstream(dir/neighbPatch().name() + postfix));
286 
287  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:" << name()
288  << " writing accumulated AMI to " << ownStr().name()
289  << " and " << neiStr().name() << endl;
290  }
291 
292  // Create another copy
293  pointField thisPoints(thisPoints0);
294  pointField nbrPoints(nbrPoints0);
295 
296  // Create patches for all the points
297 
298  // Source patch at initial location
299  const primitivePatch thisPatch0
300  (
301  SubList<face>(localFaces(), size()),
302  thisPoints0
303  );
304  // Source patch that gets moved
305  primitivePatch thisPatch
306  (
307  SubList<face>(localFaces(), size()),
308  thisPoints
309  );
310  // Target patch at initial location
311  const primitivePatch nbrPatch0
312  (
313  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
314  nbrPoints0
315  );
316  // Target patch that gets moved
317  primitivePatch nbrPatch
318  (
319  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
320  nbrPoints
321  );
322 
323  // Construct a new AMI interpolation between the initial patch locations
324  AMIPtr_.reset
325  (
327  (
328  thisPatch0,
329  nbrPatch0,
330  surfPtr(),
332  false,
334  AMILowWeightCorrection_,
335  AMIReverse_
336  )
337  );
338 
339  // Number of geometry replications
340  label iter(0);
341  label nTransformsOld(nTransforms_);
342 
343  if (ownStr.valid())
344  {
345  writeOBJ(thisPatch0, ownStr());
346  }
347  if (neiStr.valid())
348  {
349  writeOBJ(nbrPatch0, neiStr());
350  }
351 
352 
353  // Weight sum averages
354  scalar srcSum(gAverage(AMIPtr_->srcWeightsSum()));
355  scalar tgtSum(gAverage(AMIPtr_->tgtWeightsSum()));
356 
357  // Direction (or rather side of AMI : this or nbr patch) of
358  // geometry replication
359  bool direction = nTransforms_ >= 0;
360 
361  // Increase in the source weight sum for the last iteration in the
362  // opposite direction. If the current increase is less than this, the
363  // direction (= side of AMI to transform) is reversed.
364  // We switch the side to replicate instead of reversing the transform
365  // since at the coupledPolyPatch level there is no
366  // 'reverseTransformPosition' functionality.
367  scalar srcSumDiff = 0;
368 
369  if (debug)
370  {
371  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:" << name()
372  << " srcSum:" << srcSum
373  << " tgtSum:" << tgtSum
374  << " direction:" << direction
375  << endl;
376  }
377 
378  // Loop, replicating the geometry
379  while
380  (
381  (iter < maxIter_)
382  && (
383  (1 - srcSum > matchTolerance())
384  || (1 - tgtSum > matchTolerance())
385  )
386  )
387  {
388  if (direction)
389  {
390  periodicPatch.transformPosition(thisPoints);
391 
392  if (debug)
393  {
394  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:"
395  << name()
396  << " moving this side from:"
397  << gAverage(thisPatch.points())
398  << " to:" << gAverage(thisPoints) << endl;
399  }
400 
401  thisPatch.movePoints(thisPoints);
402 
403  if (debug)
404  {
405  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:"
406  << name()
407  << " appending weights with untransformed slave side"
408  << endl;
409  }
410 
411  AMIPtr_->append(thisPatch, nbrPatch0);
412 
413  if (ownStr.valid())
414  {
415  writeOBJ(thisPatch, ownStr());
416  }
417  }
418  else
419  {
420  periodicPatch.transformPosition(nbrPoints);
421 
422  if (debug)
423  {
424  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:"
425  << name()
426  << " moving neighbour side from:"
427  << gAverage(nbrPatch.points())
428  << " to:" << gAverage(nbrPoints) << endl;
429  }
430 
431  nbrPatch.movePoints(nbrPoints);
432 
433  AMIPtr_->append(thisPatch0, nbrPatch);
434 
435  if (neiStr.valid())
436  {
437  writeOBJ(nbrPatch, neiStr());
438  }
439  }
440 
441  const scalar srcSumNew = gAverage(AMIPtr_->srcWeightsSum());
442  const scalar srcSumDiffNew = srcSumNew - srcSum;
443 
444  if (srcSumDiffNew < srcSumDiff || srcSumDiffNew < SMALL)
445  {
446  direction = !direction;
447 
448  srcSumDiff = srcSumDiffNew;
449  }
450 
451  srcSum = srcSumNew;
452  tgtSum = gAverage(AMIPtr_->tgtWeightsSum());
453 
454  nTransforms_ += direction ? +1 : -1;
455 
456  ++ iter;
457 
458  if (debug)
459  {
460  Info<< "cyclicPeriodicAMIPolyPatch::resetAMI : patch:" << name()
461  << " iteration:" << iter
462  << " srcSum:" << srcSum
463  << " tgtSum:" << tgtSum
464  << " direction:" << direction
465  << endl;
466  }
467  }
468 
469 
470  // Close debug streams
471  if (ownStr.valid())
472  {
473  ownStr.clear();
474  }
475  if (neiStr.valid())
476  {
477  neiStr.clear();
478  }
479 
480 
481 
482  // Average the number of transformstions
483  nTransforms_ = (nTransforms_ + nTransformsOld)/2;
484 
485  // Check that the match is complete
486  if (iter == maxIter_)
487  {
488  // The matching algorithm has exited without getting the
489  // srcSum and tgtSum above 1. This can happen because
490  // - of an incorrect setup
491  // - or because of non-exact meshes and truncation errors
492  // (transformation, accumulation of cutting errors)
493  // so for now this situation is flagged as a SeriousError instead of
494  // a FatalError since the default matchTolerance is quite strict
495  // (0.001) and can get triggered far into the simulation.
497  << "Patches " << name() << " and " << neighbPatch().name()
498  << " do not couple to within a tolerance of "
499  << matchTolerance()
500  << " when transformed according to the periodic patch "
501  << periodicPatch.name() << "." << nl
502  << "The current sum of weights are for owner " << name()
503  << " : " << srcSum << " and for neighbour "
504  << neighbPatch().name() << " : " << tgtSum << nl
505  << "This is only acceptable during post-processing"
506  << "; not during running. Improve your mesh or increase"
507  << " the 'matchTolerance' setting in the patch specification."
508  << endl;
509  }
510 
511  // Check that both patches have replicated an integer number of times
512  if
513  (
514  mag(srcSum - floor(srcSum + 0.5)) > srcSum*matchTolerance()
515  || mag(tgtSum - floor(tgtSum + 0.5)) > tgtSum*matchTolerance()
516  )
517  {
518  // This condition is currently enforced until there is more
519  // experience with the matching algorithm and numerics.
520  // This check means that e.g. different numbers of stator and
521  // rotor partitions are not allowed.
522  // Again see the section above about tolerances.
524  << "Patches " << name() << " and " << neighbPatch().name()
525  << " do not overlap an integer number of times when transformed"
526  << " according to the periodic patch "
527  << periodicPatch.name() << "." << nl
528  << "The current matchTolerance : " << matchTolerance()
529  << ", sum of owner weights : " << srcSum
530  << ", sum of neighbour weights : " << tgtSum
531  << "." << nl
532  << "This is only acceptable during post-processing"
533  << "; not during running. Improve your mesh or increase"
534  << " the 'matchTolerance' setting in the patch specification."
535  << endl;
536  }
537 
538  // Normalise the weights. Disable printing since weights are
539  // still areas.
540  AMIPtr_->normaliseWeights(true, false);
541 
542  // Print some statistics
543  const label nFace = returnReduce(size(), sumOp<label>());
544 
545  if (nFace)
546  {
547  scalarField srcWghtSum(size(), 0);
548  scalarField tgtWghtSum(size(), 0);
549  forAll(*this, faceI)
550  {
551  srcWghtSum[faceI] = sum(AMIPtr_->srcWeights()[faceI]);
552  tgtWghtSum[faceI] = sum(AMIPtr_->tgtWeights()[faceI]);
553  }
554 
555  Info<< indent
556  << "AMI: Patch " << name()
557  << " sum(weights) min/max/average = "
558  << gMin(srcWghtSum) << ", "
559  << gMax(srcWghtSum) << ", "
560  << gAverage(srcWghtSum) << endl;
561  Info<< indent
562  << "AMI: Patch " << neighbPatch().name()
563  << " sum(weights) min/max/average = "
564  << gMin(tgtWghtSum) << ", "
565  << gMax(tgtWghtSum) << ", "
566  << gAverage(tgtWghtSum) << endl;
567  }
568  }
569 }
570 
571 
572 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
573 
575 (
576  const word& name,
577  const label size,
578  const label start,
579  const label index,
580  const polyBoundaryMesh& bm,
581  const word& patchType,
583 )
584 :
585  cyclicAMIPolyPatch(name, size, start, index, bm, patchType, transform),
586  periodicPatchName_(word::null),
587  periodicPatchID_(-1),
588  nTransforms_(0),
589  nSectors_(0),
590  maxIter_(36)
591 {}
592 
593 
595 (
596  const word& name,
597  const dictionary& dict,
598  const label index,
599  const polyBoundaryMesh& bm,
600  const word& patchType
601 )
602 :
603  cyclicAMIPolyPatch(name, dict, index, bm, patchType),
604  periodicPatchName_(dict.lookup("periodicPatch")),
605  periodicPatchID_(-1),
606  nTransforms_(dict.lookupOrDefault<label>("nTransforms", 0)),
607  nSectors_(dict.lookupOrDefault<label>("nSectors", 0)),
608  maxIter_(dict.lookupOrDefault<label>("maxIter", 36))
609 {}
610 
611 
613 (
614  const cyclicPeriodicAMIPolyPatch& pp,
615  const polyBoundaryMesh& bm
616 )
617 :
618  cyclicAMIPolyPatch(pp, bm),
619  periodicPatchName_(pp.periodicPatchName_),
620  periodicPatchID_(-1),
621  nTransforms_(pp.nTransforms_),
622  nSectors_(pp.nSectors_),
623  maxIter_(pp.maxIter_)
624 {}
625 
626 
628 (
629  const cyclicPeriodicAMIPolyPatch& pp,
630  const polyBoundaryMesh& bm,
631  const label index,
632  const label newSize,
633  const label newStart,
634  const word& nbrPatchName
635 )
636 :
637  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
638  periodicPatchName_(pp.periodicPatchName_),
639  periodicPatchID_(-1),
640  nTransforms_(pp.nTransforms_),
641  nSectors_(pp.nSectors_),
642  maxIter_(pp.maxIter_)
643 {}
644 
645 
647 (
648  const cyclicPeriodicAMIPolyPatch& pp,
649  const polyBoundaryMesh& bm,
650  const label index,
651  const labelUList& mapAddressing,
652  const label newStart
653 )
654 :
655  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
656  periodicPatchName_(pp.periodicPatchName_),
657  periodicPatchID_(-1),
658  nTransforms_(pp.nTransforms_),
659  nSectors_(pp.nSectors_),
660  maxIter_(pp.maxIter_)
661 {}
662 
663 
664 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
665 
667 {}
668 
669 
670 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
671 
673 {
674  if (periodicPatchName_ == word::null)
675  {
676  periodicPatchID_ = -1;
677 
678  return periodicPatchID_;
679  }
680 
681  if (periodicPatchID_ == -1)
682  {
683  periodicPatchID_ = this->boundaryMesh().findPatchID(periodicPatchName_);
684 
685  if (periodicPatchID_ == -1)
686  {
688  << "Illegal periodicPatch name " << periodicPatchName_
689  << nl << "Valid patch names are "
690  << this->boundaryMesh().names()
691  << exit(FatalError);
692  }
693 
694  // Check that it is a coupled patch
695  refCast<const coupledPolyPatch>
696  (
697  this->boundaryMesh()[periodicPatchID_]
698  );
699  }
700 
701  return periodicPatchID_;
702 }
703 
704 
706 {
708 
709  os.writeKeyword("periodicPatch") << periodicPatchName_
710  << token::END_STATEMENT << nl;
711 
712  if (nTransforms_ != 0)
713  {
714  os.writeKeyword("nTransforms") << nTransforms_ <<
716  }
717 
718  if (nSectors_ != 0)
719  {
720  os.writeKeyword("nSectors") << nSectors_ <<
722  }
723 
724  if (maxIter_ != 36)
725  {
726  os.writeKeyword("maxIter") << maxIter_ << token::END_STATEMENT << nl;
727  }
728 }
729 
730 
731 // ************************************************************************* //
Foam::AMIInterpolation::interpolationMethod
interpolationMethod
Enumeration specifying interpolation method.
Definition: AMIInterpolation.H:86
Foam::boundaryMesh::mesh
const bMesh & mesh() const
Definition: boundaryMesh.H:202
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::coupledPolyPatch::separation
virtual const vectorField & separation() const
If the planes are separated the separation vector.
Definition: coupledPolyPatch.H:282
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:738
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::OBJstream
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
writeOBJ
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
Definition: Test-tetTetOverlap.C:38
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:571
Foam::OBJstream::write
virtual Ostream & write(const char)
Write character.
Definition: OBJstream.C:82
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:205
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::coupledPolyPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: coupledPolyPatch.H:294
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::coupledPolyPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: coupledPolyPatch.H:300
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::coupledPolyPatch
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
Definition: coupledPolyPatch.H:51
Foam::cyclicPeriodicAMIPolyPatch::syncTransforms
void syncTransforms() const
Synchronise the periodic transformations.
Definition: cyclicPeriodicAMIPolyPatch.C:54
Foam::cyclicPeriodicAMIPolyPatch::maxIter_
const label maxIter_
Maximum number of attempts to match the AMI geometry.
Definition: cyclicPeriodicAMIPolyPatch.H:71
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::cyclicPeriodicAMIPolyPatch::~cyclicPeriodicAMIPolyPatch
virtual ~cyclicPeriodicAMIPolyPatch()
Destructor.
Definition: cyclicPeriodicAMIPolyPatch.C:666
Foam::cyclicPeriodicAMIPolyPatch::periodicPatchName_
word periodicPatchName_
Periodic patch name.
Definition: cyclicPeriodicAMIPolyPatch.H:59
Foam::AMIMethod
Base class for Arbitrary Mesh Interface (AMI) methods.
Definition: AMIMethod.H:55
Foam::cyclicPeriodicAMIPolyPatch::writeOBJ
void writeOBJ(const primitivePatch &p, OBJstream &str) const
Debug: write obj files of patch (collected on master)
Definition: cyclicPeriodicAMIPolyPatch.C:201
Foam::transform
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
Foam::coupledPolyPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: coupledPolyPatch.H:288
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::Time::rootPath
const fileName & rootPath() const
Return root path.
Definition: Time.H:269
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1040
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:302
Foam::orOp
Definition: ops.H:178
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
SeriousErrorInFunction
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Definition: messageStream.H:237
Foam::PrimitivePatch::movePoints
virtual void movePoints(const Field< PointType > &)
Correct patch after moving points.
Definition: PrimitivePatchTemplate.C:187
cyclicPeriodicAMIPolyPatch.H
Foam::cyclicPeriodicAMIPolyPatch::resetAMI
virtual void resetAMI(const AMIPatchToPatchInterpolation::interpolationMethod &AMIMethod=AMIPatchToPatchInterpolation::imFaceAreaWeight) const
Reset the AMI interpolator.
Definition: cyclicPeriodicAMIPolyPatch.C:228
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::cyclicPeriodicAMIPolyPatch::nSectors_
const label nSectors_
Number of sectors in a rotationally periodic geometry (optional)
Definition: cyclicPeriodicAMIPolyPatch.H:68
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::cyclicPeriodicAMIPolyPatch::periodicPatchID
virtual label periodicPatchID() const
Periodic patch ID.
Definition: cyclicPeriodicAMIPolyPatch.C:672
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
Foam::sumOp
Definition: ops.H:162
Foam::cyclicPeriodicAMIPolyPatch::cyclicPeriodicAMIPolyPatch
cyclicPeriodicAMIPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform=UNKNOWN)
Construct from (base couped patch) components.
Definition: cyclicPeriodicAMIPolyPatch.C:575
Foam::cyclicPeriodicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicPeriodicAMIPolyPatch.C:705
Foam::coupledPolyPatch::transformPosition
virtual void transformPosition(pointField &) const =0
Transform a patch-based position from other side to this side.
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List< vectorField >
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:77
Foam::cyclicPeriodicAMIPolyPatch::nTransforms_
label nTransforms_
Current number of transformations (+ve forward, -ve backward)
Definition: cyclicPeriodicAMIPolyPatch.H:65
Foam::autoPtr::valid
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set).
Definition: autoPtrI.H:83
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::AMIInterpolation::imPartialFaceAreaWeight
@ imPartialFaceAreaWeight
Definition: AMIInterpolation.H:91
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::autoPtr::clear
void clear()
Delete object (if the pointer is valid) and set pointer to NULL.
Definition: autoPtrI.H:126
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:59
Foam::boundaryMesh::findPatchID
label findPatchID(const polyPatchList &, const word &) const
Get index of polypatch by name.
Definition: boundaryMesh.C:254
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::coupledPolyPatch::collocated
virtual const boolList & collocated() const
Are faces collocated. Either size 0,1 or length of patch.
Definition: coupledPolyPatch.H:306
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::TimePaths::globalCaseName
const fileName & globalCaseName() const
Return global case name.
Definition: TimePaths.H:102
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:57
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::cyclicPeriodicAMIPolyPatch
Cyclic patch for periodic Arbitrary Mesh Interface (AMI)
Definition: cyclicPeriodicAMIPolyPatch.H:50
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:43
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
OBJstream.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::faceAreaIntersect::tmMesh
@ tmMesh
Definition: faceAreaIntersect.H:62
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:51