Go to the documentation of this file.
83 #ifndef sampledIsoSurface_H
84 #define sampledIsoSurface_H
99 class sampledIsoSurface
101 public sampledSurface
106 const word isoField_;
109 List<scalar> isoValues_;
112 isoSurfaceParams isoParams_;
127 mutable word exposedPatchName_;
133 mutable label prevTimeIndex_;
142 mutable autoPtr<isoSurfaceBase> isoSurfacePtr_;
148 mutable autoPtr<fvMeshSubset> subMeshPtr_;
151 mutable autoPtr<bitSet> ignoreCellsPtr_;
157 mutable autoPtr<volScalarField> storedVolFieldPtr_;
167 mutable autoPtr<volScalarField> storedVolSubFieldPtr_;
177 void getIsoFields()
const;
184 bool updateGeometry()
const;
215 return bool(isoSurfacePtr_);
268 return *isoSurfacePtr_;
278 return isoSurfacePtr_->meshCells();
323 virtual tmp<scalarField>
sample
325 const interpolation<scalar>& sampler
329 virtual tmp<vectorField>
sample
331 const interpolation<vector>& sampler
335 virtual tmp<sphericalTensorField>
sample
337 const interpolation<sphericalTensor>& sampler
List< label > labelList
A List of labels.
vectorField pointField
pointField is a vectorField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static const List< T > & 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.
bool hasIsoSurface() const
Is currently backed by an isoSurfacePtr_.
Mesh consisting of general polyhedral cells.
const labelList & meshCells() const
For each face, the original cell in mesh.
Field< vector > vectorField
Specialisation of Field<T> for vector.
const meshedSurface & surface() const
The currently created surface geometry.
virtual const labelList & zoneIds() const
Per-face zone/region information.
const scalarField & magSf() const
Face area magnitudes.
bool interpolate() const noexcept
Same as isPointData()
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
sampledIsoSurface(const isoSurfaceParams ¶ms, const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
A sampledSurface defined by a surface of iso value. It only recalculates the iso-surface if time chan...
An abstract class for surfaces with sampling.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
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,...
Preferences for controlling iso-surface algorithms.
OBJstream os(runTime.globalPath()/outputName)
virtual const scalarField & magSf() const
Face area magnitudes.
const word & name() const noexcept
Name of surface.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
MeshedSurface< face > meshedSurface
virtual void print(Ostream &os, int level=0) const
Print information.
virtual const faceList & faces() const
Faces of surface.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
A List of wordRe with additional matching capabilities.
const vectorField & Cf() const
Face centres.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
virtual ~sampledIsoSurface()
Destructor.
virtual const vectorField & Sf() const
Face area magnitudes.
TypeName("sampledIsoSurface")
Runtime type information.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField