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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 
31 // For debugging
32 #include "OBJstream.H"
33 #include "PatchTools.H"
34 #include "Time.H"
35 // Note: cannot use vtkSurfaceWriter here - circular linkage
36 // but foamVtkSurfaceWriter (vtk::surfaceWriter) would be okay.
37 //
38 // #include "foamVtkSurfaceWriter.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44  defineTypeNameAndDebug(cyclicPeriodicAMIPolyPatch, 0);
45 
46  addToRunTimeSelectionTable(polyPatch, cyclicPeriodicAMIPolyPatch, word);
48  (
49  polyPatch,
50  cyclicPeriodicAMIPolyPatch,
51  dictionary
52  );
53 }
54 
55 
56 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
57 
58 void Foam::cyclicPeriodicAMIPolyPatch::syncTransforms() const
59 {
60  if (owner())
61  {
62  // At this point we guarantee that the transformations have been
63  // updated. There is one particular case, where if the periodic patch
64  // locally has zero faces then its transformation will not be set. We
65  // need to synchronise the transforms over the zero-sized patches as
66  // well.
67  //
68  // We can't put the logic into cyclicPolyPatch as
69  // processorCyclicPolyPatch uses cyclicPolyPatch functionality.
70  // Synchronisation will fail because processor-type patches do not exist
71  // on all processors.
72  //
73  // The use in cyclicPeriodicAMI is special; we use the patch even
74  // though we have no faces in it. Ideally, in future, the transformation
75  // logic will be abstracted out, and we won't need a periodic patch
76  // here. Until then, we can just fix the behaviour of the zero-sized
77  // coupled patches here
78 
79  // Get the periodic patch
80  const coupledPolyPatch& periodicPatch
81  (
82  refCast<const coupledPolyPatch>
83  (
85  )
86  );
87 
88  // If there are any zero-sized periodic patches
89  if (returnReduceOr(size() && !periodicPatch.size()))
90  {
91  if (periodicPatch.separation().size() > 1)
92  {
94  << "Periodic patch " << periodicPatchName_
95  << " has non-uniform separation vector "
96  << periodicPatch.separation()
97  << "This is not allowed inside " << type()
98  << " patch " << name()
99  << exit(FatalError);
100  }
101 
102  if (periodicPatch.forwardT().size() > 1)
103  {
105  << "Periodic patch " << periodicPatchName_
106  << " has non-uniform transformation tensor "
107  << periodicPatch.forwardT()
108  << "This is not allowed inside " << type()
109  << " patch " << name()
110  << exit(FatalError);
111  }
112 
113  // Note that zero-sized patches will have zero-sized fields for the
114  // separation vector, forward and reverse transforms. These need
115  // replacing with the transformations from other processors.
116 
117  // Parallel in this context refers to a parallel transformation,
118  // rather than a rotational transformation.
119 
120  // Note that a cyclic with zero faces is considered parallel so
121  // explicitly check for that.
122 
123  if
124  (
126  (
127  periodicPatch.size() && periodicPatch.parallel()
128  )
129  )
130  {
131  // Sync a list of separation vectors
132  List<vectorField> sep(Pstream::nProcs());
133  sep[Pstream::myProcNo()] = periodicPatch.separation();
135 
136  List<boolList> coll(Pstream::nProcs());
137  coll[Pstream::myProcNo()] = periodicPatch.collocated();
139 
140  // If locally we have zero faces pick the first one that has a
141  // separation vector
142  if (!periodicPatch.size())
143  {
144  forAll(sep, procI)
145  {
146  if (sep[procI].size())
147  {
148  const_cast<vectorField&>
149  (
150  periodicPatch.separation()
151  ) = sep[procI];
152  const_cast<boolList&>
153  (
154  periodicPatch.collocated()
155  ) = coll[procI];
156 
157  break;
158  }
159  }
160  }
161  }
162  else
163  {
164  // Sync a list of forward and reverse transforms
165  List<tensorField> forwardT(Pstream::nProcs());
166  forwardT[Pstream::myProcNo()] = periodicPatch.forwardT();
168 
169  List<tensorField> reverseT(Pstream::nProcs());
170  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 
200 void Foam::cyclicPeriodicAMIPolyPatch::writeOBJ
201 (
202  const primitivePatch& p,
203  OBJstream& str
204 ) const
205 {
206  // Collect faces and points
208  faceList allFaces;
210  (
211  -1.0, // do not merge points
212  p,
213  allPoints,
214  allFaces
215  );
216 
217  if (Pstream::master())
218  {
219  // Write base geometry
220  str.write(allFaces, allPoints);
221  }
222 }
223 
224 
225 void Foam::cyclicPeriodicAMIPolyPatch::resetAMI() const
226 {
227  if (owner())
228  {
229  // Get the periodic patch
230  const coupledPolyPatch& periodicPatch
231  (
232  refCast<const coupledPolyPatch>
233  (
234  boundaryMesh()[periodicPatchID()]
235  )
236  );
237 
238  // Synchronise the transforms
239  syncTransforms();
240 
241  // Create copies of both patches' points, transformed to the owner
242  pointField thisPoints0(localPoints());
243  pointField nbrPoints0(neighbPatch().localPoints());
244  transformPosition(nbrPoints0);
245 
246  // Reset the stored number of periodic transformations to a lower
247  // absolute value if possible
248  if (nSectors_ > 0)
249  {
250  if (nTransforms_ > nSectors_/2)
251  {
252  nTransforms_ -= nSectors_;
253  }
254  else if (nTransforms_ < - nSectors_/2)
255  {
256  nTransforms_ += nSectors_;
257  }
258  }
259 
260  // Apply the stored number of periodic transforms
261  for (label i = 0; i < nTransforms_; ++ i)
262  {
263  periodicPatch.transformPosition(thisPoints0);
264  }
265  for (label i = 0; i > nTransforms_; -- i)
266  {
267  periodicPatch.transformPosition(nbrPoints0);
268  }
269 
270  autoPtr<OBJstream> ownStr;
271  autoPtr<OBJstream> neiStr;
272  if (debug)
273  {
274  const Time& runTime = boundaryMesh().mesh().time();
275 
276  fileName dir(runTime.globalPath());
277  fileName postfix("_" + runTime.timeName()+"_expanded.obj");
278 
279  ownStr.reset(new OBJstream(dir/name() + postfix));
280  neiStr.reset(new OBJstream(dir/neighbPatch().name() + postfix));
281 
283  << "patch:" << name()
284  << " writing accumulated AMI to " << ownStr().name()
285  << " and " << neiStr().name() << endl;
286  }
287 
288  // Create another copy
289  pointField thisPoints(thisPoints0);
290  pointField nbrPoints(nbrPoints0);
291 
292  // Create patches for all the points
293 
294  // Source patch at initial location
295  const primitivePatch thisPatch0
296  (
297  SubList<face>(localFaces(), size()),
298  thisPoints0
299  );
300  // Source patch that gets moved
301  primitivePatch thisPatch
302  (
303  SubList<face>(localFaces(), size()),
304  thisPoints
305  );
306  // Target patch at initial location
307  const primitivePatch nbrPatch0
308  (
309  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
310  nbrPoints0
311  );
312  // Target patch that gets moved
313  primitivePatch nbrPatch
314  (
315  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
316  nbrPoints
317  );
318 
319  // Construct a new AMI interpolation between the initial patch locations
320  AMIPtr_->setRequireMatch(false);
321  AMIPtr_->calculate(thisPatch0, nbrPatch0, surfPtr());
322 
323  // Number of geometry replications
324  label iter(0);
325  label nTransformsOld(nTransforms_);
326 
327  if (ownStr)
328  {
329  writeOBJ(thisPatch0, *ownStr);
330  }
331  if (neiStr)
332  {
333  writeOBJ(nbrPatch0, *neiStr);
334  }
335 
336 
337  // Weight sum averages
338  scalar srcSum(gAverage(AMIPtr_->srcWeightsSum()));
339  scalar tgtSum(gAverage(AMIPtr_->tgtWeightsSum()));
340  // Direction (or rather side of AMI : this or nbr patch) of
341  // geometry replication
342  bool direction = nTransforms_ >= 0;
343 
344  // Increase in the source weight sum for the last iteration in the
345  // opposite direction. If the current increase is less than this, the
346  // direction (= side of AMI to transform) is reversed.
347  // We switch the side to replicate instead of reversing the transform
348  // since at the coupledPolyPatch level there is no
349  // 'reverseTransformPosition' functionality.
350  scalar srcSumDiff = 0;
351 
353  << "patch:" << name()
354  << " srcSum:" << srcSum
355  << " tgtSum:" << tgtSum
356  << " direction:" << direction
357  << endl;
358 
359  // Loop, replicating the geometry
360  while
361  (
362  (iter < maxIter_)
363  && (
364  (1 - srcSum > matchTolerance())
365  || (1 - tgtSum > matchTolerance())
366  )
367  )
368  {
369  if (direction)
370  {
371  periodicPatch.transformPosition(thisPoints);
372 
374  << "patch:" << name()
375  << " moving this side from:"
376  << gAverage(thisPatch.points())
377  << " to:" << gAverage(thisPoints) << endl;
378 
379  thisPatch.movePoints(thisPoints);
380 
382  << "patch:" << name()
383  << " appending weights with untransformed slave side"
384  << endl;
385 
386  AMIPtr_->append(thisPatch, nbrPatch0);
387 
388  if (ownStr)
389  {
390  writeOBJ(thisPatch, *ownStr);
391  }
392  }
393  else
394  {
395  periodicPatch.transformPosition(nbrPoints);
396 
398  << "patch:" << name()
399  << " moving neighbour side from:"
400  << gAverage(nbrPatch.points())
401  << " to:" << gAverage(nbrPoints) << endl;
402 
403  nbrPatch.movePoints(nbrPoints);
404 
405  AMIPtr_->append(thisPatch0, nbrPatch);
406 
407  if (neiStr)
408  {
409  writeOBJ(nbrPatch, *neiStr);
410  }
411  }
412 
413  const scalar srcSumNew = gAverage(AMIPtr_->srcWeightsSum());
414  const scalar srcSumDiffNew = srcSumNew - srcSum;
415 
416  if (srcSumDiffNew < srcSumDiff || srcSumDiffNew < SMALL)
417  {
418  direction = !direction;
419 
420  srcSumDiff = srcSumDiffNew;
421  }
422 
423  srcSum = srcSumNew;
424  tgtSum = gAverage(AMIPtr_->tgtWeightsSum());
425 
426  nTransforms_ += direction ? +1 : -1;
427 
428  ++iter;
429 
431  << "patch:" << name()
432  << " iteration:" << iter
433  << " srcSum:" << srcSum
434  << " tgtSum:" << tgtSum
435  << " direction:" << direction
436  << endl;
437  }
438 
439 
440  // Close debug streams
441  ownStr.reset(nullptr);
442  neiStr.reset(nullptr);
443 
444 
445  // Average the number of transformations
446  nTransforms_ = (nTransforms_ + nTransformsOld)/2;
447 
448  // Check that the match is complete
449  if (iter == maxIter_)
450  {
451  // The matching algorithm has exited without getting the
452  // srcSum and tgtSum above 1. This can happen because
453  // - of an incorrect setup
454  // - or because of non-exact meshes and truncation errors
455  // (transformation, accumulation of cutting errors)
456  // so for now this situation is flagged as a SeriousError instead of
457  // a FatalError since the default matchTolerance is quite strict
458  // (0.001) and can get triggered far into the simulation.
460  << "Patches " << name() << " and " << neighbPatch().name()
461  << " do not couple to within a tolerance of "
462  << matchTolerance()
463  << " when transformed according to the periodic patch "
464  << periodicPatch.name() << "." << nl
465  << "The current sum of weights are for owner " << name()
466  << " : " << srcSum << " and for neighbour "
467  << neighbPatch().name() << " : " << tgtSum << nl
468  << "This is only acceptable during post-processing"
469  << "; not during running. Improve your mesh or increase"
470  << " the 'matchTolerance' setting in the patch specification."
471  << endl;
472  }
473 
474  // Check that both patches have replicated an integer number of times
475  if
476  (
477  mag(srcSum - floor(srcSum + 0.5)) > srcSum*matchTolerance()
478  || mag(tgtSum - floor(tgtSum + 0.5)) > tgtSum*matchTolerance()
479  )
480  {
481  // This condition is currently enforced until there is more
482  // experience with the matching algorithm and numerics.
483  // This check means that e.g. different numbers of stator and
484  // rotor partitions are not allowed.
485  // Again see the section above about tolerances.
487  << "Patches " << name() << " and " << neighbPatch().name()
488  << " do not overlap an integer number of times when transformed"
489  << " according to the periodic patch "
490  << periodicPatch.name() << "." << nl
491  << "The current matchTolerance : " << matchTolerance()
492  << ", sum of owner weights : " << srcSum
493  << ", sum of neighbour weights : " << tgtSum
494  << "." << nl
495  << "This is only acceptable during post-processing"
496  << "; not during running. Improve your mesh or increase"
497  << " the 'matchTolerance' setting in the patch specification."
498  << endl;
499  }
500 
501  // Print some statistics
502  if (returnReduceOr(size()))
503  {
504  scalarField srcWghtSum(size(), Zero);
505  forAll(srcWghtSum, faceI)
506  {
507  srcWghtSum[faceI] = sum(AMIPtr_->srcWeights()[faceI]);
508  }
509  scalarField tgtWghtSum(neighbPatch().size(), Zero);
510  forAll(tgtWghtSum, faceI)
511  {
512  tgtWghtSum[faceI] = sum(AMIPtr_->tgtWeights()[faceI]);
513  }
514 
515  Info<< indent
516  << "AMI: Patch " << name()
517  << " sum(weights)"
518  << " min:" << gMin(srcWghtSum)
519  << " max:" << gMax(srcWghtSum)
520  << " average:" << gAverage(srcWghtSum) << nl;
521  Info<< indent
522  << "AMI: Patch " << neighbPatch().name()
523  << " sum(weights)"
524  << " min:" << gMin(tgtWghtSum)
525  << " max:" << gMax(tgtWghtSum)
526  << " average:" << gAverage(tgtWghtSum) << nl;
527  }
528  }
529 }
530 
531 
532 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
533 
535 (
536  const word& name,
537  const label size,
538  const label start,
539  const label index,
540  const polyBoundaryMesh& bm,
541  const word& patchType,
542  const transformType transform
543 )
544 :
545  cyclicAMIPolyPatch
546  (
547  name,
548  size,
549  start,
550  index,
551  bm,
552  patchType,
553  transform,
554  faceAreaWeightAMI::typeName
555  ),
556  nTransforms_(0),
557  nSectors_(0),
558  maxIter_(36)
559 {
560  AMIPtr_->setRequireMatch(false);
561 }
562 
563 
565 (
566  const word& name,
567  const dictionary& dict,
568  const label index,
569  const polyBoundaryMesh& bm,
570  const word& patchType
571 )
572 :
574  (
575  name,
576  dict,
577  index,
578  bm,
579  patchType,
580  faceAreaWeightAMI::typeName
581  ),
582  nTransforms_(dict.getOrDefault<label>("nTransforms", 0)),
583  nSectors_(dict.getOrDefault<label>("nSectors", 0)),
584  maxIter_(dict.getOrDefault<label>("maxIter", 36))
585 {
586  AMIPtr_->setRequireMatch(false);
587 }
588 
589 
591 (
593  const polyBoundaryMesh& bm
594 )
595 :
596  cyclicAMIPolyPatch(pp, bm),
597  nTransforms_(pp.nTransforms_),
598  nSectors_(pp.nSectors_),
599  maxIter_(pp.maxIter_)
600 {
601  AMIPtr_->setRequireMatch(false);
602 }
603 
604 
606 (
608  const polyBoundaryMesh& bm,
609  const label index,
610  const label newSize,
611  const label newStart,
612  const word& nbrPatchName
613 )
614 :
615  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
616  nTransforms_(pp.nTransforms_),
617  nSectors_(pp.nSectors_),
618  maxIter_(pp.maxIter_)
619 {
620  AMIPtr_->setRequireMatch(false);
621 }
622 
623 
625 (
627  const polyBoundaryMesh& bm,
628  const label index,
629  const labelUList& mapAddressing,
630  const label newStart
631 )
632 :
633  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
634  nTransforms_(pp.nTransforms_),
635  nSectors_(pp.nSectors_),
636  maxIter_(pp.maxIter_)
637 {
638  AMIPtr_->setRequireMatch(false);
639 }
640 
641 
642 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
645 {}
646 
647 
648 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
649 
651 {
653 
654  os.writeEntryIfDifferent<label>("nTransforms", 0, nTransforms_);
655  os.writeEntryIfDifferent<label>("nSectors", 0, nSectors_);
656  os.writeEntryIfDifferent<label>("maxIter", 36, maxIter_);
657 }
658 
659 
660 // ************************************************************************* //
dictionary dict
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
uint8_t direction
Definition: direction.H:46
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
virtual const fileName & name() const
The name of the stream.
Definition: IOstream.C:33
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
coupledPolyPatch(const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm, const word &patchType, const transformType transform)
Construct from components.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
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 coupled patch) components.
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
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
Macros for easy insertion into run-time selection tables.
virtual const tensorField & forwardT() const
Return face transformation tensor.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &pp, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, globalIndex &pointAddr, globalIndex &faceAddr, labelList &pointMergeMap=const_cast< labelList &>(labelList::null()), const bool useLocal=false)
Gather points and faces onto master and merge into single patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Face area weighted Arbitrary Mesh Interface (AMI) method.
static void allGatherList(List< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses linear/tree communication. ...
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: polyPatch.C:309
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Time & time() const noexcept
Return time registry.
#define DebugInFunction
Report an information message using Foam::Info.
word periodicPatchName_
Periodic patch name.
Cyclic patch for Arbitrary Mesh Interface (AMI)
label periodicPatchID() const
Periodic patch ID (or -1)
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 polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
const word & name() const noexcept
The patch name.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field...
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Type gAverage(const FieldField< Field, Type > &f)
virtual bool owner() const
Does this side own the patch?
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Cyclic patch for periodic Arbitrary Mesh Interface (AMI)
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< bool > boolList
A List of bools.
Definition: List.H:60
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127