Loading...
Searching...
No Matches
volumeExprDriver.H
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) 2019-2025 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
26Class
27 Foam::expressions::volumeExpr::parseDriver
28
29Description
30 Driver for volume, surface, point field expressions
31
32 Additional Properties
33 \table
34 Property | Description | Required | Default
35 dimensions | Dimensions for the expression result | no |
36 \endtable
37
38 In addition to the standard mathematical functions, operations and
39 logical and relational operations, the volume expressions support the
40 following driver-specific functions:
41
42 Functions
43 \table
44 Function | Description | Number of arguments |
45 vol | The cell volumes | 0 |
46 pos | The cell centres | 0 |
47 pts | The cell points | 0 |
48 area | The face area magnitudes | 0 |
49 fpos | The face centres | 0 |
50 weightAverage| Volume or area weighted average | 1 |
51 weightSum | Volume or area weighted sum | 1 |
52 face | The face areaNormal vectors | 0 |
53 face | A surface-field face value | 1 |
54 point | A point-field point value | 1 |
55 cellToFace | Interpolate cell values onto faces | 1 |
56 cellToPoint | Interpolate cell values onto points | 1 |
57 pointToCell | Interpolate point values onto cells | 1 |
58 reconstruct | Reconstruct cell vector from surface scalar | 1 |
59 rand | Random field | 0/1 |
60 \endtable
61
62 Selections
63 \table
64 Function| Description | Number of arguments |
65 cset | Logical vol field corresponding to cellSet | 1 |
66 fset | Logical surf field corresponding to faceSet | 1 |
67 pset | Logical point field corresponding to pointSet | 1 |
68 czone | Logical vol field corresponding to cellZone | 1 |
69 fzone | Logical surf field corresponding to faceZone | 1 |
70 pzone | Logical point field corresponding to pointZone| 1 |
71 \endtable
72
73Note
74 Use namespace debug switch \c volumeExpr for scanner (2), parser (4)
75 or dictionary controls as per Foam::expressions::exprDriver.
76
77SourceFiles
78 volumeExprDriver.C
79
80\*---------------------------------------------------------------------------*/
81
82#ifndef expressions_volumeExprDriver_H
83#define expressions_volumeExprDriver_H
84
85#include "volumeExprFwd.H"
86#include "fvExprDriver.H"
88#include "volFields.H"
89#include "surfaceFields.H"
90#include "pointFields.H"
92
93// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
94
95namespace Foam
96{
97namespace expressions
98{
99namespace volumeExpr
100{
101
102/*---------------------------------------------------------------------------*\
103 Class parseDriver Declaration
104\*---------------------------------------------------------------------------*/
105
106class parseDriver
107:
108 public parsing::genericRagelLemonDriver,
109 public expressions::fvExprDriver
110{
111protected:
112
113 // Protected Data
114
115 //- The referenced mesh
116 const fvMesh& mesh_;
117
118 //- The results (volume, surface, point)
119 autoPtr<regIOobject> resultField_;
120
121 //- The result type-name.
122 // Normally volScalarField, surfaceVectorField etc,
123 // but Scalar is modified for logical as volScalarField etc
124 word resultType_;
125
126 //- A logical (bool-like) field (but actually a scalar)
127 bool isLogical_;
128
129 //- Requested use of dimensions
130 bool hasDimensions_;
131
132 //- A volume/surface/point field
134
135 //- The result dimensions
136 dimensionSet resultDimensions_;
137
138
139 // Protected Member Functions
140
141 //- Deep-copy the internalField as a result.
142 // Uses the isLogical() and isPointData() values to handle
143 // additional bookkeeping.
144 // For isLogical(), renames the resultType_ from '*Scalar*'
145 // to '*Logical*' (eg, volLogicalField)
146 template<class Type>
147 void setInternalFieldResult(const Field<Type>& fld);
148
149 //- Cell selections (as logical)
150 tmp<volScalarField> field_cellSelection
151 (
152 const word& name,
153 enum topoSetSource::sourceType setType
154 ) const;
155
156 //- Face selections (as logical)
157 tmp<surfaceScalarField> field_faceSelection
158 (
159 const word& name,
160 enum topoSetSource::sourceType setType
161 ) const;
162
163 //- Point selections (as logical)
164 tmp<pointScalarField> field_pointSelection
165 (
166 const word& name,
167 enum topoSetSource::sourceType setType
168 ) const;
169
170
171public:
172
173 ClassName("volumeExpr::driver");
174
175 // Generated Methods
176
177 // No copy copy construct
178 parseDriver(const parseDriver&) = delete;
179
180 // No copy assignment
181 void operator=(const parseDriver&) = delete;
182
183
184 // Constructors
185
186 //- Construct for specified mesh, with dictionary information
187 explicit parseDriver
188 (
189 const fvMesh& mesh,
191 );
192
193 //- Construct for specified mesh with copy of driver context
195 (
196 const fvMesh& mesh,
197 const parseDriver& driver,
198 const dictionary& dict
199 );
200
201 //- Construct with meshName for the given mesh
202 parseDriver(const word& meshName, const fvMesh& mesh);
203
204 //- Construct with patchName and region specified in dictionary
205 parseDriver(const dictionary& dict, const fvMesh& mesh);
206
207
208 // Not generally clonable
209
210
211 //- Destructor
212 virtual ~parseDriver() = default;
213
214
215 // Public Member Functions
216
217 //- The mesh we are attached to
218 virtual const fvMesh& mesh() const
219 {
220 return mesh_;
221 }
222
223 //- The natural field size for the expression
224 virtual label size() const
225 {
226 return mesh_.nCells();
227 }
228
229 //- The point field size for the expression
230 virtual label pointSize() const
231 {
232 return mesh_.nPoints();
233 }
234
235 //- Field size associated with different geometric field types
236 inline label size(const FieldAssociation geoType) const;
237
238 //- Apply dimensions() to geometric fields
239 inline bool hasDimensions() const noexcept;
240
241 //- The preferred result dimensions (if any)
242 inline const dimensionSet& dimensions() const noexcept;
244
245 //- Clear out local copies of the field
246 void clearField();
247
248
249 // Reading
250
251 //- Read variables, tables etc.
252 // Adds support for "dimensions"
253 virtual bool readDict(const dictionary& dict);
254
255
256 // Evaluation
257
258 //- Perform parsing on (sub) string
260
261 //- Execute the parser.
262 // The return value currently has no meaning.
263 virtual unsigned parse
264 (
265 const std::string& expr,
266 size_t pos = 0,
267 size_t len = std::string::npos
268 );
269
270
271 // Field Information
272
273 //- The result type-name.
274 // Normally volScalarField, surfaceVectorField etc,
275 // but Scalar is modified for logical as volScalarField etc
276 const word& resultType() const noexcept
277 {
278 return resultType_;
279 }
280
281 //- The geometric field association
283 {
284 return fieldGeoType_;
285 }
286
287 //- A logical (bool-like) field. Actually stored as a scalar.
288 bool isLogical() const noexcept
289 {
290 return isLogical_;
291 }
292
293 //- A volume field
294 bool isVolumeData() const noexcept
295 {
297 }
298
299 //- A surface field
300 bool isFaceData() const noexcept
301 {
303 }
304
305 //- A point field
306 bool isPointData() const noexcept
307 {
309 }
310
311 //- Test if stored result pointer is the specified type
312 template<class GeoField>
313 const GeoField* isResultType() const;
314
315 //- Test if stored result pointer is the specified type
316 //- and matches the specified logical type
317 template<class GeoField>
318 const GeoField* isResultType(bool logical, bool dieOnNull=false) const;
319
320 //- A zero-initialized field with the same type as the result field.
323
324 // Set Fields
326 //- Set result (vol field)
327 template<class Type>
328 void setResult(VolumeField<Type>* ptr, bool logical = false);
329
330 //- Set result (surface field)
331 template<class Type>
332 void setResult(SurfaceField<Type>* ptr, bool logical = false);
333
334 //- Set result (point field)
335 template<class Type>
336 void setResult(PointField<Type>* ptr, bool logical = false);
337
338
339 // New Fields
340
341 //- Return a new volume field with the mesh size
342 template<class Type>
344 newVolField(const Type& val = pTraits<Type>::zero) const;
345
346 //- Return a new surface field with the mesh nInternalFaces size
347 template<class Type>
349 newSurfaceField(const Type& val = pTraits<Type>::zero) const;
350
351 //- Return a new point field with the mesh nPoints size
352 template<class Type>
354 newPointField(const Type& val = pTraits<Type>::zero) const;
355
356
357 //- Retrieve field (vol field)
358 template<class Type>
360 getVolField(const word& fldName, bool getOldTime=false);
361
362 //- Retrieve field (surface field)
363 template<class Type>
365 getSurfaceField(const word& fldName, bool getOldTime=false);
367 //- Retrieve field (surface field)
368 template<class Type>
370 getPointField(const word& fldName, bool getOldTime=false);
371
372
373 // Field "shape" conversions
375 //- Interpolate cell to face values
376 template<class Type>
378 cellToFace(const VolumeField<Type>& field) const;
379
380 //- Interpolate cell to point values
381 template<class Type>
383 cellToPoint(const VolumeField<Type>& field) const;
384
385 //- Interpolate point to cell values
386 template<class Type>
388 pointToCell(const PointField<Type>& field) const;
389
391 // Custom Field Functions
392
393 //- The volume-weighted average of a field
394 template<class Type>
395 Type volAverage(VolumeField<Type>& fld) const
396 {
397 return gWeightedAverage(fld.mesh().V(), fld.primitiveField());
398 }
399
400 //- The volume-weighted sum of a field
401 template<class Type>
402 Type volSum(VolumeField<Type>& fld) const
403 {
404 return gWeightedSum(fld.mesh().V(), fld.primitiveField());
405 }
406
407 //- The area-weighted average of a field
408 template<class Type>
410 {
411 return gWeightedAverage
412 (
413 fld.mesh().magSf().primitiveField(),
414 fld.primitiveField()
415 );
416 }
417
418 //- The area-weighted sum of a field
419 template<class Type>
420 Type areaSum(SurfaceField<Type>& fld) const
421 {
422 return gWeightedSum
423 (
424 fld.mesh().magSf().primitiveField(),
425 fld.primitiveField()
426 );
427 }
428
429
430 //- The cell volumes - (swak = vol)
431 tmp<volScalarField> field_cellVolume() const;
432
433 //- The cell centres - (swak = pos)
434 tmp<volVectorField> field_cellCentre() const;
435
436 //- The face area magnitudes [magSf] - (swak = area)
437 tmp<surfaceScalarField> field_faceArea() const;
438
439 //- The face centres - (swak = fpos)
440 tmp<surfaceVectorField> field_faceCentre() const;
441
442 //- The face areas with their vector direction [Sf] - (swak = face)
443 tmp<surfaceVectorField> field_areaNormal() const;
444
445 //- The mesh point locations - (swak = pts)
446 tmp<pointVectorField> field_pointField() const;
447
448
449 //- Cell selection (set)
450 inline tmp<volScalarField> field_cellSet(const word& name) const;
451
452 //- Cell selection (zone)
453 inline tmp<volScalarField> field_cellZone(const word& name) const;
454
455 //- Face selection (set)
456 inline tmp<surfaceScalarField> field_faceSet(const word& name) const;
457
458 //- Face selection (zone)
459 inline tmp<surfaceScalarField> field_faceZone(const word& name) const;
460
461 //- Point selection (set)
462 inline tmp<pointScalarField> field_pointSet(const word& name) const;
464 //- Point selection (zone)
465 inline tmp<pointScalarField> field_pointZone(const word& name) const;
466
467 //- A uniform random field
468 tmp<volScalarField> field_rand(label seed=0, bool gaussian=false) const;
469
470 //- A Gaussian random field
472 {
473 return field_rand(seed, true);
474 }
475};
476
477
478// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480} // End namespace volumeExpr
481} // End namespace expressions
482} // End namespace Foam
483
484
485// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
486
488
489#ifdef NoRepository
491#endif
492
493// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
494
495#endif
496
497// ************************************************************************* //
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))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const dictionary & dict() const noexcept
The dictionary with all input data/specification.
Definition exprDriver.H:455
Base driver for parsing value expressions associated with an fvMesh.
Driver for volume, surface, point field expressions.
virtual label size() const
The natural field size for the expression.
Type volAverage(VolumeField< Type > &fld) const
The volume-weighted average of a field.
const fvMesh & mesh_
The referenced mesh.
const word & resultType() const noexcept
The result type-name.
tmp< volScalarField > field_cellVolume() const
The cell volumes - (swak = vol).
tmp< volScalarField > field_randGaussian(label seed=0) const
A Gaussian random field.
tmp< VolumeField< Type > > newVolField(const Type &val=pTraits< Type >::zero) const
Return a new volume field with the mesh size.
autoPtr< regIOobject > dupZeroField() const
A zero-initialized field with the same type as the result field.
FieldAssociation fieldAssociation() const noexcept
The geometric field association.
expressions::FieldAssociation fieldGeoType_
A volume/surface/point field.
tmp< surfaceScalarField > field_faceArea() const
The face area magnitudes [magSf] - (swak = area).
bool isPointData() const noexcept
A point field.
tmp< VolumeField< Type > > pointToCell(const PointField< Type > &field) const
Interpolate point to cell values.
virtual label pointSize() const
The point field size for the expression.
tmp< PointField< Type > > getPointField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field).
virtual unsigned parse(const std::string &expr, size_t pos=0, size_t len=std::string::npos)
Execute the parser.
tmp< pointVectorField > field_pointField() const
The mesh point locations - (swak = pts).
void operator=(const parseDriver &)=delete
tmp< PointField< Type > > newPointField(const Type &val=pTraits< Type >::zero) const
Return a new point field with the mesh nPoints size.
tmp< surfaceScalarField > field_faceZone(const word &name) const
Face selection (zone).
tmp< pointScalarField > field_pointSet(const word &name) const
Point selection (set).
const GeoField * isResultType(bool logical, bool dieOnNull=false) const
Test if stored result pointer is the specified type and matches the specified logical type.
bool isLogical() const noexcept
A logical (bool-like) field. Actually stored as a scalar.
virtual bool readDict(const dictionary &dict)
Read variables, tables etc.
Type volSum(VolumeField< Type > &fld) const
The volume-weighted sum of a field.
parseDriver(const parseDriver &)=delete
tmp< pointScalarField > field_pointZone(const word &name) const
Point selection (zone).
const GeoField * isResultType() const
Test if stored result pointer is the specified type.
tmp< VolumeField< Type > > getVolField(const word &fldName, bool getOldTime=false)
Retrieve field (vol field).
tmp< surfaceScalarField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
Face selections (as logical).
tmp< SurfaceField< Type > > cellToFace(const VolumeField< Type > &field) const
Interpolate cell to face values.
bool hasDimensions_
Requested use of dimensions.
tmp< surfaceVectorField > field_areaNormal() const
The face areas with their vector direction [Sf] - (swak = face).
bool hasDimensions() const noexcept
Apply dimensions() to geometric fields.
tmp< SurfaceField< Type > > newSurfaceField(const Type &val=pTraits< Type >::zero) const
Return a new surface field with the mesh nInternalFaces size.
tmp< PointField< Type > > cellToPoint(const VolumeField< Type > &field) const
Interpolate cell to point values.
tmp< surfaceScalarField > field_faceSet(const word &name) const
Face selection (set).
tmp< volScalarField > field_cellSet(const word &name) const
Cell selection (set).
tmp< surfaceVectorField > field_faceCentre() const
The face centres - (swak = fpos).
tmp< volScalarField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
Cell selections (as logical).
tmp< volScalarField > field_cellZone(const word &name) const
Cell selection (zone).
tmp< volScalarField > field_rand(label seed=0, bool gaussian=false) const
A uniform random field.
tmp< SurfaceField< Type > > getSurfaceField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field).
void setInternalFieldResult(const Field< Type > &fld)
Deep-copy the internalField as a result.
Type areaSum(SurfaceField< Type > &fld) const
The area-weighted sum of a field.
tmp< volVectorField > field_cellCentre() const
The cell centres - (swak = pos).
virtual const fvMesh & mesh() const
The mesh we are attached to.
bool isFaceData() const noexcept
A surface field.
bool isVolumeData() const noexcept
A volume field.
autoPtr< regIOobject > resultField_
The results (volume, surface, point).
bool isLogical_
A logical (bool-like) field (but actually a scalar).
const dimensionSet & dimensions() const noexcept
The preferred result dimensions (if any).
Type areaAverage(SurfaceField< Type > &fld) const
The area-weighted average of a field.
virtual ~parseDriver()=default
Destructor.
tmp< pointScalarField > field_pointSelection(const word &name, enum topoSetSource::sourceType setType) const
Point selections (as logical).
void setResult(VolumeField< Type > *ptr, bool logical=false)
Set result (vol field).
void clearField()
Clear out local copies of the field.
dimensionSet resultDimensions_
The result dimensions.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Generic interface code for Ragel/Lemon combination Subclasses should implement one or more process() ...
const std::string & content() const
Get reference to the input buffer content.
label nPoints() const noexcept
Number of mesh points.
A class for managing temporary objects.
Definition tmp.H:75
sourceType
Enumeration defining the types of sources.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
rDeltaTY field()
auto & name
Namespace for volume field expressions parsing and evaluation.
A namespace for expression-related classes/traits etc.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
GeometricField< Type, fvPatchField, volMesh > VolumeField
A volume field for a given type.
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
const direction noexcept
Definition scalarImpl.H:265
GeometricField< Type, pointPatchField, pointMesh > PointField
A point field for a given type.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
GeometricField< Type, fvsPatchField, surfaceMesh > SurfaceField
A (volume) surface field for a given type.
dictionary dict
Foam::surfaceFields.