Go to the documentation of this file.
70 #ifndef sampledIsoSurface_H
71 #define sampledIsoSurface_H
86 class sampledIsoSurface
99 const scalar mergeTol_;
108 const boundBox bounds_;
114 mutable word exposedPatchName_;
116 mutable autoPtr<isoSurface> surfPtr_;
122 mutable label prevTimeIndex_;
125 mutable autoPtr<volScalarField> storedVolFieldPtr_;
134 mutable autoPtr<fvMeshSubset> subMeshPtr_;
137 mutable autoPtr<volScalarField> storedVolSubFieldPtr_;
147 void getIsoFields()
const;
151 bool updateGeometry()
const;
List< label > labelList
A List of labels.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const List< label > & null()
Return a null List.
A class for handling words, derived from Foam::string.
virtual bool needsUpdate() const
Does the surface need an update?
A class for managing temporary objects.
virtual const vectorField & Cf() const
Face centres.
const vectorField & Sf() const
Face area vectors (normals)
const List< Face > & surfFaces() const
Return const access to the faces.
virtual const pointField & points() const
Points of surface.
virtual bool update()
Update the surface as required.
virtual void print(Ostream &) const
Write.
Mesh consisting of general polyhedral cells.
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual const labelList & zoneIds() const
Per-face zone/region information.
A surface formed by the iso value. After "Regularised Marching Tetrahedra: improved iso-surface extra...
const scalarField & magSf() const
Face area magnitudes.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
A sampledSurface defined by a surface of iso value. Always triangulated. To be used in sampleSurfaces...
An abstract class for surfaces with sampling.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Abstract base class for interpolation.
virtual bool expire()
Mark the surface as needing an update.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual const scalarField & magSf() const
Face area magnitudes.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
void getIsoField()
Lookup or read isoField.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
List< face > faceList
A List of faces.
const isoSurface & surface() const
virtual const faceList & faces() const
Faces of surface.
const word & name() const
Name of surface.
A List of wordRe with additional matching capabilities.
A bounding box defined in terms of min/max extrema points.
const vectorField & Cf() const
Face centres.
const polyMesh & mesh() const
Access to the underlying mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
bool interpolate() const
Interpolation to nodes requested for surface.
virtual ~sampledIsoSurface()
Destructor.
sampledIsoSurface(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
filterType
The filtering (regularization) to apply.
virtual const vectorField & Sf() const
Face area magnitudes.
TypeName("sampledIsoSurface")
Runtime type information.