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,
443 fld.mesh().magSf().primitiveField(),
457 fld.mesh().magSf().primitiveField(),
dimensionSet resultDimensions_
tmp< pointScalarField > field_pointZone(const word &name) const
const word & resultType() const noexcept
virtual label pointSize() const
genericRagelLemonDriver()
autoPtr< regIOobject > dupZeroField() const
const dimensionSet & dimensions() const noexcept
virtual label size() const
A class for handling words, derived from Foam::string.
tmp< GeometricField< Type, pointPatchField, pointMesh > > cellToPoint(const GeometricField< Type, fvPatchField, volMesh > &field) const
const GeoField * isResultType() const
A class for managing temporary objects.
Driver for volume, surface, point field expressions.
virtual const fvMesh & mesh() const
void operator=(const parseDriver &)=delete
parseDriver(const parseDriver &)=delete
const std::string & content() const
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
tmp< volScalarField > field_randGaussian(label seed=0) const
@ VOLUME_DATA
Volume data.
label nPoints() const noexcept
tmp< volScalarField > field_cellVolume() const
static const dictionary null
virtual ~parseDriver()=default
label nCells() const noexcept
Type areaSum(GeometricField< Type, fvsPatchField, surfaceMesh > &fld) const
Type volAverage(GeometricField< Type, fvPatchField, volMesh > &fld) const
Type volSum(GeometricField< Type, fvPatchField, volMesh > &fld) const
tmp< GeometricField< Type, pointPatchField, pointMesh > > getPointField(const word &fldName, bool getOldTime=false)
static Type weightedAverage(const scalarField &weights, const Field< Type > &fld)
virtual bool readDict(const dictionary &dict)
Generic interface code for Ragel/Lemon combination Subclasses should implement one or more process() ...
autoPtr< regIOobject > resultField_
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > getSurfaceField(const word &fldName, bool getOldTime=false)
tmp< surfaceVectorField > field_faceCentre() const
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))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
expressions::FieldAssociation fieldGeoType_
Mesh data needed to do the Finite Volume discretisation.
void setInternalFieldResult(const Field< Type > &fld)
tmp< pointVectorField > field_pointField() const
Type areaAverage(GeometricField< Type, fvsPatchField, surfaceMesh > &fld) const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
tmp< surfaceScalarField > field_faceArea() const
tmp< pointScalarField > field_pointSet(const word &name) const
tmp< surfaceScalarField > field_faceZone(const word &name) const
tmp< volScalarField > field_cellSet(const word &name) const
const dictionary & dict() const noexcept
tmp< surfaceVectorField > field_areaNormal() const
A traits class, which is primarily used for primitives.
tmp< GeometricField< Type, fvPatchField, volMesh > > getVolField(const word &fldName, bool getOldTime=false)
tmp< surfaceScalarField > field_faceSet(const word &name) const
tmp< GeometricField< Type, fvPatchField, volMesh > > pointToCell(const GeometricField< Type, pointPatchField, pointMesh > &field) const
bool hasDimensions() const noexcept
tmp< volScalarField > field_cellSelection(const word &name, enum topoSetSource::sourceType setType) const
word name(const expressions::valueTypeCode typeCode)
tmp< volVectorField > field_cellCentre() const
tmp< GeometricField< Type, fvPatchField, volMesh > > newVolField(const Type &val=pTraits< Type >::zero) const
void setResult(GeometricField< Type, fvPatchField, volMesh > *ptr, bool logical=false)
tmp< pointScalarField > field_pointSelection(const word &name, enum topoSetSource::sourceType setType) const
tmp< volScalarField > field_cellZone(const word &name) const
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > newSurfaceField(const Type &val=pTraits< Type >::zero) const
bool isLogical() const noexcept
bool isPointData() const noexcept
Generic GeometricField class.
static Type weightedSum(const scalarField &weights, const Field< Type > &fld)
virtual unsigned parse(const std::string &expr, size_t pos=0, size_t len=std::string::npos)
tmp< surfaceScalarField > field_faceSelection(const word &name, enum topoSetSource::sourceType setType) const
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > cellToFace(const GeometricField< Type, fvPatchField, volMesh > &field) const
ClassName("volumeExpr::driver")
bool isVolumeData() const noexcept
tmp< GeometricField< Type, pointPatchField, pointMesh > > newPointField(const Type &val=pTraits< Type >::zero) const
FieldAssociation fieldAssociation() const noexcept
bool isFaceData() const noexcept
dimensionedScalar pos(const dimensionedScalar &ds)