Go to the documentation of this file.
82 #ifndef expressions_volumeExprDriver_H
83 #define expressions_volumeExprDriver_H
108 public parsing::genericRagelLemonDriver,
109 public expressions::fvExprDriver
198 const dictionary&
dict
218 virtual const fvMesh&
mesh()
const
224 virtual label
size()
const
263 virtual
unsigned parse
265 const std::
string& expr,
267 size_t len = std::
string::npos
312 template<
class GeoField>
317 template<
class GeoField>
318 const GeoField*
isResultType(
bool logical,
bool dieOnNull=
false)
const;
330 GeometricField<Type, fvPatchField, volMesh>* ptr,
338 GeometricField<Type, fvsPatchField, surfaceMesh>* ptr,
443 fld.mesh().magSf().primitiveField(),
457 fld.mesh().magSf().primitiveField(),
dimensionSet resultDimensions_
The result dimensions.
tmp< pointScalarField > field_pointZone(const word &name) const
Point selection (zone)
sourceType
Enumeration defining the types of sources.
const word & resultType() const noexcept
The result type-name.
virtual label pointSize() const
The point field size for the expression.
genericRagelLemonDriver()
Construct null.
autoPtr< regIOobject > dupZeroField() const
A zero-initialized field with the same type as the result field.
const dimensionSet & dimensions() const noexcept
The preferred result dimensions (if any)
virtual label size() const
The natural field size for the expression.
A class for handling words, derived from Foam::string.
tmp< GeometricField< Type, pointPatchField, pointMesh > > cellToPoint(const GeometricField< Type, fvPatchField, volMesh > &field) const
Interpolate cell to point values.
const GeoField * isResultType() const
Test if stored result pointer is the specified type.
A class for managing temporary objects.
Driver for volume, surface, point field expressions.
virtual const fvMesh & mesh() const
The mesh we are attached to.
void operator=(const parseDriver &)=delete
void clearField()
Clear out local copies of the field.
parseDriver(const parseDriver &)=delete
const std::string & content() const
Get reference to the input buffer content.
Base driver for parsing value expressions associated with an fvMesh.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
tmp< volScalarField > field_rand(label seed=0, bool gaussian=false) const
A uniform random field.
const fvMesh & mesh_
The referenced mesh.
tmp< volScalarField > field_randGaussian(label seed=0) const
A Gaussian random field.
label nPoints() const noexcept
Number of mesh points.
tmp< volScalarField > field_cellVolume() const
The cell volumes - (swak = vol)
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
virtual ~parseDriver()=default
Destructor.
word resultType_
The result type-name.
label nCells() const noexcept
Number of mesh cells.
Type areaSum(GeometricField< Type, fvsPatchField, surfaceMesh > &fld) const
The area-weighted sum of a field.
Type volAverage(GeometricField< Type, fvPatchField, volMesh > &fld) const
The volume-weighted average of a field.
Generic templated field type.
Type volSum(GeometricField< Type, fvPatchField, volMesh > &fld) const
The volume-weighted sum of a field.
tmp< GeometricField< Type, pointPatchField, pointMesh > > getPointField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
static Type weightedAverage(const scalarField &weights, const Field< Type > &fld)
The (global) weighted average of a field, with stabilisation.
virtual bool readDict(const dictionary &dict)
Read variables, tables etc.
Generic interface code for Ragel/Lemon combination Subclasses should implement one or more process() ...
autoPtr< regIOobject > resultField_
The results (volume, surface, point)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > getSurfaceField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
tmp< surfaceVectorField > field_faceCentre() const
The face centres - (swak = fpos)
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool hasDimensions_
Requested use of dimensions.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
expressions::FieldAssociation fieldGeoType_
A volume/surface/point field.
Mesh data needed to do the Finite Volume discretisation.
void setInternalFieldResult(const Field< Type > &fld)
Deep-copy the internalField as a result.
tmp< pointVectorField > field_pointField() const
The mesh point locations - (swak = pts)
bool isLogical_
A logical (bool-like) field (but actually a scalar)
Type areaAverage(GeometricField< Type, fvsPatchField, surfaceMesh > &fld) const
The area-weighted average of a field.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
tmp< surfaceScalarField > field_faceArea() const
The face area magnitudes [magSf] - (swak = area)
tmp< pointScalarField > field_pointSet(const word &name) const
Point selection (set)
tmp< surfaceScalarField > field_faceZone(const word &name) const
Face selection (zone)
tmp< volScalarField > field_cellSet(const word &name) const
Cell selection (set)
const dictionary & dict() const noexcept
The dictionary with all input data/specification.
tmp< surfaceVectorField > field_areaNormal() const
The face areas with their vector direction [Sf] - (swak = face)
A traits class, which is primarily used for primitives.
tmp< GeometricField< Type, fvPatchField, volMesh > > getVolField(const word &fldName, bool getOldTime=false)
Retrieve field (vol field)
tmp< surfaceScalarField > field_faceSet(const word &name) const
Face selection (set)
tmp< GeometricField< Type, fvPatchField, volMesh > > pointToCell(const GeometricField< Type, pointPatchField, pointMesh > &field) const
Interpolate point to cell values.
bool hasDimensions() const noexcept
Apply dimensions() to geometric fields.
tmp< volScalarField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
Cell selections (as logical)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
tmp< volVectorField > field_cellCentre() const
The cell centres - (swak = pos)
tmp< GeometricField< Type, fvPatchField, volMesh > > newVolField(const Type &val=pTraits< Type >::zero) const
Return a new volume field with the mesh size.
void setResult(GeometricField< Type, fvPatchField, volMesh > *ptr, bool logical=false)
Set result (vol field)
tmp< pointScalarField > field_pointSelection(const word &name, enum topoSetSource::sourceType setType) const
Point selections (as logical)
tmp< volScalarField > field_cellZone(const word &name) const
Cell selection (zone)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > newSurfaceField(const Type &val=pTraits< Type >::zero) const
Return a new surface field with the mesh nInternalFaces size.
bool isLogical() const noexcept
A logical (bool-like) field. Actually stored as a scalar.
bool isPointData() const noexcept
A point field.
static Type weightedSum(const scalarField &weights, const Field< Type > &fld)
The (global) weighted sum (integral) of a field.
virtual unsigned parse(const std::string &expr, size_t pos=0, size_t len=std::string::npos)
Execute the parser.
tmp< surfaceScalarField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
Face selections (as logical)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > cellToFace(const GeometricField< Type, fvPatchField, volMesh > &field) const
Interpolate cell to face values.
ClassName("volumeExpr::driver")
bool isVolumeData() const noexcept
A volume field.
tmp< GeometricField< Type, pointPatchField, pointMesh > > newPointField(const Type &val=pTraits< Type >::zero) const
Return a new point field with the mesh nPoints size.
FieldAssociation fieldAssociation() const noexcept
The geometric field association.
bool isFaceData() const noexcept
A surface field.
dimensionedScalar pos(const dimensionedScalar &ds)