48 Foam::distanceSurface::topologyFilterType
50Foam::distanceSurface::topoFilterNames_
52 { topologyFilterType::NONE,
"none" },
53 { topologyFilterType::LARGEST_REGION,
"largestRegion" },
54 { topologyFilterType::NEAREST_POINTS,
"nearestPoints" },
55 { topologyFilterType::PROXIMITY_REGIONS,
"proximityRegions" },
56 { topologyFilterType::PROXIMITY_FACES,
"proximityFaces" },
57 { topologyFilterType::PROXIMITY_FACES,
"proximity" },
81 <<
"Had " << notHit <<
" faces/cells from "
82 << nearest.
size() <<
" without a point hit." <<
nl
83 <<
"May be caused by a severely degenerate input surface" <<
nl
113 const scalar normDist = (
diff & norm);
124 const List<pointIndexHit>& nearest,
164 const bitSet& ignoreLocation,
172 if (ignoreLocation.
test(i))
195template<
bool WantPo
intFilter = false>
201 const scalar boundBoxInflate = 0.1
216 const point& pt = nearest[celli].point();
219 cellBb.
inflate(boundBoxInflate);
223 ignoreCells.
set(celli);
225 else if (WantPointFilter)
228 ignorePoints.
unset(
mesh.cellPoints(celli));
243 const word& defaultSurfaceName,
256 dict.getOrDefault(
"surfaceName", defaultSurfaceName),
266 distance_(
dict.getOrDefault<scalar>(
"distance", 0)),
267 withZeroDistance_(
equal(distance_, 0)),
272 ||
dict.getOrDefault<bool>(
"signed", true)
283 topoFilterNames_.getOrDefault
287 topologyFilterType::NONE
291 maxDistanceSqr_(
Foam::
sqr(GREAT)),
292 absProximity_(
dict.getOrDefault<scalar>(
"absProximity", 1
e-5)),
293 cellDistancePtr_(nullptr),
297 isoSurfacePtr_(nullptr)
299 if (topologyFilterType::NEAREST_POINTS == topoFilter_)
301 dict.readEntry(
"nearestPoints", nearestPoints_);
304 if (
dict.readIfPresent(
"maxDistance", maxDistanceSqr_))
306 maxDistanceSqr_ = Foam::sqr(maxDistanceSqr_);
314 const word& surfaceType,
315 const word& surfaceName,
337 const word& surfaceType,
338 const word& surfaceName,
340 const bool useSignedDistance,
375 const bool useSignedDistance,
380 geometryPtr_(surface),
382 withZeroDistance_(
equal(distance_, 0)),
391 topoFilter_(topologyFilterType::NONE),
393 maxDistanceSqr_(
Foam::
sqr(GREAT)),
395 cellDistancePtr_(nullptr),
399 isoSurfacePtr_(nullptr)
409 Pout<<
"distanceSurface::createGeometry updating geometry." <<
endl;
413 isoSurfacePtr_.reset(
nullptr);
418 const searchableSurface& geom = geometryPtr_();
420 const fvMesh& fvmesh =
static_cast<const fvMesh&
>(mesh_);
425 cellDistancePtr_.reset
431 "distanceSurface.cellDistance",
432 fvmesh.time().timeName(),
442 auto& cellDistance = *cellDistancePtr_;
450 bitSet ignoreCells, ignoreCellPoints;
452 const bool filterCells =
464 List<pointIndexHit> nearest;
494 topologyFilterType::PROXIMITY_REGIONS == topoFilter_
501 if (withSignDistance_)
504 geom.getNormal(nearest, norms);
518 else if (withZeroDistance_)
529 calcAbsoluteDistance(
fld, cc, nearest);
536 forAll(fvmesh.C().boundaryField(), patchi)
538 const pointField& cc = fvmesh.C().boundaryField()[patchi];
550 if (withSignDistance_)
553 geom.getNormal(nearest, norms);
555 if (withZeroDistance_)
571 calcAbsoluteDistance(
fld, cc, nearest);
582 pointDistance_.resize(fvmesh.nPoints());
583 pointDistance_ = GREAT;
598 if (withSignDistance_)
601 geom.getNormal(nearest, norms);
615 else if (withZeroDistance_)
626 calcAbsoluteDistance(
fld,
pts, nearest);
632 if (ignoreCells.none())
634 ignoreCells.clearStorage();
636 else if (filterCells && topologyFilterType::NONE != topoFilter_)
647 if (topologyFilterType::LARGEST_REGION == topoFilter_)
649 refineBlockedCells(ignoreCells, isoCutter);
650 filterKeepLargestRegion(ignoreCells);
652 else if (topologyFilterType::NEAREST_POINTS == topoFilter_)
654 refineBlockedCells(ignoreCells, isoCutter);
655 filterKeepNearestRegions(ignoreCells);
662 ignoreCellPoints.clearStorage();
667 Pout<<
"Writing cell distance:" << cellDistance.objectPath() <<
endl;
668 cellDistance.write();
674 "distanceSurface.pointDistance",
675 fvmesh.time().timeName(),
684 pDist.primitiveFieldRef() = pointDistance_;
686 Pout<<
"Writing point distance:" << pDist.objectPath() <<
endl;
704 if (filterCells && topologyFilterType::PROXIMITY_REGIONS == topoFilter_)
714 refineBlockedCells(ignoreCells, isoCutter);
725 || topologyFilterType::PROXIMITY_FACES == topoFilter_
726 || topologyFilterType::PROXIMITY_REGIONS == topoFilter_
729 surface_.transfer(
static_cast<meshedSurface&
>(*isoSurfacePtr_));
730 meshCells_.transfer(isoSurfacePtr_->meshCells());
732 isoSurfacePtr_.reset(
nullptr);
733 cellDistancePtr_.reset(
nullptr);
734 pointDistance_.clear();
737 if (topologyFilterType::PROXIMITY_FACES == topoFilter_)
739 filterFaceProximity();
741 else if (topologyFilterType::PROXIMITY_REGIONS == topoFilter_)
743 filterRegionProximity(ignoreCells);
756 os <<
" surface:" << surfaceName()
758 <<
" topology:" << topoFilterNames_[topoFilter_];
760 isoParams_.print(
os);
764 os <<
" faces:" << surface().surfFaces().
size()
765 <<
" points:" << surface().points().
size();
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
fileName objectPath() const
The complete path + object name.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
void resize(const label numElem, const unsigned int val=0u)
Reset addressable list size, does not shrink the allocated size.
void clearStorage()
Clear the list and delete storage.
const point_type & point() const noexcept
Return point, no checks.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
void size(const label n)
Older name for setAddressableSize.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
void set(const bitSet &bitset)
Set specified bits from another bitset.
bool none() const
True if no bits in this bitset are set.
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
bitSet & unset(const bitSet &other)
Unset (subtract) the bits specified in the other bitset, which is a set difference corresponds to the...
A bounding box defined in terms of min/max extrema points.
bool contains(const point &pt) const
Contains point? (inside or on edge).
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A surface defined by a distance from an input searchable surface. Uses an iso-surface algorithm (cell...
void filterFaceProximity()
Adjust extracted iso-surface to remove far faces.
void createGeometry()
Create/recreate the distance surface.
void filterRegionProximity(bitSet &ignoreCells) const
Remove region(s) that have far faces.
scalar distance() const noexcept
The distance to the underlying searchableSurface.
distanceSurface(const word &defaultSurfaceName, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
void print(Ostream &os, int level=0) const
Print information.
void filterKeepNearestRegions(bitSet &ignoreCells) const
Keep region(s) closest to the nearest points.
const meshedSurface & surface() const
The underlying surface.
bool refineBlockedCells(bitSet &ignoreCells, const isoSurfaceBase &isoContext) const
Re-filter the blocked cells based on the anticipated cuts.
void filterKeepLargestRegion(bitSet &ignoreCells) const
Keep region with the most cuts (after region split).
const word & surfaceName() const
The name of the underlying searchableSurface.
Mesh data needed to do the Finite Volume discretisation.
const volVectorField & C() const
Return cell centres as volVectorField.
const Time & time() const
Return the top-level database.
Low-level components common to various iso-surface algorithms.
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.
Preferences for controlling iso-surface algorithms.
Mesh consisting of general polyhedral cells.
virtual const pointField & points() const
Return raw points.
label nPoints() const noexcept
Number of mesh points.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const =0
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const =0
From a set of points and indices get the normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Different types of constants.
Namespace for handling debugging switches.
static void calcNormalDistance_zero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
scalar distance(const vector &p1, const vector &p2)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
static scalar normalDistance_zero(const point &pt, const pointIndexHit &pHit, const vector &norm)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
bool equal(const T &a, const T &b)
Compare two values for equality.
static void calcNormalDistance_filtered(scalarField &distance, const bitSet &ignoreLocation, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
static scalar normalDistance_nonzero(const point &pt, const pointIndexHit &pHit, const vector &norm)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
MeshedSurface< face > meshedSurface
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
static void calcNormalDistance_nonzero(scalarField &distance, const pointField &points, const List< pointIndexHit > &nearest, const vectorField &normals)
static void checkAllHits(const UList< pointIndexHit > &nearest)
static bitSet simpleGeometricFilter(bitSet &ignoreCells, const List< pointIndexHit > &nearest, const polyMesh &mesh, const scalar boundBoxInflate=0.1)
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.