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-2021 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 (returnReduce((size() && !periodicPatch.size()), orOp<bool>()))
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  bool isParallel =
123  (
124  periodicPatch.size()
125  && periodicPatch.parallel()
126  );
127  reduce(isParallel, orOp<bool>());
128 
129  if (isParallel)
130  {
131  // Sync a list of separation vectors
132  List<vectorField> sep(Pstream::nProcs());
133  sep[Pstream::myProcNo()] = periodicPatch.separation();
134  Pstream::gatherList(sep);
136 
137  List<boolList> coll(Pstream::nProcs());
138  coll[Pstream::myProcNo()] = periodicPatch.collocated();
139  Pstream::gatherList(coll);
140  Pstream::scatterList(coll);
141 
142  // If locally we have zero faces pick the first one that has a
143  // separation vector
144  if (!periodicPatch.size())
145  {
146  forAll(sep, procI)
147  {
148  if (sep[procI].size())
149  {
150  const_cast<vectorField&>
151  (
152  periodicPatch.separation()
153  ) = sep[procI];
154  const_cast<boolList&>
155  (
156  periodicPatch.collocated()
157  ) = coll[procI];
158 
159  break;
160  }
161  }
162  }
163  }
164  else
165  {
166  // Sync a list of forward and reverse transforms
167  List<tensorField> forwardT(Pstream::nProcs());
168  forwardT[Pstream::myProcNo()] = periodicPatch.forwardT();
171 
172  List<tensorField> reverseT(Pstream::nProcs());
173  reverseT[Pstream::myProcNo()] = periodicPatch.reverseT();
176 
177  // If locally we have zero faces pick the first one that has a
178  // transformation vector
179  if (!periodicPatch.size())
180  {
181  forAll(forwardT, procI)
182  {
183  if (forwardT[procI].size())
184  {
185  const_cast<tensorField&>
186  (
187  periodicPatch.forwardT()
188  ) = forwardT[procI];
189  const_cast<tensorField&>
190  (
191  periodicPatch.reverseT()
192  ) = reverseT[procI];
193 
194  break;
195  }
196  }
197  }
198  }
199  }
200  }
201 }
202 
203 
204 void Foam::cyclicPeriodicAMIPolyPatch::writeOBJ
205 (
206  const primitivePatch& p,
207  OBJstream& str
208 ) const
209 {
210  // Collect faces and points
212  faceList allFaces;
213  labelList pointMergeMap;
215  (
216  -1.0, // do not merge points
217  p,
218  allPoints,
219  allFaces,
220  pointMergeMap
221  );
222 
223  if (Pstream::master())
224  {
225  // Write base geometry
226  str.write(allFaces, allPoints);
227  }
228 }
229 
230 
231 void Foam::cyclicPeriodicAMIPolyPatch::resetAMI() const
232 {
233  if (owner())
234  {
235  // Get the periodic patch
236  const coupledPolyPatch& periodicPatch
237  (
238  refCast<const coupledPolyPatch>
239  (
240  boundaryMesh()[periodicPatchID()]
241  )
242  );
243 
244  // Synchronise the transforms
245  syncTransforms();
246 
247  // Create copies of both patches' points, transformed to the owner
248  pointField thisPoints0(localPoints());
249  pointField nbrPoints0(neighbPatch().localPoints());
250  transformPosition(nbrPoints0);
251 
252  // Reset the stored number of periodic transformations to a lower
253  // absolute value if possible
254  if (nSectors_ > 0)
255  {
256  if (nTransforms_ > nSectors_/2)
257  {
258  nTransforms_ -= nSectors_;
259  }
260  else if (nTransforms_ < - nSectors_/2)
261  {
262  nTransforms_ += nSectors_;
263  }
264  }
265 
266  // Apply the stored number of periodic transforms
267  for (label i = 0; i < nTransforms_; ++ i)
268  {
269  periodicPatch.transformPosition(thisPoints0);
270  }
271  for (label i = 0; i > nTransforms_; -- i)
272  {
273  periodicPatch.transformPosition(nbrPoints0);
274  }
275 
276  autoPtr<OBJstream> ownStr;
277  autoPtr<OBJstream> neiStr;
278  if (debug)
279  {
280  const Time& runTime = boundaryMesh().mesh().time();
281 
282  fileName dir(runTime.globalPath());
283  fileName postfix("_" + runTime.timeName()+"_expanded.obj");
284 
285  ownStr.reset(new OBJstream(dir/name() + postfix));
286  neiStr.reset(new OBJstream(dir/neighbPatch().name() + postfix));
287 
289  << "patch:" << name()
290  << " writing accumulated AMI to " << ownStr().name()
291  << " and " << neiStr().name() << endl;
292  }
293 
294  // Create another copy
295  pointField thisPoints(thisPoints0);
296  pointField nbrPoints(nbrPoints0);
297 
298  // Create patches for all the points
299 
300  // Source patch at initial location
301  const primitivePatch thisPatch0
302  (
303  SubList<face>(localFaces(), size()),
304  thisPoints0
305  );
306  // Source patch that gets moved
307  primitivePatch thisPatch
308  (
309  SubList<face>(localFaces(), size()),
310  thisPoints
311  );
312  // Target patch at initial location
313  const primitivePatch nbrPatch0
314  (
315  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
316  nbrPoints0
317  );
318  // Target patch that gets moved
319  primitivePatch nbrPatch
320  (
321  SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
322  nbrPoints
323  );
324 
325  // Construct a new AMI interpolation between the initial patch locations
326  AMIPtr_->setRequireMatch(false);
327  AMIPtr_->calculate(thisPatch0, nbrPatch0, surfPtr());
328 
329  // Number of geometry replications
330  label iter(0);
331  label nTransformsOld(nTransforms_);
332 
333  if (ownStr)
334  {
335  writeOBJ(thisPatch0, *ownStr);
336  }
337  if (neiStr)
338  {
339  writeOBJ(nbrPatch0, *neiStr);
340  }
341 
342 
343  // Weight sum averages
344  scalar srcSum(gAverage(AMIPtr_->srcWeightsSum()));
345  scalar tgtSum(gAverage(AMIPtr_->tgtWeightsSum()));
346  // Direction (or rather side of AMI : this or nbr patch) of
347  // geometry replication
348  bool direction = nTransforms_ >= 0;
349 
350  // Increase in the source weight sum for the last iteration in the
351  // opposite direction. If the current increase is less than this, the
352  // direction (= side of AMI to transform) is reversed.
353  // We switch the side to replicate instead of reversing the transform
354  // since at the coupledPolyPatch level there is no
355  // 'reverseTransformPosition' functionality.
356  scalar srcSumDiff = 0;
357 
359  << "patch:" << name()
360  << " srcSum:" << srcSum
361  << " tgtSum:" << tgtSum
362  << " direction:" << direction
363  << endl;
364 
365  // Loop, replicating the geometry
366  while
367  (
368  (iter < maxIter_)
369  && (
370  (1 - srcSum > matchTolerance())
371  || (1 - tgtSum > matchTolerance())
372  )
373  )
374  {
375  if (direction)
376  {
377  periodicPatch.transformPosition(thisPoints);
378 
380  << "patch:" << name()
381  << " moving this side from:"
382  << gAverage(thisPatch.points())
383  << " to:" << gAverage(thisPoints) << endl;
384 
385  thisPatch.movePoints(thisPoints);
386 
388  << "patch:" << name()
389  << " appending weights with untransformed slave side"
390  << endl;
391 
392  AMIPtr_->append(thisPatch, nbrPatch0);
393 
394  if (ownStr)
395  {
396  writeOBJ(thisPatch, *ownStr);
397  }
398  }
399  else
400  {
401  periodicPatch.transformPosition(nbrPoints);
402 
404  << "patch:" << name()
405  << " moving neighbour side from:"
406  << gAverage(nbrPatch.points())
407  << " to:" << gAverage(nbrPoints) << endl;
408 
409  nbrPatch.movePoints(nbrPoints);
410 
411  AMIPtr_->append(thisPatch0, nbrPatch);
412 
413  if (neiStr)
414  {
415  writeOBJ(nbrPatch, *neiStr);
416  }
417  }
418 
419  const scalar srcSumNew = gAverage(AMIPtr_->srcWeightsSum());
420  const scalar srcSumDiffNew = srcSumNew - srcSum;
421 
422  if (srcSumDiffNew < srcSumDiff || srcSumDiffNew < SMALL)
423  {
424  direction = !direction;
425 
426  srcSumDiff = srcSumDiffNew;
427  }
428 
429  srcSum = srcSumNew;
430  tgtSum = gAverage(AMIPtr_->tgtWeightsSum());
431 
432  nTransforms_ += direction ? +1 : -1;
433 
434  ++iter;
435 
437  << "patch:" << name()
438  << " iteration:" << iter
439  << " srcSum:" << srcSum
440  << " tgtSum:" << tgtSum
441  << " direction:" << direction
442  << endl;
443  }
444 
445 
446  // Close debug streams
447  ownStr.reset(nullptr);
448  neiStr.reset(nullptr);
449 
450 
451  // Average the number of transformations
452  nTransforms_ = (nTransforms_ + nTransformsOld)/2;
453 
454  // Check that the match is complete
455  if (iter == maxIter_)
456  {
457  // The matching algorithm has exited without getting the
458  // srcSum and tgtSum above 1. This can happen because
459  // - of an incorrect setup
460  // - or because of non-exact meshes and truncation errors
461  // (transformation, accumulation of cutting errors)
462  // so for now this situation is flagged as a SeriousError instead of
463  // a FatalError since the default matchTolerance is quite strict
464  // (0.001) and can get triggered far into the simulation.
466  << "Patches " << name() << " and " << neighbPatch().name()
467  << " do not couple to within a tolerance of "
468  << matchTolerance()
469  << " when transformed according to the periodic patch "
470  << periodicPatch.name() << "." << nl
471  << "The current sum of weights are for owner " << name()
472  << " : " << srcSum << " and for neighbour "
473  << neighbPatch().name() << " : " << tgtSum << nl
474  << "This is only acceptable during post-processing"
475  << "; not during running. Improve your mesh or increase"
476  << " the 'matchTolerance' setting in the patch specification."
477  << endl;
478  }
479 
480  // Check that both patches have replicated an integer number of times
481  if
482  (
483  mag(srcSum - floor(srcSum + 0.5)) > srcSum*matchTolerance()
484  || mag(tgtSum - floor(tgtSum + 0.5)) > tgtSum*matchTolerance()
485  )
486  {
487  // This condition is currently enforced until there is more
488  // experience with the matching algorithm and numerics.
489  // This check means that e.g. different numbers of stator and
490  // rotor partitions are not allowed.
491  // Again see the section above about tolerances.
493  << "Patches " << name() << " and " << neighbPatch().name()
494  << " do not overlap an integer number of times when transformed"
495  << " according to the periodic patch "
496  << periodicPatch.name() << "." << nl
497  << "The current matchTolerance : " << matchTolerance()
498  << ", sum of owner weights : " << srcSum
499  << ", sum of neighbour weights : " << tgtSum
500  << "." << nl
501  << "This is only acceptable during post-processing"
502  << "; not during running. Improve your mesh or increase"
503  << " the 'matchTolerance' setting in the patch specification."
504  << endl;
505  }
506 
507  // Print some statistics
508  const label nFace = returnReduce(size(), sumOp<label>());
509 
510  if (nFace)
511  {
512  scalarField srcWghtSum(size(), Zero);
513  forAll(srcWghtSum, faceI)
514  {
515  srcWghtSum[faceI] = sum(AMIPtr_->srcWeights()[faceI]);
516  }
517  scalarField tgtWghtSum(neighbPatch().size(), Zero);
518  forAll(tgtWghtSum, faceI)
519  {
520  tgtWghtSum[faceI] = sum(AMIPtr_->tgtWeights()[faceI]);
521  }
522 
523  Info<< indent
524  << "AMI: Patch " << name()
525  << " sum(weights)"
526  << " min:" << gMin(srcWghtSum)
527  << " max:" << gMax(srcWghtSum)
528  << " average:" << gAverage(srcWghtSum) << nl;
529  Info<< indent
530  << "AMI: Patch " << neighbPatch().name()
531  << " sum(weights)"
532  << " min:" << gMin(tgtWghtSum)
533  << " max:" << gMax(tgtWghtSum)
534  << " average:" << gAverage(tgtWghtSum) << nl;
535  }
536  }
537 }
538 
539 
540 // * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
541 
543 (
544  const word& name,
545  const label size,
546  const label start,
547  const label index,
548  const polyBoundaryMesh& bm,
549  const word& patchType,
551 )
552 :
554  (
555  name,
556  size,
557  start,
558  index,
559  bm,
560  patchType,
561  transform,
562  faceAreaWeightAMI::typeName
563  ),
564  nTransforms_(0),
565  nSectors_(0),
566  maxIter_(36)
567 {
568  AMIPtr_->setRequireMatch(false);
569 }
570 
571 
573 (
574  const word& name,
575  const dictionary& dict,
576  const label index,
577  const polyBoundaryMesh& bm,
578  const word& patchType
579 )
580 :
582  (
583  name,
584  dict,
585  index,
586  bm,
587  patchType,
588  faceAreaWeightAMI::typeName
589  ),
590  nTransforms_(dict.getOrDefault<label>("nTransforms", 0)),
591  nSectors_(dict.getOrDefault<label>("nSectors", 0)),
592  maxIter_(dict.getOrDefault<label>("maxIter", 36))
593 {
594  AMIPtr_->setRequireMatch(false);
595 }
596 
597 
599 (
600  const cyclicPeriodicAMIPolyPatch& pp,
601  const polyBoundaryMesh& bm
602 )
603 :
604  cyclicAMIPolyPatch(pp, bm),
605  nTransforms_(pp.nTransforms_),
606  nSectors_(pp.nSectors_),
607  maxIter_(pp.maxIter_)
608 {
609  AMIPtr_->setRequireMatch(false);
610 }
611 
612 
614 (
615  const cyclicPeriodicAMIPolyPatch& pp,
616  const polyBoundaryMesh& bm,
617  const label index,
618  const label newSize,
619  const label newStart,
620  const word& nbrPatchName
621 )
622 :
623  cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
624  nTransforms_(pp.nTransforms_),
625  nSectors_(pp.nSectors_),
626  maxIter_(pp.maxIter_)
627 {
628  AMIPtr_->setRequireMatch(false);
629 }
630 
631 
633 (
634  const cyclicPeriodicAMIPolyPatch& pp,
635  const polyBoundaryMesh& bm,
636  const label index,
637  const labelUList& mapAddressing,
638  const label newStart
639 )
640 :
641  cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
642  nTransforms_(pp.nTransforms_),
643  nSectors_(pp.nSectors_),
644  maxIter_(pp.maxIter_)
645 {
646  AMIPtr_->setRequireMatch(false);
647 }
648 
649 
650 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
651 
653 {}
654 
655 
656 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
657 
659 {
661 
662  os.writeEntryIfDifferent<label>("nTransforms", 0, nTransforms_);
663  os.writeEntryIfDifferent<label>("nSectors", 0, nSectors_);
664  os.writeEntryIfDifferent<label>("maxIter", 36, maxIter_);
665 }
666 
667 
668 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::cyclicAMIPolyPatch::owner
virtual bool owner() const
Does this side own the patch?
Definition: cyclicAMIPolyPatch.C:964
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PatchTools.H
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::cyclicAMIPolyPatch::periodicPatchName_
word periodicPatchName_
Periodic patch name.
Definition: cyclicAMIPolyPatch.H:130
Foam::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:215
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:38
Foam::coupledPolyPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: coupledPolyPatch.H:301
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::coupledPolyPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: coupledPolyPatch.H:307
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::cyclicPeriodicAMIPolyPatch::~cyclicPeriodicAMIPolyPatch
virtual ~cyclicPeriodicAMIPolyPatch()
Destructor.
Definition: cyclicPeriodicAMIPolyPatch.C:652
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::cyclicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicAMIPolyPatch.C:1225
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
SeriousErrorInFunction
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
Definition: messageStream.H:306
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.C:40
cyclicPeriodicAMIPolyPatch.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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:52
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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 coupled patch) components.
Definition: cyclicPeriodicAMIPolyPatch.C:543
Foam::cyclicPeriodicAMIPolyPatch::write
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
Definition: cyclicPeriodicAMIPolyPatch.C:658
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::UList< label >
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::TimePaths::globalPath
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:80
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::primitivePatch
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Definition: primitivePatch.H:51
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::coupledPolyPatch::transformType
transformType
Definition: coupledPolyPatch.H:60
Foam::cyclicPeriodicAMIPolyPatch
Cyclic patch for periodic Arbitrary Mesh Interface (AMI)
Definition: cyclicPeriodicAMIPolyPatch.H:52
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::cyclicAMIPolyPatch::periodicPatchID
label periodicPatchID() const
Periodic patch ID (or -1)
Definition: cyclicAMIPolyPatch.C:938
Foam::coupledPolyPatch::coupledPolyPatch
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.
Definition: coupledPolyPatch.C:476
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
OBJstream.H
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:68