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;
328 void setResult(VolumeField<Type>* ptr,
bool logical =
false);
332 void setResult(SurfaceField<Type>* ptr,
bool logical =
false);
336 void setResult(PointField<Type>* ptr,
bool logical =
false);
343 tmp<VolumeField<Type>>
344 newVolField(
const Type& val = pTraits<Type>::zero)
const;
348 tmp<SurfaceField<Type>>
353 tmp<PointField<Type>>
359 tmp<VolumeField<Type>>
360 getVolField(
const word& fldName,
bool getOldTime=
false);
364 tmp<SurfaceField<Type>>
369 tmp<PointField<Type>>
413 fld.mesh().magSf().primitiveField(),
424 fld.mesh().magSf().primitiveField(),
tmp< volScalarField > field_cellVolume() const
The cell volumes - (swak = vol)
tmp< PointField< Type > > newPointField(const Type &val=pTraits< Type >::zero) const
Return a new point field with the mesh nPoints size.
Type volAverage(VolumeField< Type > &fld) const
The volume-weighted average of a field.
tmp< surfaceScalarField > field_faceArea() const
The face area magnitudes [magSf] - (swak = area)
FieldAssociation fieldAssociation() const noexcept
The geometric field association.
void setInternalFieldResult(const Field< Type > &fld)
Deep-copy the internalField as a result.
tmp< PointField< Type > > getPointField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
label nPoints() const noexcept
Number of mesh points.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
The field association for mesh (patch/volume) values.
bool hasDimensions() const noexcept
Apply dimensions() to geometric fields.
ClassName("volumeExpr::driver")
autoPtr< regIOobject > dupZeroField() const
A zero-initialized field with the same type as the result field.
void setResult(VolumeField< Type > *ptr, bool logical=false)
Set result (vol field)
tmp< VolumeField< Type > > getVolField(const word &fldName, bool getOldTime=false)
Retrieve field (vol field)
Generic GeometricField class.
parseDriver(const parseDriver &)=delete
autoPtr< regIOobject > resultField_
The results (volume, surface, point)
virtual ~parseDriver()=default
Destructor.
virtual const fvMesh & mesh() const
The mesh we are attached to.
tmp< pointVectorField > field_pointField() const
The mesh point locations - (swak = pts)
tmp< SurfaceField< Type > > getSurfaceField(const word &fldName, bool getOldTime=false)
Retrieve field (surface field)
tmp< surfaceScalarField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
Face selections (as logical)
bool hasDimensions_
Requested use of dimensions.
void clearField()
Clear out local copies of the field.
dimensionedScalar pos(const dimensionedScalar &ds)
void operator=(const parseDriver &)=delete
tmp< SurfaceField< Type > > cellToFace(const VolumeField< Type > &field) const
Interpolate cell to face values.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
bool isVolumeData() const noexcept
A volume field.
const word & resultType() const noexcept
The result type-name.
bool isPointData() const noexcept
A point field.
const fvMesh & mesh_
The referenced mesh.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
tmp< pointScalarField > field_pointSet(const word &name) const
Point selection (set)
bool isFaceData() const noexcept
A surface field.
tmp< surfaceVectorField > field_faceCentre() const
The face centres - (swak = fpos)
A class for handling words, derived from Foam::string.
dimensionSet resultDimensions_
The result dimensions.
tmp< VolumeField< Type > > newVolField(const Type &val=pTraits< Type >::zero) const
Return a new volume field with the mesh size.
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Type areaAverage(SurfaceField< Type > &fld) const
The area-weighted average of a field.
virtual unsigned parse(const std::string &expr, size_t pos=0, size_t len=std::string::npos)
Execute the parser.
const GeoField * isResultType() const
Test if stored result pointer is the specified type.
tmp< surfaceScalarField > field_faceZone(const word &name) const
Face selection (zone)
virtual bool readDict(const dictionary &dict)
Read variables, tables etc.
expressions::FieldAssociation fieldGeoType_
A volume/surface/point field.
static Type weightedAverage(const scalarField &weights, const Field< Type > &fld)
The (global) weighted average of a field, with stabilisation.
const std::string & content() const
Get reference to the input buffer content.
genericRagelLemonDriver()
Construct null.
tmp< volScalarField > field_cellZone(const word &name) const
Cell selection (zone)
tmp< pointScalarField > field_pointZone(const word &name) const
Point selection (zone)
bool isLogical() const noexcept
A logical (bool-like) field. Actually stored as a scalar.
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))
tmp< surfaceScalarField > field_faceSet(const word &name) const
Face selection (set)
tmp< PointField< Type > > cellToPoint(const VolumeField< Type > &field) const
Interpolate cell to point values.
Type volSum(VolumeField< Type > &fld) const
The volume-weighted sum of a field.
static Type weightedSum(const scalarField &weights, const Field< Type > &fld)
The (global) weighted sum (integral) of a field.
tmp< VolumeField< Type > > pointToCell(const PointField< Type > &field) const
Interpolate point to cell values.
tmp< volScalarField > field_randGaussian(label seed=0) const
A Gaussian random field.
label nCells() const noexcept
Number of mesh cells.
tmp< surfaceVectorField > field_areaNormal() const
The face areas with their vector direction [Sf] - (swak = face)
tmp< SurfaceField< Type > > newSurfaceField(const Type &val=pTraits< Type >::zero) const
Return a new surface field with the mesh nInternalFaces size.
sourceType
Enumeration defining the types of sources.
const dictionary & dict() const noexcept
The dictionary with all input data/specification.
word resultType_
The result type-name.
Type areaSum(SurfaceField< Type > &fld) const
The area-weighted sum of a field.
tmp< volScalarField > field_cellSet(const word &name) const
Cell selection (set)
tmp< volVectorField > field_cellCentre() const
The cell centres - (swak = pos)
A class for managing temporary objects.
tmp< volScalarField > field_rand(label seed=0, bool gaussian=false) const
A uniform random field.
const dimensionSet & dimensions() const noexcept
The preferred result dimensions (if any)
bool isLogical_
A logical (bool-like) field (but actually a scalar)
virtual label pointSize() const
The point field size for the expression.
tmp< pointScalarField > field_pointSelection(const word &name, enum topoSetSource::sourceType setType) const
Point selections (as logical)
virtual label size() const
The natural field size for the expression.
tmp< volScalarField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
Cell selections (as logical)