Go to the documentation of this file.
43 #ifndef isoSurfaceBase_H
44 #define isoSurfaceBase_H
135 const uint8_t maskValue
271 #undef declareIsoSurfaceInterpolateMethod
272 #define declareIsoSurfaceInterpolateMethod(Type) \
274 virtual tmp<Field<Type>> \
277 const GeometricField<Type, fvPatchField, volMesh>& cellValues, \
278 const Field<Type>& pointValues \
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
type
Volume classification types.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams ¶ms, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
Low-level components common to various iso-surface algorithms.
void ignoreCyclics()
Set ignoreBoundaryFaces to ignore cyclics (cyclicACMI)
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A class for managing temporary objects.
const scalarField & cVals_
Cell values.
cutType getCellCutType(const label celli) const
#define declareIsoSurfaceInterpolateMethod(Type)
algorithmType
The algorithm types.
labelList & meshCells() noexcept
For each face, the original cell in mesh.
Mesh consisting of general polyhedral cells.
static label countCutType(const UList< cutType > &cuts, const uint8_t maskValue)
Count the number of cuts matching the mask type.
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
isoSurfaceBase(const isoSurfaceBase &)=delete
No copy construct.
const polyMesh & mesh_
Reference to mesh.
isoSurfaceParams(const algorithmType algo=algorithmType::ALGO_DEFAULT, const filterType filter=filterType::DIAGCELL) noexcept
Default construct, or with specified algorithm.
labelList meshCells_
For every face, the original cell in mesh.
const scalarField & pVals_
Point values.
cutType getFaceCutType(const label facei) const
Determine face cut for an individual face.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
const scalar iso_
Iso value.
void operator=(const isoSurfaceBase &)=delete
No copy assignment.
All edges to cell centre cut.
label calcCellCuts(List< cutType > &cuts) const
Populate a list of candidate cell cuts using getCellCutType()
Preferences for controlling iso-surface algorithms.
const labelList & meshCells() const noexcept
For each face, the original cell in mesh.
static void resetCuts(UList< cutType > &cuts)
Restore non-BLOCKED state to an UNVISITED state.
Vector< scalar > vector
A scalar version of the templated Vector.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
tmp< Field< Type > > interpolateTemplate(const GeometricField< Type, fvPatchField, volMesh > &cellValues, const Field< Type > &pointValues) const
Dummy templated interpolate method.
scalar isoValue() const noexcept
The iso-value associated with the surface.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
cutType
The type of cell/face cuts.
A bounding box defined in terms of min/max extrema points.
bitSet ignoreBoundaryFaces_
Optional boundary faces to ignore.
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
filterType
The filtering (regularization) to apply.