MeshLib C++ Docs
Loading...
Searching...
No Matches
MRMesh/MRMesh.h
Go to the documentation of this file.
1#pragma once
2
4#include "MRMeshTopology.h"
5#include "MRMeshProject.h"
6#include "MREdgePoint.h"
7#include "MRLineSegm.h"
9#include "MRWriter.h"
10#include "MRConstants.h"
11#include "MRProgressCallback.h"
12#include <cfloat>
13
14namespace MR
15{
16
18
22struct [[nodiscard]] Mesh
23{
25 VertCoords points;
26
28 [[nodiscard]] MRMESH_API static Mesh fromTriangles(
29 VertCoords vertexCoordinates,
30 const Triangulation& t, const MeshBuilder::BuildSettings& settings = {}, ProgressCallback cb = {} );
31
33 [[nodiscard]] MRMESH_API static Mesh fromTriMesh(
34 TriMesh && triMesh,
35 const MeshBuilder::BuildSettings& settings = {}, ProgressCallback cb = {} );
36
40 VertCoords vertexCoordinates,
41 Triangulation & t,
42 std::vector<MeshBuilder::VertDuplication> * dups = nullptr,
43 const MeshBuilder::BuildSettings & settings = {} );
44
48 [[nodiscard]] MRMESH_API static Mesh fromFaceSoup(
49 VertCoords vertexCoordinates,
50 const std::vector<VertId> & verts, const Vector<MeshBuilder::VertSpan, FaceId> & faces,
51 const MeshBuilder::BuildSettings& settings = {}, ProgressCallback cb = {} );
52
56 [[nodiscard]] MRMESH_API static Mesh fromPointTriples( const std::vector<Triangle3f> & posTriples, bool duplicateNonManifoldVertices );
57
59 [[nodiscard]] MRMESH_API bool operator ==( const Mesh & b ) const;
60
62 [[nodiscard]] Vector3f orgPnt( EdgeId e ) const { return points[ topology.org( e ) ]; }
63
65 [[nodiscard]] Vector3f destPnt( EdgeId e ) const { return points[ topology.dest( e ) ]; }
66
68 [[nodiscard]] Vector3f edgeVector( EdgeId e ) const { return destPnt( e ) - orgPnt( e ); }
69
71 [[nodiscard]] LineSegm3f edgeSegment( EdgeId e ) const { return { orgPnt( e ), destPnt( e ) }; }
72
74 [[nodiscard]] Vector3f edgePoint( EdgeId e, float f ) const { return f * destPnt( e ) + ( 1 - f ) * orgPnt( e ); }
75
77 [[nodiscard]] Vector3f edgePoint( const MeshEdgePoint & ep ) const { return edgePoint( ep.e, ep.a ); }
78
80 [[nodiscard]] Vector3f edgeCenter( UndirectedEdgeId e ) const { return edgePoint( e, 0.5f ); }
81
83 MRMESH_API void getLeftTriPoints( EdgeId e, Vector3f & v0, Vector3f & v1, Vector3f & v2 ) const;
84
86 void getLeftTriPoints( EdgeId e, Vector3f (&v)[3] ) const { getLeftTriPoints( e, v[0], v[1], v[2] ); }
87
89 [[nodiscard]] Triangle3f getLeftTriPoints( EdgeId e ) const { Triangle3f res; getLeftTriPoints( e, res[0], res[1], res[2] ); return res; }
90
92 void getTriPoints( FaceId f, Vector3f & v0, Vector3f & v1, Vector3f & v2 ) const { getLeftTriPoints( topology.edgeWithLeft( f ), v0, v1, v2 ); }
93
95 void getTriPoints( FaceId f, Vector3f (&v)[3] ) const { getTriPoints( f, v[0], v[1], v[2] ); }
96
98 [[nodiscard]] Triangle3f getTriPoints( FaceId f ) const { Triangle3f res; getTriPoints( f, res[0], res[1], res[2] ); return res; }
99
101 [[nodiscard]] MRMESH_API Vector3f triPoint( const MeshTriPoint & p ) const;
102
104 [[nodiscard]] MRMESH_API Vector3f triCenter( FaceId f ) const;
105
107 [[nodiscard]] MRMESH_API float triangleAspectRatio( FaceId f ) const;
108
110 [[nodiscard]] MRMESH_API float circumcircleDiameterSq( FaceId f ) const;
111
113 [[nodiscard]] MRMESH_API float circumcircleDiameter( FaceId f ) const;
114
116 [[nodiscard]] MRMESH_API MeshTriPoint toTriPoint( VertId v ) const;
117
119 [[nodiscard]] MRMESH_API MeshTriPoint toTriPoint( FaceId f, const Vector3f & p ) const;
120
122 [[nodiscard]] MRMESH_API MeshTriPoint toTriPoint( const PointOnFace& p ) const;
123
125 [[nodiscard]] MRMESH_API MeshEdgePoint toEdgePoint( VertId v ) const;
126
128 [[nodiscard]] MRMESH_API MeshEdgePoint toEdgePoint( EdgeId e, const Vector3f & p ) const;
129
131 [[nodiscard]] MRMESH_API VertId getClosestVertex( const PointOnFace & p ) const;
132
134 [[nodiscard]] VertId getClosestVertex( const MeshTriPoint & p ) const { return getClosestVertex( PointOnFace{ topology.left( p.e ), triPoint( p ) } ); }
135
137 [[nodiscard]] MRMESH_API UndirectedEdgeId getClosestEdge( const PointOnFace & p ) const;
138
140 [[nodiscard]] UndirectedEdgeId getClosestEdge( const MeshTriPoint & p ) const { return getClosestEdge( PointOnFace{ topology.left( p.e ), triPoint( p ) } ); }
141
143 [[nodiscard]] float edgeLength( UndirectedEdgeId e ) const { return edgeVector( e ).length(); }
144
146 [[nodiscard]] float edgeLengthSq( UndirectedEdgeId e ) const { return edgeVector( e ).lengthSq(); }
147
149 [[nodiscard]] MRMESH_API Vector3f leftDirDblArea( EdgeId e ) const;
150
152 [[nodiscard]] Vector3f dirDblArea( FaceId f ) const { return leftDirDblArea( topology.edgeWithLeft( f ) ); }
153
155 [[nodiscard]] float dblArea( FaceId f ) const { return dirDblArea( f ).length(); }
156
158 [[nodiscard]] float area( FaceId f ) const { return 0.5f * dblArea( f ); }
159
161 [[nodiscard]] MRMESH_API double area( const FaceBitSet & fs ) const;
162
164 [[nodiscard]] double area( const FaceBitSet * fs = nullptr ) const { return area( topology.getFaceIds( fs ) ); }
165
167 [[nodiscard]] MRMESH_API Vector3d dirArea( const FaceBitSet & fs ) const;
168
170 [[nodiscard]] Vector3d dirArea( const FaceBitSet * fs = nullptr ) const { return dirArea( topology.getFaceIds( fs ) ); }
171
173 [[nodiscard]] MRMESH_API double projArea( const Vector3f & dir, const FaceBitSet & fs ) const;
174
176 [[nodiscard]] double projArea( const Vector3f & dir, const FaceBitSet * fs = nullptr ) const { return projArea( dir, topology.getFaceIds( fs ) ); }
177
180 [[nodiscard]] MRMESH_API double volume( const FaceBitSet* region = nullptr ) const;
181
183 [[nodiscard]] MRMESH_API double holePerimiter( EdgeId e ) const;
184
187 [[nodiscard]] MRMESH_API Vector3d holeDirArea( EdgeId e ) const;
188
190 [[nodiscard]] MRMESH_API Vector3f leftTangent( EdgeId e ) const;
191
193 [[nodiscard]] Vector3f leftNormal( EdgeId e ) const { return leftDirDblArea( e ).normalized(); }
194
196 [[nodiscard]] Vector3f normal( FaceId f ) const { return dirDblArea( f ).normalized(); }
197
199 [[nodiscard]] MRMESH_API Vector3f dirDblArea( VertId v ) const;
200
202 [[nodiscard]] float dblArea( VertId v ) const { return dirDblArea( v ).length(); }
203
205 [[nodiscard]] Vector3f normal( VertId v ) const { return dirDblArea( v ).normalized(); }
206
209 [[nodiscard]] MRMESH_API Vector3f normal( const MeshTriPoint & p ) const;
210
213 [[nodiscard]] MRMESH_API Vector3f pseudonormal( VertId v, const FaceBitSet * region = nullptr ) const;
214
216 [[nodiscard]] MRMESH_API Vector3f pseudonormal( UndirectedEdgeId e, const FaceBitSet * region = nullptr ) const;
217
222 [[nodiscard]] MRMESH_API Vector3f pseudonormal( const MeshTriPoint & p, const FaceBitSet * region = nullptr ) const;
223
227 [[nodiscard]] MRMESH_API float signedDistance( const Vector3f & pt, const MeshProjectionResult & proj, const FaceBitSet * region = nullptr ) const;
228 [[deprecated]] MRMESH_API float signedDistance( const Vector3f & pt, const MeshTriPoint & proj, const FaceBitSet * region = nullptr ) const;
229
233 [[nodiscard]] MRMESH_API float signedDistance( const Vector3f & pt ) const;
234
239 [[nodiscard]] MRMESH_API std::optional<float> signedDistance( const Vector3f & pt, float maxDistSq, const FaceBitSet * region = nullptr ) const;
240
246 [[nodiscard]] MRMESH_API float calcFastWindingNumber( const Vector3f & pt, float beta = 2 ) const;
247
250 [[nodiscard]] bool isOutside( const Vector3f & pt, float windingNumberThreshold = 0.5f, float beta = 2 ) const { return calcFastWindingNumber( pt, beta ) <= windingNumberThreshold; }
251
255 [[nodiscard]] MRMESH_API bool isOutsideByProjNorm( const Vector3f & pt, const MeshProjectionResult & proj, const FaceBitSet * region = nullptr ) const;
256
258 [[nodiscard]] MRMESH_API float sumAngles( VertId v, bool * outBoundaryVert = nullptr ) const;
259
261 [[nodiscard]] MRMESH_API Expected<VertBitSet> findSpikeVertices( float minSumAngle, const VertBitSet* region = nullptr, ProgressCallback cb = {} ) const;
262
267 [[nodiscard]] MRMESH_API float dihedralAngleSin( UndirectedEdgeId e ) const;
268
273 [[nodiscard]] MRMESH_API float dihedralAngleCos( UndirectedEdgeId e ) const;
274
280 [[nodiscard]] MRMESH_API float dihedralAngle( UndirectedEdgeId e ) const;
281
284 [[nodiscard]] MRMESH_API float discreteMeanCurvature( VertId v ) const;
285
288 [[nodiscard]] MRMESH_API float discreteMeanCurvature( UndirectedEdgeId e ) const;
289
293 [[nodiscard]] float discreteGaussianCurvature( VertId v, bool * outBoundaryVert = nullptr ) const { return 2 * PI_F - sumAngles( v, outBoundaryVert ); }
294
296 [[nodiscard]] MRMESH_API UndirectedEdgeBitSet findCreaseEdges( float angleFromPlanar ) const;
297
300 [[nodiscard]] MRMESH_API float leftCotan( EdgeId e ) const;
301
304 [[nodiscard]] float cotan( UndirectedEdgeId ue ) const { EdgeId e{ ue }; return leftCotan( e ) + leftCotan( e.sym() ); }
305
309 [[nodiscard]] MRMESH_API QuadraticForm3f quadraticForm( VertId v, bool angleWeigted,
310 const FaceBitSet * region = nullptr, const UndirectedEdgeBitSet * creases = nullptr ) const;
311
314 [[nodiscard]] MRMESH_API Box3f computeBoundingBox( const AffineXf3f * toWorld = nullptr ) const;
315
318 [[nodiscard]] MRMESH_API Box3f getBoundingBox() const;
319
322 [[nodiscard]] MRMESH_API Box3f computeBoundingBox( const FaceBitSet* region, const AffineXf3f* toWorld = nullptr ) const;
323
325 [[nodiscard]] MRMESH_API float averageEdgeLength() const;
326
328 [[nodiscard]] MRMESH_API Vector3f findCenterFromPoints() const;
329
331 [[nodiscard]] MRMESH_API Vector3f findCenterFromFaces() const;
332
334 [[nodiscard]] MRMESH_API Vector3f findCenterFromBBox() const;
335
338
341 MRMESH_API void transform( const AffineXf3f& xf, const VertBitSet* region = nullptr );
342
344 MRMESH_API VertId addPoint( const Vector3f & pos );
345
348 MRMESH_API EdgeId addSeparateEdgeLoop(const std::vector<Vector3f>& contourPoints);
349
352 MRMESH_API EdgeId addSeparateContours( const Contours3f& contours, const AffineXf3f* xf = nullptr );
353
358 MRMESH_API void attachEdgeLoopPart( EdgeId first, EdgeId last, const std::vector<Vector3f>& contourPoints );
359
368 MRMESH_API EdgeId splitEdge( EdgeId e, const Vector3f & newVertPos, FaceBitSet * region = nullptr, FaceHashMap * new2Old = nullptr );
369 // same, but split given edge on two equal parts
370 EdgeId splitEdge( EdgeId e, FaceBitSet * region = nullptr, FaceHashMap * new2Old = nullptr ) { return splitEdge( e, edgeCenter( e ), region, new2Old ); }
371
375 MRMESH_API VertId splitFace( FaceId f, const Vector3f & newVertPos, FaceBitSet * region = nullptr, FaceHashMap * new2Old = nullptr );
376 // same, putting new vertex in the centroid of original triangle
377 VertId splitFace( FaceId f, FaceBitSet * region = nullptr, FaceHashMap * new2Old = nullptr ) { return splitFace( f, triCenter( f ), region, new2Old ); }
378
380 MRMESH_API void addMesh( const Mesh & from,
381 // optionally returns mappings: from.id -> this.id
382 FaceMap * outFmap = nullptr, VertMap * outVmap = nullptr, WholeEdgeMap * outEmap = nullptr, bool rearrangeTriangles = false );
383 [[deprecated]] void addPart( const Mesh & from, FaceMap * outFmap = nullptr, VertMap * outVmap = nullptr, WholeEdgeMap * outEmap = nullptr, bool rearrangeTriangles = false )
384 { addMesh( from, outFmap, outVmap, outEmap, rearrangeTriangles ); }
385
387 MRMESH_API void addMeshPart( const MeshPart & from, const PartMapping & map );
388 [[deprecated]] void addPartByMask( const Mesh & from, const FaceBitSet & fromFaces, const PartMapping & map ) { addMeshPart( { from, &fromFaces }, map ); }
389
392 MRMESH_API void addMeshPart( const MeshPart & from, bool flipOrientation = false,
393 const std::vector<EdgePath> & thisContours = {}, // contours on this mesh that have to be stitched with
394 const std::vector<EdgePath> & fromContours = {}, // contours on from mesh during addition
395 // optionally returns mappings: from.id -> this.id
396 const PartMapping & map = {} );
397 [[deprecated]] void addPartByMask( const Mesh & from, const FaceBitSet & fromFaces, bool flipOrientation = false,
398 const std::vector<EdgePath> & thisContours = {}, const std::vector<EdgePath> & fromContours = {}, const PartMapping & map = {} )
399 { addMeshPart( { from, &fromFaces }, flipOrientation, thisContours, fromContours, map ); }
400
402 MRMESH_API void addPartByFaceMap( const Mesh & from, const FaceMap & fromFaces, bool flipOrientation = false,
403 const std::vector<EdgePath> & thisContours = {}, // contours on this mesh that have to be stitched with
404 const std::vector<EdgePath> & fromContours = {}, // contours on from mesh during addition
405 // optionally returns mappings: from.id -> this.id
406 const PartMapping & map = {} );
407
409 template<typename I>
410 MRMESH_API void addPartBy( const Mesh & from, I fbegin, I fend, size_t fcount, bool flipOrientation = false,
411 const std::vector<EdgePath> & thisContours = {},
412 const std::vector<EdgePath> & fromContours = {},
413 PartMapping map = {} );
414
416 MRMESH_API Mesh cloneRegion( const FaceBitSet & region, bool flipOrientation = false, const PartMapping & map = {} ) const;
417
420 MRMESH_API void pack( FaceMap * outFmap = nullptr, VertMap * outVmap = nullptr, WholeEdgeMap * outEmap = nullptr, bool rearrangeTriangles = false );
421
424 MRMESH_API PackMapping packOptimally( bool preserveAABBTree = true );
426
428 MRMESH_API void deleteFaces( const FaceBitSet & fs, const UndirectedEdgeBitSet * keepEdges = nullptr );
429
436 [[nodiscard]] MRMESH_API bool projectPoint( const Vector3f& point, PointOnFace& res, float maxDistSq = FLT_MAX, const FaceBitSet* region = nullptr, const AffineXf3f * xf = nullptr ) const;
437
444 [[nodiscard]] MRMESH_API bool projectPoint( const Vector3f& point, MeshProjectionResult& res, float maxDistSq = FLT_MAX, const FaceBitSet* region = nullptr, const AffineXf3f * xf = nullptr ) const;
445 [[nodiscard]] bool findClosestPoint( const Vector3f& point, MeshProjectionResult& res, float maxDistSq = FLT_MAX, const FaceBitSet* region = nullptr, const AffineXf3f * xf = nullptr ) const { return projectPoint( point, res, maxDistSq, region, xf ); }
446
453 [[nodiscard]] MRMESH_API std::optional<MeshProjectionResult> projectPoint( const Vector3f& point, float maxDistSq = FLT_MAX, const FaceBitSet * region = nullptr, const AffineXf3f * xf = nullptr ) const;
454 [[nodiscard]] std::optional<MeshProjectionResult> findClosestPoint( const Vector3f& point, float maxDistSq = FLT_MAX, const FaceBitSet * region = nullptr, const AffineXf3f * xf = nullptr ) const { return projectPoint( point, maxDistSq, region, xf ); }
455
458
460 [[nodiscard]] const AABBTree * getAABBTreeNotCreate() const { return AABBTreeOwner_.get(); }
461
464
466 [[nodiscard]] const AABBTreePoints * getAABBTreePointsNotCreate() const { return AABBTreePointsOwner_.get(); }
467
469 MRMESH_API const Dipoles & getDipoles() const;
470
472 [[nodiscard]] const Dipoles * getDipolesNotCreate() const { return dipolesOwner_.get(); }
473
476 MRMESH_API void invalidateCaches( bool pointsChanged = true );
477
481 MRMESH_API void updateCaches( const VertBitSet & changedVerts );
482
483 // returns the amount of memory this object occupies on heap
484 [[nodiscard]] MRMESH_API size_t heapBytes() const;
485
488
490 MRMESH_API void mirror( const Plane3f& plane );
491
492private:
493 mutable SharedThreadSafeOwner<AABBTree> AABBTreeOwner_;
494 mutable SharedThreadSafeOwner<AABBTreePoints> AABBTreePointsOwner_;
495 mutable SharedThreadSafeOwner<Dipoles> dipolesOwner_;
496};
497
498} //namespace MR
constexpr bool operator==(ImVec2 a, ImVec2 b)
Definition MRImGuiVectorOperators.h:117
#define MRMESH_API
Definition MRMesh/MRMeshFwd.h:79
bounding volume hierarchy for point cloud structure
Definition MRAABBTreePoints.h:16
Definition MRAABBTree.h:16
Definition MRMesh/MRMeshTopology.h:18
FaceId left(EdgeId he) const
returns left face of half-edge
Definition MRMesh/MRMeshTopology.h:92
VertId dest(EdgeId he) const
returns destination vertex of half-edge
Definition MRMesh/MRMeshTopology.h:89
VertId org(EdgeId he) const
returns origin vertex of half-edge
Definition MRMesh/MRMeshTopology.h:86
EdgeId edgeWithLeft(FaceId a) const
returns valid edge if given vertex is present in the mesh
Definition MRMesh/MRMeshTopology.h:217
const FaceBitSet & getFaceIds(const FaceBitSet *region) const
if region pointer is not null then converts it in reference, otherwise returns all valid faces in the...
Definition MRMesh/MRMeshTopology.h:271
Definition MRSharedThreadSafeOwner.h:21
std::vector<T>-like container that requires specific indexing type,
Definition MRMesh/MRVector.h:20
std::function< bool(float)> ProgressCallback
Definition MRMesh/MRMeshFwd.h:626
tl::expected< T, E > Expected
Definition MRExpected.h:59
T dblArea(const Vector3< T > &p, const Vector3< T > &q, const Vector3< T > &r)
computes twice the area of given triangle
Definition MRTriMath.h:153
HashMap< FaceId, FaceId > FaceHashMap
Definition MRMesh/MRMeshFwd.h:509
Vector3< T > dirDblArea(const Triangle3< T > &t)
computes directed double area of given triangle
Definition MRTriMath.h:125
Triangle3< float > Triangle3f
Definition MRMesh/MRMeshFwd.h:382
Contours3< float > Contours3f
Definition MRMesh/MRMeshFwd.h:318
encodes a point on an edge of mesh or of polyline
Definition MREdgePoint.h:11
SegmPointf a
a in [0,1], a=0 => point is in org( e ), a=1 => point is in dest( e )
Definition MREdgePoint.h:13
EdgeId e
Definition MREdgePoint.h:12
Definition MRMeshBuilderTypes.h:30
Definition MRMesh/MRMeshProject.h:18
Definition MRMesh/MRMeshTriPoint.h:23
EdgeId e
Definition MRMesh/MRMeshTriPoint.h:24
Definition MRMesh/MRMesh.h:23
MRMESH_API float signedDistance(const Vector3f &pt, const MeshProjectionResult &proj, const FaceBitSet *region=nullptr) const
Vector3d dirArea(const FaceBitSet *fs=nullptr) const
computes the sum of directed areas for faces from given region (or whole mesh)
Definition MRMesh/MRMesh.h:170
float area(FaceId f) const
returns the area of given face
Definition MRMesh/MRMesh.h:158
MRMESH_API void addMesh(const Mesh &from, FaceMap *outFmap=nullptr, VertMap *outVmap=nullptr, WholeEdgeMap *outEmap=nullptr, bool rearrangeTriangles=false)
appends another mesh as separate connected component(s) to this
MRMESH_API float dihedralAngleSin(UndirectedEdgeId e) const
VertId getClosestVertex(const MeshTriPoint &p) const
returns one of three face vertices, closest to given point
Definition MRMesh/MRMesh.h:134
MRMESH_API Box3f computeBoundingBox(const FaceBitSet *region, const AffineXf3f *toWorld=nullptr) const
MRMESH_API Vector3f normal(const MeshTriPoint &p) const
MRMESH_API double area(const FaceBitSet &fs) const
computes the area of given face-region
MRMESH_API float calcFastWindingNumber(const Vector3f &pt, float beta=2) const
MRMESH_API bool isOutsideByProjNorm(const Vector3f &pt, const MeshProjectionResult &proj, const FaceBitSet *region=nullptr) const
MRMESH_API void shrinkToFit()
requests the removal of unused capacity
MRMESH_API Vector3f findCenterFromFaces() const
computes center of mass considering that density of all triangles is the same
MRMESH_API const AABBTree & getAABBTree() const
returns cached aabb-tree for this mesh, creating it if it did not exist in a thread-safe manner
float discreteGaussianCurvature(VertId v, bool *outBoundaryVert=nullptr) const
Definition MRMesh/MRMesh.h:293
Vector3f edgeVector(EdgeId e) const
returns vector equal to edge destination point minus edge origin point
Definition MRMesh/MRMesh.h:68
MRMESH_API QuadraticForm3f quadraticForm(VertId v, bool angleWeigted, const FaceBitSet *region=nullptr, const UndirectedEdgeBitSet *creases=nullptr) const
MRMESH_API std::optional< MeshProjectionResult > projectPoint(const Vector3f &point, float maxDistSq=FLT_MAX, const FaceBitSet *region=nullptr, const AffineXf3f *xf=nullptr) const
MRMESH_API UndirectedEdgeId getClosestEdge(const PointOnFace &p) const
returns one of three face edges, closest to given point
MRMESH_API Vector3d dirArea(const FaceBitSet &fs) const
computes the sum of directed areas for faces from given region
float dblArea(FaceId f) const
returns twice the area of given face
Definition MRMesh/MRMesh.h:155
MRMESH_API void deleteFaces(const FaceBitSet &fs, const UndirectedEdgeBitSet *keepEdges=nullptr)
deletes multiple given faces, also deletes adjacent edges and vertices if they were not shared by rem...
MRMESH_API float signedDistance(const Vector3f &pt, const MeshTriPoint &proj, const FaceBitSet *region=nullptr) const
MRMESH_API Vector3f triPoint(const MeshTriPoint &p) const
computes coordinates of point given as face and barycentric representation
MRMESH_API MeshTriPoint toTriPoint(VertId v) const
converts vertex into barycentric representation
MRMESH_API void addMeshPart(const MeshPart &from, const PartMapping &map)
appends whole or part of another mesh as separate connected component(s) to this
MRMESH_API Mesh cloneRegion(const FaceBitSet &region, bool flipOrientation=false, const PartMapping &map={}) const
creates new mesh from given triangles of this mesh
MRMESH_API Box3f computeBoundingBox(const AffineXf3f *toWorld=nullptr) const
MRMESH_API Vector3f pseudonormal(UndirectedEdgeId e, const FaceBitSet *region=nullptr) const
computes normalized half sum of face normals sharing given edge (only (region) faces will be consider...
MRMESH_API UndirectedEdgeBitSet findCreaseEdges(float angleFromPlanar) const
finds all mesh edges where dihedral angle is distinct from planar PI angle on at least given value
const AABBTree * getAABBTreeNotCreate() const
returns cached aabb-tree for this mesh, but does not create it if it did not exist
Definition MRMesh/MRMesh.h:460
MRMESH_API void addPartBy(const Mesh &from, I fbegin, I fend, size_t fcount, bool flipOrientation=false, const std::vector< EdgePath > &thisContours={}, const std::vector< EdgePath > &fromContours={}, PartMapping map={})
both addPartByMask and addPartByFaceMap call this general implementation
MRMESH_API Vector3f pseudonormal(VertId v, const FaceBitSet *region=nullptr) const
MRMESH_API std::optional< float > signedDistance(const Vector3f &pt, float maxDistSq, const FaceBitSet *region=nullptr) const
MRMESH_API Vector3f dirDblArea(VertId v) const
computes sum of directed double areas of all triangles around given vertex
MRMESH_API const Dipoles & getDipoles() const
returns cached dipoles of aabb-tree nodes for this mesh, creating it if it did not exist in a thread-...
static MRMESH_API Mesh fromTriangles(VertCoords vertexCoordinates, const Triangulation &t, const MeshBuilder::BuildSettings &settings={}, ProgressCallback cb={})
construct mesh from vertex coordinates and a set of triangles with given ids
double area(const FaceBitSet *fs=nullptr) const
computes the area of given face-region (or whole mesh)
Definition MRMesh/MRMesh.h:164
MRMESH_API VertId splitFace(FaceId f, const Vector3f &newVertPos, FaceBitSet *region=nullptr, FaceHashMap *new2Old=nullptr)
MRMESH_API float triangleAspectRatio(FaceId f) const
returns aspect ratio of given mesh triangle equal to the ratio of the circum-radius to twice its in-r...
MRMESH_API Vector3f findCenterFromPoints() const
computes average position of all valid mesh vertices
float edgeLengthSq(UndirectedEdgeId e) const
returns squared Euclidean length of the edge (faster to compute than length)
Definition MRMesh/MRMesh.h:146
MRMESH_API VertId getClosestVertex(const PointOnFace &p) const
returns one of three face vertices, closest to given point
MRMESH_API void zeroUnusedPoints()
for all points not in topology.getValidVerts() sets coordinates to (0,0,0)
static MRMESH_API Mesh fromPointTriples(const std::vector< Triangle3f > &posTriples, bool duplicateNonManifoldVertices)
Triangle3f getTriPoints(FaceId f) const
returns three points of given face
Definition MRMesh/MRMesh.h:98
EdgeId splitEdge(EdgeId e, FaceBitSet *region=nullptr, FaceHashMap *new2Old=nullptr)
Definition MRMesh/MRMesh.h:370
Vector3f edgePoint(EdgeId e, float f) const
returns a point on the edge: origin point for f=0 and destination point for f=1
Definition MRMesh/MRMesh.h:74
static MRMESH_API Mesh fromTrianglesDuplicatingNonManifoldVertices(VertCoords vertexCoordinates, Triangulation &t, std::vector< MeshBuilder::VertDuplication > *dups=nullptr, const MeshBuilder::BuildSettings &settings={})
Vector3f orgPnt(EdgeId e) const
returns coordinates of the edge origin
Definition MRMesh/MRMesh.h:62
Vector3f destPnt(EdgeId e) const
returns coordinates of the edge destination
Definition MRMesh/MRMesh.h:65
MRMESH_API Expected< VertBitSet > findSpikeVertices(float minSumAngle, const VertBitSet *region=nullptr, ProgressCallback cb={}) const
returns vertices where the sum of triangle angles is below given threshold
MRMESH_API double projArea(const Vector3f &dir, const FaceBitSet &fs) const
computes the sum of absolute projected area of faces from given region as visible if look from given ...
void addPartByMask(const Mesh &from, const FaceBitSet &fromFaces, bool flipOrientation=false, const std::vector< EdgePath > &thisContours={}, const std::vector< EdgePath > &fromContours={}, const PartMapping &map={})
Definition MRMesh/MRMesh.h:397
MRMESH_API Box3f getBoundingBox() const
bool findClosestPoint(const Vector3f &point, MeshProjectionResult &res, float maxDistSq=FLT_MAX, const FaceBitSet *region=nullptr, const AffineXf3f *xf=nullptr) const
Definition MRMesh/MRMesh.h:445
MRMESH_API void pack(FaceMap *outFmap=nullptr, VertMap *outVmap=nullptr, WholeEdgeMap *outEmap=nullptr, bool rearrangeTriangles=false)
Triangle3f getLeftTriPoints(EdgeId e) const
returns three points of left face of e
Definition MRMesh/MRMesh.h:89
MRMESH_API bool projectPoint(const Vector3f &point, PointOnFace &res, float maxDistSq=FLT_MAX, const FaceBitSet *region=nullptr, const AffineXf3f *xf=nullptr) const
MRMESH_API float circumcircleDiameterSq(FaceId f) const
returns squared circumcircle diameter of given mesh triangle
MRMESH_API EdgeId splitEdge(EdgeId e, const Vector3f &newVertPos, FaceBitSet *region=nullptr, FaceHashMap *new2Old=nullptr)
Vector3f leftNormal(EdgeId e) const
computes triangular face normal from its vertices
Definition MRMesh/MRMesh.h:193
MRMESH_API void getLeftTriPoints(EdgeId e, Vector3f &v0, Vector3f &v1, Vector3f &v2) const
returns three points of left face of e
MRMESH_API Vector3f findCenterFromBBox() const
computes bounding box and returns its center
bool isOutside(const Vector3f &pt, float windingNumberThreshold=0.5f, float beta=2) const
Definition MRMesh/MRMesh.h:250
LineSegm3f edgeSegment(EdgeId e) const
returns line segment of given edge
Definition MRMesh/MRMesh.h:71
MRMESH_API float signedDistance(const Vector3f &pt) const
std::optional< MeshProjectionResult > findClosestPoint(const Vector3f &point, float maxDistSq=FLT_MAX, const FaceBitSet *region=nullptr, const AffineXf3f *xf=nullptr) const
Definition MRMesh/MRMesh.h:454
const AABBTreePoints * getAABBTreePointsNotCreate() const
returns cached aabb-tree for points of this mesh, but does not create it if it did not exist
Definition MRMesh/MRMesh.h:466
void addPartByMask(const Mesh &from, const FaceBitSet &fromFaces, const PartMapping &map)
Definition MRMesh/MRMesh.h:388
VertId splitFace(FaceId f, FaceBitSet *region=nullptr, FaceHashMap *new2Old=nullptr)
Definition MRMesh/MRMesh.h:377
MRMESH_API void mirror(const Plane3f &plane)
reflects the mesh from a given plane
MRMESH_API float sumAngles(VertId v, bool *outBoundaryVert=nullptr) const
computes the sum of triangle angles at given vertex; optionally returns whether the vertex is on boun...
MeshTopology topology
Definition MRMesh/MRMesh.h:24
UndirectedEdgeId getClosestEdge(const MeshTriPoint &p) const
returns one of three face edges, closest to given point
Definition MRMesh/MRMesh.h:140
MRMESH_API void transform(const AffineXf3f &xf, const VertBitSet *region=nullptr)
MRMESH_API float dihedralAngle(UndirectedEdgeId e) const
MRMESH_API double volume(const FaceBitSet *region=nullptr) const
MRMESH_API bool projectPoint(const Vector3f &point, MeshProjectionResult &res, float maxDistSq=FLT_MAX, const FaceBitSet *region=nullptr, const AffineXf3f *xf=nullptr) const
MRMESH_API MeshTriPoint toTriPoint(const PointOnFace &p) const
converts face id and 3d point into barycentric representation
VertCoords points
Definition MRMesh/MRMesh.h:25
MRMESH_API float circumcircleDiameter(FaceId f) const
returns circumcircle diameter of given mesh triangle
MRMESH_API Vector3f leftTangent(EdgeId e) const
computes unit vector that is both orthogonal to given edge and to the normal of its left triangle,...
MRMESH_API Vector3f triCenter(FaceId f) const
returns the centroid of given triangle
const Dipoles * getDipolesNotCreate() const
returns cached dipoles of aabb-tree nodes for this mesh, but does not create it if it did not exist
Definition MRMesh/MRMesh.h:472
MRMESH_API float leftCotan(EdgeId e) const
MRMESH_API MeshEdgePoint toEdgePoint(EdgeId e, const Vector3f &p) const
converts edge and 3d point into edge-point representation
float cotan(UndirectedEdgeId ue) const
Definition MRMesh/MRMesh.h:304
MRMESH_API MeshTriPoint toTriPoint(FaceId f, const Vector3f &p) const
converts face id and 3d point into barycentric representation
MRMESH_API float discreteMeanCurvature(VertId v) const
float dblArea(VertId v) const
computes the length of summed directed double areas of all triangles around given vertex
Definition MRMesh/MRMesh.h:202
MRMESH_API void attachEdgeLoopPart(EdgeId first, EdgeId last, const std::vector< Vector3f > &contourPoints)
MRMESH_API Expected< PackMapping > packOptimally(bool preserveAABBTree, ProgressCallback cb)
MRMESH_API EdgeId addSeparateContours(const Contours3f &contours, const AffineXf3f *xf=nullptr)
Vector3f edgePoint(const MeshEdgePoint &ep) const
computes coordinates of point given as edge and relative position on it
Definition MRMesh/MRMesh.h:77
void getTriPoints(FaceId f, Vector3f &v0, Vector3f &v1, Vector3f &v2) const
returns three points of given face
Definition MRMesh/MRMesh.h:92
void getTriPoints(FaceId f, Vector3f(&v)[3]) const
returns three points of given face
Definition MRMesh/MRMesh.h:95
MRMESH_API void addPartByFaceMap(const Mesh &from, const FaceMap &fromFaces, bool flipOrientation=false, const std::vector< EdgePath > &thisContours={}, const std::vector< EdgePath > &fromContours={}, const PartMapping &map={})
fromFaces contains mapping from this-mesh (considering it is empty) to from-mesh
MRMESH_API size_t heapBytes() const
MRMESH_API void invalidateCaches(bool pointsChanged=true)
static MRMESH_API Mesh fromTriMesh(TriMesh &&triMesh, const MeshBuilder::BuildSettings &settings={}, ProgressCallback cb={})
construct mesh from TriMesh representation
MRMESH_API PackMapping packOptimally(bool preserveAABBTree=true)
static MRMESH_API Mesh fromFaceSoup(VertCoords vertexCoordinates, const std::vector< VertId > &verts, const Vector< MeshBuilder::VertSpan, FaceId > &faces, const MeshBuilder::BuildSettings &settings={}, ProgressCallback cb={})
Vector3f normal(FaceId f) const
computes triangular face normal from its vertices
Definition MRMesh/MRMesh.h:196
MRMESH_API float averageEdgeLength() const
computes average length of an edge in this mesh
MRMESH_API const AABBTreePoints & getAABBTreePoints() const
returns cached aabb-tree for points of this mesh, creating it if it did not exist in a thread-safe ma...
MRMESH_API Vector3f leftDirDblArea(EdgeId e) const
computes directed double area of left triangular face of given edge
Vector3f edgeCenter(UndirectedEdgeId e) const
computes the center of given edge
Definition MRMesh/MRMesh.h:80
MRMESH_API void updateCaches(const VertBitSet &changedVerts)
float edgeLength(UndirectedEdgeId e) const
returns Euclidean length of the edge
Definition MRMesh/MRMesh.h:143
void addPart(const Mesh &from, FaceMap *outFmap=nullptr, VertMap *outVmap=nullptr, WholeEdgeMap *outEmap=nullptr, bool rearrangeTriangles=false)
Definition MRMesh/MRMesh.h:383
MRMESH_API float discreteMeanCurvature(UndirectedEdgeId e) const
MRMESH_API Vector3f pseudonormal(const MeshTriPoint &p, const FaceBitSet *region=nullptr) const
MRMESH_API MeshEdgePoint toEdgePoint(VertId v) const
converts vertex into edge-point representation
MRMESH_API Vector3d holeDirArea(EdgeId e) const
MRMESH_API void addMeshPart(const MeshPart &from, bool flipOrientation=false, const std::vector< EdgePath > &thisContours={}, const std::vector< EdgePath > &fromContours={}, const PartMapping &map={})
MRMESH_API VertId addPoint(const Vector3f &pos)
creates new point and assigns given position to it
void getLeftTriPoints(EdgeId e, Vector3f(&v)[3]) const
returns three points of left face of e
Definition MRMesh/MRMesh.h:86
MRMESH_API double holePerimiter(EdgeId e) const
computes the perimeter of the hole specified by one of its edges with no valid left face (left is hol...
double projArea(const Vector3f &dir, const FaceBitSet *fs=nullptr) const
computes the sum of absolute projected area of faces from given region (or whole mesh) as visible if ...
Definition MRMesh/MRMesh.h:176
MRMESH_API float dihedralAngleCos(UndirectedEdgeId e) const
Vector3f dirDblArea(FaceId f) const
computes directed double area for a triangular face from its vertices
Definition MRMesh/MRMesh.h:152
MRMESH_API EdgeId addSeparateEdgeLoop(const std::vector< Vector3f > &contourPoints)
Vector3f normal(VertId v) const
computes normal in a vertex using sum of directed areas of neighboring triangles
Definition MRMesh/MRMesh.h:205
Definition MRBuffer.h:151
Definition MRPartMapping.h:10
Definition MRMesh/MRPointOnFace.h:11
Definition MRTriMesh.h:13