Go to the documentation of this file.
122 #ifndef distanceSurface_H
123 #define distanceSurface_H
139 class distanceSurface
144 enum class topologyFilterType : uint8_t
153 static const Enum<topologyFilterType> topoFilterNames_;
159 const polyMesh& mesh_;
162 const autoPtr<searchableSurface> geometryPtr_;
165 const scalar distance_;
168 const bool withZeroDistance_;
171 const bool withSignDistance_;
174 isoSurfaceParams isoParams_;
177 topologyFilterType topoFilter_;
183 scalar maxDistanceSqr_;
186 scalar absProximity_;
189 autoPtr<volScalarField> cellDistancePtr_;
211 static inline void calcAbsoluteDistance
232 return bool(isoSurfacePtr_);
239 const GeometricField<Type, fvPatchField, volMesh>& cellValues,
240 const Field<Type>& pointValues
245 return isoSurfacePtr_->interpolate(cellValues, pointValues);
259 const isoSurfaceBase& isoContext
286 const word& defaultSurfaceName,
287 const polyMesh&
mesh,
288 const dictionary&
dict
295 const word& surfaceType,
306 const word& surfaceType,
309 const bool useSignedDistance,
326 return geometryPtr_->name();
340 return *isoSurfacePtr_;
350 return *isoSurfacePtr_;
360 return isoSurfacePtr_->meshCells();
370 return isoSurfacePtr_->meshCells();
379 void print(Ostream& os)
const;
List< label > labelList
A List of labels.
vectorField pointField
pointField is a vectorField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A class for handling words, derived from Foam::string.
Low-level components common to various iso-surface algorithms.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
A class for managing temporary objects.
bool hasIsoSurface() const
Is currently backed by an isoSurfacePtr_.
distanceSurface(const word &defaultSurfaceName, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
void filterKeepNearestRegions(bitSet &ignoreCells) const
Keep region(s) closest to the nearest points.
scalar distance() const noexcept
The distance to the underlying searchableSurface.
Mesh consisting of general polyhedral cells.
#define forAll(list, i)
Loop across all elements in list.
bool refineBlockedCells(bitSet &ignoreCells, const isoSurfaceBase &isoContext) const
Re-filter the blocked cells based on the anticipated cuts.
tmp< Field< Type > > isoSurfaceInterpolate(const GeometricField< Type, fvPatchField, volMesh > &cellValues, const Field< Type > &pointValues) const
Interpolate volume field onto surface points.
A surface defined by a distance from an input searchable surface. Uses an iso-surface algorithm (cell...
const labelList & meshCells() const
For each face, the original cell in mesh.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
const word & surfaceName() const
The name of the underlying searchableSurface.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Preferences for controlling iso-surface algorithms.
void print(Ostream &os) const
Print information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
virtual ~distanceSurface()=default
Destructor.
MeshedSurface< face > meshedSurface
const meshedSurface & surface() const
The underlying surface.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void filterKeepLargestRegion(bitSet &ignoreCells) const
Keep region with the most cuts (after region split)
labelList & meshCells()
For each face, the original cell in mesh.
void createGeometry()
Create/recreate the distance surface.
TypeName("distanceSurface")
Runtime type information.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
bitSet filterPrepareRegionSplit(const bitSet &ignoreCells) const
Prepare blockedFaces for region split.
vector point
Point is a vector.
void filterByProximity()
Adjust extracted iso-surface to remove far faces.
meshedSurface & surface()
The underlying surface.