Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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"
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
45
48 (
52 );
53}
54
55
56// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
57
58void 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 (
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
200void 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
226{
227 if (owner())
228 {
229 // Get the periodic patch
230 const coupledPolyPatch& periodicPatch
231 (
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 const fileName dir(runTime.globalPath());
277 const 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 "
285 << ownStr().name().relative(dir)
286 << " and " << neiStr().name().relative(dir) << endl;
287 }
288
289 // Create another copy
290 pointField thisPoints(thisPoints0);
291 pointField nbrPoints(nbrPoints0);
292
293 // Create patches for all the points
294
295 // Source patch at initial location
296 const primitivePatch thisPatch0
297 (
298 SubList<face>(localFaces(), size()),
299 thisPoints0
300 );
301 // Source patch that gets moved
302 primitivePatch thisPatch
303 (
304 SubList<face>(localFaces(), size()),
305 thisPoints
306 );
307 // Target patch at initial location
308 const primitivePatch nbrPatch0
309 (
310 SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
311 nbrPoints0
312 );
313 // Target patch that gets moved
314 primitivePatch nbrPatch
315 (
316 SubList<face>(neighbPatch().localFaces(), neighbPatch().size()),
317 nbrPoints
318 );
319
320 // Construct a new AMI interpolation between the initial patch locations
321 AMIPtr_->setRequireMatch(false);
322 AMIPtr_->calculate(thisPatch0, nbrPatch0, surfPtr());
323
324 // Number of geometry replications
325 label iter(0);
326 label nTransformsOld(nTransforms_);
327
328 if (ownStr)
329 {
330 writeOBJ(thisPatch0, *ownStr);
331 }
332 if (neiStr)
333 {
334 writeOBJ(nbrPatch0, *neiStr);
335 }
336
337
338 // Weight sum averages
339 scalar srcSum(gAverage(AMIPtr_->srcWeightsSum()));
340 scalar tgtSum(gAverage(AMIPtr_->tgtWeightsSum()));
341 // Direction (or rather side of AMI : this or nbr patch) of
342 // geometry replication
343 bool direction = nTransforms_ >= 0;
344
345 // Increase in the source weight sum for the last iteration in the
346 // opposite direction. If the current increase is less than this, the
347 // direction (= side of AMI to transform) is reversed.
348 // We switch the side to replicate instead of reversing the transform
349 // since at the coupledPolyPatch level there is no
350 // 'reverseTransformPosition' functionality.
351 scalar srcSumDiff = 0;
352
354 << "patch:" << name()
355 << " srcSum:" << srcSum
356 << " tgtSum:" << tgtSum
357 << " direction:" << direction
358 << endl;
359
360 // Loop, replicating the geometry
361 while
362 (
363 (iter < maxIter_)
364 && (
365 (1 - srcSum > matchTolerance())
366 || (1 - tgtSum > matchTolerance())
367 )
368 )
369 {
370 if (direction)
371 {
372 periodicPatch.transformPosition(thisPoints);
373
375 << "patch:" << name()
376 << " moving this side from:"
377 << gAverage(thisPatch.points())
378 << " to:" << gAverage(thisPoints) << endl;
379
380 thisPatch.movePoints(thisPoints);
381
383 << "patch:" << name()
384 << " appending weights with untransformed slave side"
385 << endl;
386
387 AMIPtr_->append(thisPatch, nbrPatch0);
388
389 if (ownStr)
390 {
391 writeOBJ(thisPatch, *ownStr);
392 }
393 }
394 else
395 {
396 periodicPatch.transformPosition(nbrPoints);
397
399 << "patch:" << name()
400 << " moving neighbour side from:"
401 << gAverage(nbrPatch.points())
402 << " to:" << gAverage(nbrPoints) << endl;
403
404 nbrPatch.movePoints(nbrPoints);
405
406 AMIPtr_->append(thisPatch0, nbrPatch);
407
408 if (neiStr)
409 {
410 writeOBJ(nbrPatch, *neiStr);
411 }
412 }
413
414 const scalar srcSumNew = gAverage(AMIPtr_->srcWeightsSum());
415 const scalar srcSumDiffNew = srcSumNew - srcSum;
416
417 if (srcSumDiffNew < srcSumDiff || srcSumDiffNew < SMALL)
418 {
420
421 srcSumDiff = srcSumDiffNew;
422 }
423
424 srcSum = srcSumNew;
425 tgtSum = gAverage(AMIPtr_->tgtWeightsSum());
426
427 nTransforms_ += direction ? +1 : -1;
428
429 ++iter;
430
432 << "patch:" << name()
433 << " iteration:" << iter
434 << " srcSum:" << srcSum
435 << " tgtSum:" << tgtSum
436 << " direction:" << direction
437 << endl;
438 }
439
440
441 // Close debug streams
442 ownStr.reset(nullptr);
443 neiStr.reset(nullptr);
444
445
446 // Average the number of transformations
447 nTransforms_ = (nTransforms_ + nTransformsOld)/2;
448
449 // Check that the match is complete
450 if (iter == maxIter_)
451 {
452 // The matching algorithm has exited without getting the
453 // srcSum and tgtSum above 1. This can happen because
454 // - of an incorrect setup
455 // - or because of non-exact meshes and truncation errors
456 // (transformation, accumulation of cutting errors)
457 // so for now this situation is flagged as a SeriousError instead of
458 // a FatalError since the default matchTolerance is quite strict
459 // (0.001) and can get triggered far into the simulation.
461 << "Patches " << name() << " and " << neighbPatch().name()
462 << " do not couple to within a tolerance of "
463 << matchTolerance()
464 << " when transformed according to the periodic patch "
465 << periodicPatch.name() << "." << nl
466 << "The current sum of weights are for owner " << name()
467 << " : " << srcSum << " and for neighbour "
468 << neighbPatch().name() << " : " << tgtSum << nl
469 << "This is only acceptable during post-processing"
470 << "; not during running. Improve your mesh or increase"
471 << " the 'matchTolerance' setting in the patch specification."
472 << endl;
473 }
474
475 // Check that both patches have replicated an integer number of times
476 if
477 (
478 mag(srcSum - floor(srcSum + 0.5)) > srcSum*matchTolerance()
479 || mag(tgtSum - floor(tgtSum + 0.5)) > tgtSum*matchTolerance()
480 )
481 {
482 // This condition is currently enforced until there is more
483 // experience with the matching algorithm and numerics.
484 // This check means that e.g. different numbers of stator and
485 // rotor partitions are not allowed.
486 // Again see the section above about tolerances.
488 << "Patches " << name() << " and " << neighbPatch().name()
489 << " do not overlap an integer number of times when transformed"
490 << " according to the periodic patch "
491 << periodicPatch.name() << "." << nl
492 << "The current matchTolerance : " << matchTolerance()
493 << ", sum of owner weights : " << srcSum
494 << ", sum of neighbour weights : " << tgtSum
495 << "." << nl
496 << "This is only acceptable during post-processing"
497 << "; not during running. Improve your mesh or increase"
498 << " the 'matchTolerance' setting in the patch specification."
499 << endl;
500 }
501
502 // Print some statistics
503 if (returnReduceOr(size()))
504 {
505 scalarField srcWghtSum(size(), Zero);
506 forAll(srcWghtSum, faceI)
507 {
508 srcWghtSum[faceI] = sum(AMIPtr_->srcWeights()[faceI]);
509 }
510 scalarField tgtWghtSum(neighbPatch().size(), Zero);
511 forAll(tgtWghtSum, faceI)
512 {
513 tgtWghtSum[faceI] = sum(AMIPtr_->tgtWeights()[faceI]);
514 }
515
516 {
517 auto limits = gMinMax(srcWghtSum);
518 auto avg = gAverage(srcWghtSum);
519
520 Info<< indent
521 << "AMI: Patch " << name()
522 << " sum(weights)"
523 << " min:" << limits.min()
524 << " max:" << limits.max()
525 << " average:" << avg << nl;
526 }
527
528 {
529 auto limits = gMinMax(tgtWghtSum);
530 auto avg = gAverage(tgtWghtSum);
531
532 Info<< indent
533 << "AMI: Patch " << neighbPatch().name()
534 << " sum(weights)"
535 << " min:" << limits.min()
536 << " max:" << limits.max()
537 << " average:" << avg << nl;
538 }
540 }
541}
542
543
544// * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
545
547(
548 const word& name,
549 const label size,
550 const label start,
551 const label index,
552 const polyBoundaryMesh& bm,
553 const word& patchType,
554 const transformType transform
555)
556:
558 (
559 name,
560 size,
561 start,
562 index,
563 bm,
564 patchType,
565 transform,
567 ),
568 nTransforms_(0),
569 nSectors_(0),
570 maxIter_(36)
571{
572 AMIPtr_->setRequireMatch(false);
573}
574
575
577(
578 const word& name,
579 const dictionary& dict,
580 const label index,
581 const polyBoundaryMesh& bm,
582 const word& patchType
583)
584:
586 (
587 name,
588 dict,
589 index,
590 bm,
591 patchType,
593 ),
594 nTransforms_(dict.getOrDefault<label>("nTransforms", 0)),
595 nSectors_(dict.getOrDefault<label>("nSectors", 0)),
596 maxIter_(dict.getOrDefault<label>("maxIter", 36))
597{
598 AMIPtr_->setRequireMatch(false);
599}
600
601
603(
605 const polyBoundaryMesh& bm
606)
607:
609 nTransforms_(pp.nTransforms_),
610 nSectors_(pp.nSectors_),
611 maxIter_(pp.maxIter_)
612{
613 AMIPtr_->setRequireMatch(false);
614}
615
616
618(
620 const polyBoundaryMesh& bm,
621 const label index,
622 const label newSize,
623 const label newStart,
624 const word& nbrPatchName
625)
626:
627 cyclicAMIPolyPatch(pp, bm, index, newSize, newStart, nbrPatchName),
628 nTransforms_(pp.nTransforms_),
629 nSectors_(pp.nSectors_),
630 maxIter_(pp.maxIter_)
631{
632 AMIPtr_->setRequireMatch(false);
633}
634
635
637(
639 const polyBoundaryMesh& bm,
640 const label index,
641 const labelUList& mapAddressing,
642 const label newStart
643)
644:
645 cyclicAMIPolyPatch(pp, bm, index, mapAddressing, newStart),
646 nTransforms_(pp.nTransforms_),
647 nSectors_(pp.nSectors_),
648 maxIter_(pp.maxIter_)
650 AMIPtr_->setRequireMatch(false);
651}
652
653
654// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
657{}
658
659
660// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
661
663{
665
666 os.writeEntryIfDifferent<label>("nTransforms", 0, nTransforms_);
667 os.writeEntryIfDifferent<label>("nSectors", 0, nSectors_);
668 os.writeEntryIfDifferent<label>("maxIter", 36, maxIter_);
669}
670
671
672// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
fileName globalPath() const
The global path for the case = rootPath/globalCaseName.
Definition TimePathsI.H:108
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
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:1697
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
const bMesh & mesh() const
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
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.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrary Mesh Interface (AMI).
virtual void resetAMI() const
Reset the AMI interpolator, use current patch points.
virtual bool owner() const
Does this side own the patch?
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
word periodicPatchName_
Periodic patch name.
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< AMIPatchToPatchInterpolation > AMIPtr_
AMI interpolation class.
label periodicPatchID() const
Periodic patch ID (or -1).
Cyclic patch for periodic Arbitrary Mesh Interface (AMI).
virtual void write(Ostream &) const
Write the polyPatch data as a dictionary.
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Face area weighted Arbitrary Mesh Interface (AMI) method.
A class for handling file names.
Definition fileName.H:75
label index() const noexcept
The index of this patch in the boundaryMesh.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition polyPatch.C:295
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
auto limits
Definition setRDeltaT.H:186
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
#define SeriousErrorInFunction
Report an error message using Foam::SeriousError.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Namespace for handling debugging switches.
Definition debug.C:45
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
List< bool > boolList
A List of bools.
Definition List.H:60
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299