MeshLib Documentation
Loading...
Searching...
No Matches
MRFastWindingNumber.h
Go to the documentation of this file.
1#pragma once
2
4#include "MRExpected.h"
5#include "MRId.h"
6
7namespace MR
8{
9
12{
13public:
14 virtual ~IFastWindingNumber() = default;
15
23 virtual void calcFromVector( std::vector<float>& res, const std::vector<Vector3f>& points, float beta, FaceId skipFace = {} ) = 0;
24
31 virtual bool calcSelfIntersections( FaceBitSet& res, float beta, ProgressCallback cb = {} ) = 0;
32
40 virtual Expected<void> calcFromGrid( std::vector<float>& res, const Vector3i& dims, const AffineXf3f& gridToMeshXf, float beta, ProgressCallback cb = {} ) = 0;
41
50 virtual Expected<void> calcFromGridWithDistances( std::vector<float>& res, const Vector3i& dims, const AffineXf3f& gridToMeshXf, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq, ProgressCallback cb ) = 0;
51};
52
58{
59public:
62 [[nodiscard]] MRMESH_API FastWindingNumber( const Mesh & mesh );
63
64 // see methods' descriptions in IFastWindingNumber
65 MRMESH_API void calcFromVector( std::vector<float>& res, const std::vector<Vector3f>& points, float beta, FaceId skipFace = {} ) override;
66 MRMESH_API bool calcSelfIntersections( FaceBitSet& res, float beta, ProgressCallback cb ) override;
67 MRMESH_API Expected<void> calcFromGrid( std::vector<float>& res, const Vector3i& dims, const AffineXf3f& gridToMeshXf, float beta, ProgressCallback cb ) override;
68 MRMESH_API float calcWithDistances( const Vector3f& p, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq );
69 MRMESH_API Expected<void> calcFromGridWithDistances( std::vector<float>& res, const Vector3i& dims, const AffineXf3f& gridToMeshXf, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq, ProgressCallback cb ) override;
70
71private:
72 [[nodiscard]] float calc_( const Vector3f & q, float beta, FaceId skipFace = {} ) const;
73
74 const Mesh & mesh_;
75 const AABBTree & tree_;
76 const Dipoles & dipoles_;
77};
78
79} // namespace MR
#define MRMESH_API
Definition MRMesh/MRMeshFwd.h:46
#define MRMESH_CLASS
Definition MRMesh/MRMeshFwd.h:50
Definition MRAABBTree.h:16
Definition MRFastWindingNumber.h:58
MRMESH_API Expected< void > calcFromGridWithDistances(std::vector< float > &res, const Vector3i &dims, const AffineXf3f &gridToMeshXf, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq, ProgressCallback cb) override
calculates distances with the sign obtained from winding number in each point from a three-dimensiona...
MRMESH_API float calcWithDistances(const Vector3f &p, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq)
MRMESH_API FastWindingNumber(const Mesh &mesh)
MRMESH_API bool calcSelfIntersections(FaceBitSet &res, float beta, ProgressCallback cb) override
calculates winding numbers for all centers of mesh's triangles. if winding number is less than 0 or g...
MRMESH_API void calcFromVector(std::vector< float > &res, const std::vector< Vector3f > &points, float beta, FaceId skipFace={}) override
calculates winding numbers in the points from given vector
MRMESH_API Expected< void > calcFromGrid(std::vector< float > &res, const Vector3i &dims, const AffineXf3f &gridToMeshXf, float beta, ProgressCallback cb) override
calculates winding numbers in each point from a three-dimensional grid
Abstract class for fast approximate computation of generalized winding number for a mesh (using its A...
Definition MRFastWindingNumber.h:12
virtual void calcFromVector(std::vector< float > &res, const std::vector< Vector3f > &points, float beta, FaceId skipFace={})=0
calculates winding numbers in the points from given vector
virtual Expected< void > calcFromGridWithDistances(std::vector< float > &res, const Vector3i &dims, const AffineXf3f &gridToMeshXf, float windingNumberThreshold, float beta, float maxDistSq, float minDistSq, ProgressCallback cb)=0
calculates distances with the sign obtained from winding number in each point from a three-dimensiona...
virtual Expected< void > calcFromGrid(std::vector< float > &res, const Vector3i &dims, const AffineXf3f &gridToMeshXf, float beta, ProgressCallback cb={})=0
calculates winding numbers in each point from a three-dimensional grid
virtual bool calcSelfIntersections(FaceBitSet &res, float beta, ProgressCallback cb={})=0
calculates winding numbers for all centers of mesh's triangles. if winding number is less than 0 or g...
virtual ~IFastWindingNumber()=default
std::function< bool(float)> ProgressCallback
Definition MRMesh/MRMeshFwd.h:576
Definition MRCameraOrientationPlugin.h:8
tl::expected< T, E > Expected
Definition MRExpected.h:58
Definition MRMesh/MRMesh.h:23