Loading...
Searching...
No Matches
displacementLayeredMotionMotionSolver.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
32#include "pointFields.H"
33#include "PointEdgeWave.H"
34#include "syncTools.H"
35#include "interpolationTable.H"
36#include "mapPolyMesh.H"
37#include "pointConstraints.H"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
44
46 (
50 );
51
53 (
56 displacement
57 );
58}
59
60
61// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62
64Foam::displacementLayeredMotionMotionSolver::getCylindrical
65(
66 const label cellZoneI,
67 const dictionary& zoneDict
68)
69{
70 auto iter = cylSystems_.cfind(cellZoneI);
71
72 if (iter.good())
73 {
74 return *(iter.val());
75 }
76
77 cylSystems_.set(cellZoneI, new coordSystem::cylindrical(zoneDict));
78
79 return *cylSystems_[cellZoneI];
80}
81
82
83void Foam::displacementLayeredMotionMotionSolver::calcZoneMask
84(
85 const label cellZoneI,
86 bitSet& isZonePoint,
87 bitSet& isZoneEdge
88) const
89{
90 isZonePoint.resize(mesh().nPoints());
91 isZoneEdge.resize(mesh().nEdges());
92
93 if (cellZoneI == -1)
94 {
95 isZonePoint = true;
96 isZoneEdge = true;
97 return;
98 }
99
100
101 isZonePoint.reset();
102 isZoneEdge.reset();
103
104 const cellZone& cz = mesh().cellZones()[cellZoneI];
105
106 // Mark points, edges inside cellZone
107 for (const label celli : cz)
108 {
109 isZonePoint.set(mesh().cellPoints(celli));
110 isZoneEdge.set(mesh().cellEdges(celli));
111 }
112
114 (
115 mesh(),
116 isZonePoint,
118 0
119 );
120
122 (
123 mesh(),
124 isZoneEdge,
126 0
127 );
128
130 << "On cellZone " << cz.name() << " marked "
131 << returnReduce(isZonePoint.count(), sumOp<label>()) << " points and "
132 << returnReduce(isZoneEdge.count(), sumOp<label>()) << " edges" << nl;
133}
134
135
136// Find distance to starting point
137void Foam::displacementLayeredMotionMotionSolver::walkStructured
138(
139 const label cellZoneI,
140 const bitSet& isZonePoint,
141 const bitSet& isZoneEdge,
142 const labelList& seedPoints,
143 const vectorField& seedData,
145 vectorField& data,
146 labelField& patchPoints
147) const
148{
149 List<pointEdgeStructuredWalk> seedInfo(seedPoints.size());
150
151 forAll(seedPoints, i)
152 {
153 const label pointi = seedPoints[i];
154
155 seedInfo[i] = pointEdgeStructuredWalk
156 (
157 points0()[pointi], // location of data
158 points0()[pointi], // previous location
159 0.0,
160 seedData[i],
161 pointi // pass thru seed point id
162 );
163 }
164
165 // Current info on points
167
168 // Mark points inside cellZone.
169 // Note that we use points0, not mesh.points()
170 // so as not to accumulate errors.
171 for (const label pointi : isZonePoint)
172 {
173 allPointInfo[pointi] = pointEdgeStructuredWalk
174 (
175 points0()[pointi], // location of data
176 point::max, // not valid
177 0.0,
178 Zero // passive data
179 );
180 }
181
182
183 // Current info on edges
184 List<pointEdgeStructuredWalk> allEdgeInfo(mesh().nEdges());
185
186 // Mark edges inside cellZone
187 for (const label edgei : isZoneEdge)
188 {
189 allEdgeInfo[edgei] = pointEdgeStructuredWalk
190 (
191 mesh().edges()[edgei].centre(points0()), // location of data
192 point::max, // not valid
193 0.0,
194 Zero
195 );
196 }
197
198 // Walk
200 (
201 mesh(),
202 seedPoints,
203 seedInfo,
204
205 allPointInfo,
206 allEdgeInfo,
207 mesh().globalData().nTotalPoints() // max iterations
208 );
209
210 // Extract distance and passive data
211 for (const label pointi : isZonePoint)
212 {
213 const auto& pointInfo = allPointInfo[pointi];
214
215 distance[pointi] = pointInfo.dist();
216 data[pointi] = pointInfo.data();
217
218 // Optional information
219 if (patchPoints.size())
220 {
221 patchPoints[pointi] = pointInfo.index();
222 }
223 }
224}
225
226
227// Evaluate faceZone patch
228Foam::tmp<Foam::vectorField>
229Foam::displacementLayeredMotionMotionSolver::faceZoneEvaluate
230(
231 const faceZone& fz,
232 const labelList& meshPoints,
233 const dictionary& dict,
234 const PtrList<pointVectorField>& patchDisp,
235 const label patchi
236) const
237{
238 auto tfld = tmp<vectorField>::New(meshPoints.size());
239 auto& fld = tfld.ref();
240
241 const word type(dict.get<word>("type"));
242
243 if (type == "fixedValue")
244 {
245 fld.assign("value", dict, meshPoints.size());
246 }
247 else if (type == "timeVaryingUniformFixedValue")
248 {
250
251 fld = timeSeries(mesh().time().timeOutputValue());
252 }
253 else if (type == "slip")
254 {
255 if ((patchi % 2) != 1)
256 {
258 << "FaceZone:" << fz.name()
259 << exit(FatalIOError);
260 }
261 // Use field set by previous bc
262 fld = vectorField(patchDisp[patchi - 1], meshPoints);
263 }
264 else if (type == "follow")
265 {
266 // Only on boundary faces - follow boundary conditions
267 fld = vectorField(pointDisplacement_, meshPoints);
268 }
269 else if (type == "uniformFollow")
270 {
271 // Reads name of name of patch. Then get average point displacement on
272 // patch. That becomes the value of fld.
273 const word patchName(dict.get<word>("patch"));
274 const label patchID = mesh().boundaryMesh().findPatchID(patchName);
275 pointField pdf
276 (
277 pointDisplacement_.boundaryField()[patchID].patchInternalField()
278 );
279 fld = gAverage(pdf);
280 }
281 else
282 {
284 << "Unknown faceZonePatch type " << type
285 << " for faceZone " << fz.name() << nl
286 << exit(FatalIOError);
287 }
288
289 return tfld;
290}
291
292
293void Foam::displacementLayeredMotionMotionSolver::cellZoneSolve
294(
295 const label cellZoneI,
296 const dictionary& zoneDict
297)
298{
299 bitSet isZonePoint, isZoneEdge;
300 calcZoneMask(cellZoneI, isZonePoint, isZoneEdge);
301
302 const dictionary& patchesDict = zoneDict.subDict("boundaryField");
303
304 if (patchesDict.size() != 2)
305 {
307 << "Two faceZones (patches) must be specified per cellZone. "
308 << " cellZone:" << cellZoneI
309 << " patches:" << patchesDict.toc()
310 << exit(FatalIOError);
311 }
312
313 PtrList<labelField> patchPoints(patchesDict.size());
314 PtrList<scalarField> patchDist(patchesDict.size());
315 PtrList<pointVectorField> patchDisp(patchesDict.size());
316
317 // Allocate the fields
318 label patchi = 0;
319 for (const entry& dEntry : patchesDict)
320 {
321 const word& faceZoneName = dEntry.keyword();
322 const label zoneI = mesh().faceZones().findZoneID(faceZoneName);
323 if (zoneI == -1)
324 {
326 << "Cannot find faceZone " << faceZoneName
327 << endl << "Valid zones are " << mesh().faceZones().names()
328 << exit(FatalIOError);
329 }
330
331 // Determine the points of the faceZone within the cellZone
332 const faceZone& fz = mesh().faceZones()[zoneI];
333
334 patchPoints.set(patchi, new labelField(mesh().nPoints(), label(-1)));
335 patchDist.set(patchi, new scalarField(mesh().nPoints()));
336 patchDisp.set
337 (
338 patchi,
340 (
341 mesh().newIOobject
342 (
343 mesh().cellZones()[cellZoneI].name() + "_" + fz.name()
344 ),
345 pointDisplacement_ // to inherit the boundary conditions
346 )
347 );
348
349 ++patchi;
350 }
351
352
353
354 // 'correctBoundaryConditions'
355 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
356 // Loops over all the faceZones and walks their boundary values
357
358 // Make sure we can pick up bc values from field
359 pointDisplacement_.correctBoundaryConditions();
360
361 patchi = 0;
362 for (const entry& dEntry : patchesDict)
363 {
364 const word& faceZoneName = dEntry.keyword();
365 const dictionary& faceZoneDict = dEntry.dict();
366
367 // Determine the points of the faceZone within the cellZone
368 const faceZone& fz = mesh().faceZones()[faceZoneName];
369 const labelList& fzMeshPoints = fz().meshPoints();
370 DynamicList<label> meshPoints(fzMeshPoints.size());
371 forAll(fzMeshPoints, i)
372 {
373 if (isZonePoint[fzMeshPoints[i]])
374 {
375 meshPoints.append(fzMeshPoints[i]);
376 }
377 }
378
379 // Get initial value for all the faceZone points
380 tmp<vectorField> tseed = faceZoneEvaluate
381 (
382 fz,
383 meshPoints,
384 faceZoneDict,
385 patchDisp,
386 patchi
387 );
388
389 if (debug)
390 {
391 auto limits = gMinMax(tseed());
392 auto avg = gAverage(tseed());
393
394 Info
395 << "For cellZone:" << cellZoneI
396 << " for faceZone:" << fz.name()
397 << " nPoints:" << tseed().size()
398 << " have patchField:"
399 << " min:" << limits.min()
400 << " max:" << limits.max()
401 << " avg:" << avg
402 << endl;
403 }
404
405
406 // Set distance and transported value
407 walkStructured
408 (
409 cellZoneI,
410 isZonePoint,
411 isZoneEdge,
412
413 meshPoints,
414 tseed,
415 patchDist[patchi],
416 patchDisp[patchi],
417 patchPoints[patchi]
418 );
419
420 // Implement real bc.
421 patchDisp[patchi].correctBoundaryConditions();
422
423 ++patchi;
424 }
425
426
427 // Solve
428 // ~~~~~
429
430 if (debug)
431 {
432 // Normalised distance
433 auto tdistance = pointScalarField::New
434 (
436 (
437 mesh().cellZones()[cellZoneI].name(),
438 "distance"
439 ),
443 );
444 auto& distance = tdistance.ref();
445
446 for (const label pointi : isZonePoint)
447 {
448 const scalar d1 = patchDist[0][pointi];
449 const scalar d2 = patchDist[1][pointi];
450 if (d1 + d2 > SMALL)
451 {
452 distance[pointi] = d1/(d1 + d2);
453 }
454 }
455
456 Info<< "Writing " << pointScalarField::typeName
457 << distance.name() << " to "
458 << mesh().time().timeName() << endl;
459 distance.write();
460 }
461
462
463 const word interpolationScheme(zoneDict.get<word>("interpolationScheme"));
464
465 if (interpolationScheme == "oneSided")
466 {
467 for (const label pointi : isZonePoint)
468 {
469 pointDisplacement_[pointi] = patchDisp[0][pointi];
470 }
471 }
472 else if (interpolationScheme == "linear")
473 {
474 for (const label pointi : isZonePoint)
475 {
476 const scalar d1 = patchDist[0][pointi];
477 const scalar d2 = patchDist[1][pointi];
478 const scalar s = d1/(d1 + d2 + VSMALL);
479
480 const vector& pd1 = patchDisp[0][pointi];
481 const vector& pd2 = patchDisp[1][pointi];
482
483 pointDisplacement_[pointi] = (1 - s)*pd1 + s*pd2;
484 }
485 }
486 else if (interpolationScheme == "cylindrical")
487 {
488 const coordSystem::cylindrical& cs =
489 this->getCylindrical(cellZoneI, zoneDict);
490
491 // May wish to implement alternative distance calculation here
492
493
495 << "cylindrical interpolation not yet available" << nl
496 << "coordinate system " << cs << nl
497 << exit(FatalError);
498 }
499 else
500 {
502 << "Invalid interpolationScheme: " << interpolationScheme
503 << ". Valid schemes: (oneSided linear cylindrical)" << nl
505 }
506}
507
508
509// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
510
511Foam::displacementLayeredMotionMotionSolver::
512displacementLayeredMotionMotionSolver
513(
514 const polyMesh& mesh,
516)
517:
519{}
520
521
522Foam::displacementLayeredMotionMotionSolver::
523displacementLayeredMotionMotionSolver
524(
525 const polyMesh& mesh,
526 const IOdictionary& dict,
527 const pointVectorField& pointDisplacement,
528 const pointIOField& points0
529)
530:
532{}
533
534
535// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
536
539{
540 tmp<pointField> tcurPoints
541 (
542 points0() + pointDisplacement_.primitiveField()
543 );
544
545 return tcurPoints;
546}
547
548
550{
551 // The points have moved so before interpolation update the motionSolver
552 movePoints(mesh().points());
553
554 // Apply boundary conditions
555 pointDisplacement_.boundaryFieldRef().updateCoeffs();
556
557 // Solve motion on all regions (=cellZones)
558 for (const entry& dEntry : coeffDict().subDict("regions"))
559 {
560 const word& cellZoneName = dEntry.keyword();
561 const dictionary& regionDict = dEntry.dict();
562
563 const label zoneI = mesh().cellZones().findZoneID(cellZoneName);
564
565 Info<< "solving for zone: " << cellZoneName << endl;
566
567 if (zoneI == -1)
568 {
570 << "Cannot find cellZone " << cellZoneName
571 << endl << "Valid zones are " << mesh().cellZones().names()
572 << exit(FatalIOError);
573 }
574
575 cellZoneSolve(zoneI, regionDict);
576 }
577
578 // Update pointDisplacement for solved values
579 const pointConstraints& pcs =
582}
583
584
586(
587 const mapPolyMesh& mpm
588)
589{
591
592 const vectorField displacement(this->newPoints() - points0_);
593
594 forAll(points0_, pointi)
595 {
596 const label oldPointi = mpm.pointMap()[pointi];
597
598 if (oldPointi >= 0)
599 {
600 const label masterPointi = mpm.reversePointMap()[oldPointi];
601
602 if ((masterPointi != pointi))
603 {
604 // newly inserted point in this cellZone
605
606 // need to set point0 so that it represents the position that
607 // it would have had if it had existed for all time
608 points0_[pointi] -= displacement[pointi];
609 }
610 }
611 }
612}
613
614
615// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
static const char *const typeName
Typename for Field.
Definition Field.H:93
static tmp< GeometricField< scalar, pointPatchField, pointMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=pointPatchField< scalar >::calculatedType())
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Wave propagation of information through grid. Every iteration information goes through one layer of e...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
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
static const Form max
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition ZoneMesh.C:757
wordList names() const
A list of the zone names.
Definition ZoneMesh.C:463
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
A subset of mesh cells.
Definition cellZone.H:61
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
Mesh motion solver for an (multi-block) extruded fvMesh. Gets given the structure of the mesh blocks ...
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
virtual void updateMesh(const mapPolyMesh &mpm)
Update topology.
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
pointVectorField pointDisplacement_
Point motion field.
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const word & name() const
Return reference to name.
Definition fvMesh.H:387
An interpolation/look-up table of scalar vs <Type> values. The reference scalar values must be monoto...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const noexcept
Reverse point map.
const labelList & pointMap() const noexcept
Old point map.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
virtual tmp< pointField > newPoints()
Provide new points for motion. Solves for motion.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Application of (multi-)patch point constraints.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Determines length of string of edges walked to point.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
pointIOField points0_
Starting points.
virtual void movePoints(const pointField &)
Update local data for geometry changes.
virtual void updateMesh(const mapPolyMesh &)
Update local data for topology changes.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition polyMesh.H:679
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
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
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define DebugInfo
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
scalar distance(const vector &p1, const vector &p2)
Definition curveTools.C:12
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
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.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
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
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))