MeshLib C++ Docs
Loading...
Searching...
No Matches
MRTriMath.h
Go to the documentation of this file.
1#pragma once
3
4#include "MRVector2.h"
5#include "MRVector3.h"
6#include <algorithm>
7#include <cassert>
8#include <cmath>
9#include <limits>
10#include <optional>
11
12namespace MR
13{
16
17
20template <typename T>
21[[nodiscard]] T circumcircleDiameterSq( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
22{
23 const auto ab = ( b - a ).lengthSq();
24 const auto ca = ( a - c ).lengthSq();
25 const auto bc = ( c - b ).lengthSq();
26 if ( ab <= 0 )
27 return ca;
28 if ( ca <= 0 )
29 return bc;
30 if ( bc <= 0 )
31 return ab;
32 const auto f = cross( b - a, c - a ).lengthSq();
33 if ( f <= 0 )
34 return std::numeric_limits<T>::infinity();
35 return ab * ca * bc / f;
36}
37
40template <typename T>
41[[nodiscard]] inline T circumcircleDiameter( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
42{
43 return std::sqrt( circumcircleDiameterSq( a, b, c ) );
44}
45
47template <typename T>
48[[nodiscard]] Vector3<T> circumcircleCenter( const Vector3<T> & a, const Vector3<T> & b )
49{
50 const auto xabSq = cross( a, b ).lengthSq();
51 const auto aa = a.lengthSq();
52 const auto bb = b.lengthSq();
53 if ( xabSq <= 0 )
54 {
55 if ( aa <= 0 )
56 return b / T(2);
58 return a / T(2);
59 }
60 const auto ab = dot( a, b );
61 return ( bb * ( aa - ab ) * a + aa * ( bb - ab ) * b ) / ( 2 * xabSq );
62}
63
65template <typename T>
66[[nodiscard]] inline Vector3<T> circumcircleCenter( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
67{
68 return circumcircleCenter( a - c, b - c ) + c;
69}
70
73template <typename T>
74[[nodiscard]] bool circumballCenters( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T radius,
75 Vector3<T> & centerPos,
76 Vector3<T> & centerNeg )
77{
78 const auto rr = sqr( radius );
79 const auto circRadSq = circumcircleDiameterSq( a, b, c ) / T( 4 );
80 if ( rr < circRadSq )
81 return false;
82
83 const auto x = std::sqrt( rr - circRadSq );
84 const auto xn = x * normal( a, b, c );
85 const auto circCenter = circumcircleCenter( a, b, c );
86 centerPos = circCenter + xn;
87 centerNeg = circCenter - xn;
88
89 return true;
90}
91
94template <typename T>
95[[nodiscard]] T minTriangleAngleSin( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
96{
97 const auto ab = ( b - a ).length();
98 const auto ca = ( a - c ).length();
99 const auto bc = ( c - b ).length();
100 if ( ab <= 0 || ca <= 0 || bc <= 0 )
101 return 0;
102 const auto f = cross( b - a, c - a ).length();
103 return f * std::min( { ab, ca, bc } ) / ( ab * ca * bc );
104}
105
106template <typename T>
107[[nodiscard]] T minTriangleAngle( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
108{
109 return std::asin( minTriangleAngleSin( a, b, c ) );
110}
111
114template<typename T>
115[[nodiscard]] T triangleAspectRatio( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
116{
117 const auto bc = ( c - b ).length();
118 const auto ca = ( a - c ).length();
119 const auto ab = ( b - a ).length();
120 auto halfPerimeter = ( bc + ca + ab ) / 2;
121 auto den = 8 * ( halfPerimeter - bc ) * ( halfPerimeter - ca ) * ( halfPerimeter - ab );
122 if ( den <= 0 )
123 return std::numeric_limits<T>::max();
124
125 return bc * ca * ab / den;
126}
127
129template<typename T>
130[[nodiscard]] inline Vector3<T> dirDblArea( const Triangle3<T> & t )
131{
132 return cross( t[1] - t[0], t[2] - t[0] );
133}
134
136template<typename T>
137[[nodiscard]] inline Vector3<T> dirDblArea( const Vector3<T> & q, const Vector3<T> & r )
138{
139 return cross( q, r );
140}
141
143template<typename T>
144[[nodiscard]] inline Vector3<T> dirDblArea( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
145{
146 return cross( q - p, r - p );
147}
148
150template<typename T>
151[[nodiscard]] inline Vector3<T> normal( const Vector3<T> & q, const Vector3<T> & r )
152{
153 return dirDblArea( q, r ).normalized();
154}
155
157template<typename T>
158[[nodiscard]] inline Vector3<T> normal( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
159{
160 return normal( q - p, r - p );
161}
162
164template<typename T>
165[[nodiscard]] inline Vector3<T> normal( const Triangle3<T> & t )
166{
167 return normal( t[1] - t[0], t[2] - t[0] );
168}
169
171template<typename T>
172[[nodiscard]] inline T dblAreaSq( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
173{
174 return dirDblArea( p, q, r ).lengthSq();
175}
176
178template<typename T>
179[[nodiscard]] inline T dblArea( const Triangle3<T> & t )
180{
181 return dirDblArea( t ).length();
182}
183
185template<typename T>
186[[nodiscard]] inline T dblArea( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
187{
188 return dirDblArea( p, q, r ).length();
189}
190
192template<typename T>
193[[nodiscard]] inline T area( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
194{
195 return dblArea( p, q, r ) / 2;
196}
197
199template<typename T>
200[[nodiscard]] inline T dblArea( const Vector2<T> & p, const Vector2<T> & q, const Vector2<T> & r )
201{
202 return std::abs( cross( q - p, r - p ) );
203}
204
206template<typename T>
207[[nodiscard]] inline T area( const Vector2<T> & p, const Vector2<T> & q, const Vector2<T> & r )
208{
209 return dblArea( p, q, r ) / 2;
210}
211
213template <typename T>
214[[nodiscard]] Triangle3<T> makeDegenerate( const Triangle3<T> & t )
215{
216 const auto c = ( t[0] + t[1] + t[2] ) / T(3);
217 int longest = 0;
218 T longestSq = 0;
219 for ( int i = 0; i < 3; ++i )
220 {
221 const auto sq = ( t[i] - c ).lengthSq();
222 if ( longestSq >= sq )
223 continue;
224 longest = i;
225 longestSq = sq;
226 }
227 const auto d = ( t[longest] - c ).normalized();
228
230 Triangle3<T> res;
231 for ( int i = 0; i < 3; ++i )
232 res[i] = c + d * dot( d, t[i] - c );
233 return res;
234}
235
238template <typename T>
239[[nodiscard]] Triangle3<T> triangleWithNormal( const Triangle3<T> & t, const Vector3<T> & n )
240{
241 const auto c = ( t[0] + t[1] + t[2] ) / T(3);
242 Triangle3<T> res;
243 for ( int i = 0; i < 3; ++i )
244 res[i] = t[i] - n * dot( n, t[i] - c );
245
246 if ( dot( n, dirDblArea( res ) ) < 0 )
247 res = makeDegenerate( res );
248 return res;
249}
250
256template <typename T>
257[[nodiscard]] T dihedralAngleSin( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm, const Vector3<T>& edgeVec )
258{
259 auto edgeDir = edgeVec.normalized();
260 return dot( edgeDir, cross( leftNorm, rightNorm ) );
261}
262
268template <typename T>
269[[nodiscard]] T dihedralAngleCos( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm )
270{
271 return dot( leftNorm, rightNorm );
272}
273
280template <typename T>
281[[nodiscard]] T dihedralAngle( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm, const Vector3<T>& edgeVec )
282{
283 auto sin = dihedralAngleSin( leftNorm, rightNorm, edgeVec );
284 auto cos = dihedralAngleCos( leftNorm, rightNorm );
285 return std::atan2( sin, cos );
286}
287
291template <typename T>
292[[nodiscard]] std::optional<Vector2<T>> posFromTriEdgeLengths( T a, T b, T c )
293{
294 if ( c == 0 )
295 {
296 if ( a == b )
297 return Vector2<T>{ a, 0 };
298 else
299 return {};
300 }
301 const auto aa = sqr( a );
302 const auto y = ( aa - sqr( b ) + sqr( c ) ) / ( 2 * c );
303 const auto yy = sqr( y );
304 if ( yy > aa )
305 return {};
306 const auto x = std::sqrt( aa - yy );
307 return Vector2<T>{ x, y };
308}
309
313template <typename T>
314[[nodiscard]] std::optional<T> quadrangleOtherDiagonal( T a, T b, T c, T a1, T b1 )
315{
316 const auto p = posFromTriEdgeLengths( a, b, c );
317 if ( !p )
318 return {};
319 auto p1 = posFromTriEdgeLengths( a1, b1, c );
320 if ( !p1 )
321 return {};
322 p1->x = -p1->x;
324 auto y = ( p->x * p1->y - p1->x * p->y ) / ( p->x - p1->x );
325 if ( y < 0 || y > c )
326 return {};
327 return ( *p - *p1 ).length();
328}
329
333template <typename T>
334[[nodiscard]] inline T tanSqOfHalfAngle( T a, T b, T c )
335{
336 const T den = ( a + b + c ) * ( b + c - a );
337 if ( den <= 0 )
338 return std::numeric_limits<T>::infinity();
339 const T num = ( a + c - b ) * ( a + b - c );
340 if ( num <= 0 )
341 return 0;
342 return num / den;
343}
344
347template <typename T>
348[[nodiscard]] inline T cotan( const Triangle3<T> & t, T absMaxVal = std::numeric_limits<T>::max() )
349{
350 auto a = t[0] - t[2];
351 auto b = t[1] - t[2];
352 auto nom = dot( a, b );
353 auto den = cross( a, b ).length();
354 if ( fabs( nom ) >= absMaxVal * den )
355 return absMaxVal * sgn( nom );
356 return nom / den;
357}
358
362template <typename T>
363[[nodiscard]] inline T cotan( T a, T b, T c )
364{
365 const T den = ( a + b + c ) * ( b + c - a );
366 if ( den <= 0 )
367 return -std::numeric_limits<T>::infinity();
368 const T num = ( a + c - b ) * ( a + b - c );
369 if ( num <= 0 )
370 return std::numeric_limits<T>::infinity();
371 const auto tanSq = num / den;
372 return ( 1 - tanSq ) / ( 2 * std::sqrt( tanSq ) );
373}
374
377template <typename T>
378[[nodiscard]] std::optional<Vector3<T>> gradientInTri( const Vector3<T> & b, const Vector3<T> & c, T vb, T vc )
379{
380 const auto bb = dot( b, b );
381 const auto bc = dot( b, c );
382 const auto cc = dot( c, c );
383 const auto det = bb * cc - bc * bc;
384 if ( det <= 0 )
385 return {};
386 const auto kb = ( 1 / det ) * ( cc * vb - bc * vc );
387 const auto kc = ( 1 / det ) * (-bc * vb + bb * vc );
388 return kb * b + kc * c;
389}
390
393template <typename T>
394[[nodiscard]] std::optional<Vector3<T>> gradientInTri( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T va, T vb, T vc )
395{
396 return gradientInTri( b - a, c - a, vb - va, vc - va );
397}
398
401template <typename T>
402[[nodiscard]] std::optional<T> findTriExitPos( const Vector3<T> & b, const Vector3<T> & c, const Vector3<T> & grad )
403{
404 const auto gradSq = grad.lengthSq();
405 if ( gradSq <= 0 )
406 return {};
407 const auto d = c - b;
409 const auto gort = d - ( dot( d, grad ) / gradSq ) * grad;
410 const auto god = dot( gort, d );
411 if ( god <= 0 )
412 return {};
413 const auto gob = -dot( gort, b );
414 if ( gob <= 0 || gob >= god )
415 return {};
416 const auto a = gob / god;
417 assert( a < std::numeric_limits<T>::max() );
418 const auto ip = a * c + ( 1 - a ) * b;
419 if ( dot( grad, ip ) >= 0 )
420 return {};
421 return a;
422}
423
430template <typename T>
431[[nodiscard]] std::optional<Vector3<T>> tangentPlaneNormalToSpheres( const Vector3<T> & b, const Vector3<T> & c, T rb, T rc )
432{
433 auto grad = gradientInTri( b, c, rb, rc );
434 if ( !grad.has_value() )
435 return {};
436 auto gradSq = grad->lengthSq();
437 if ( gradSq >= 1 )
438 return {};
439 return sqrt( 1 - gradSq ) * normal( b, c ) - *grad;
440}
441
448template <typename T>
449[[nodiscard]] std::optional<Plane3<T>> tangentPlaneToSpheres( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T ra, T rb, T rc )
450{
451 if ( auto n = tangentPlaneNormalToSpheres( b - a, c - a, rb - ra, rc - ra ) )
452 return Plane3<T>( *n, dot( *n, a ) + ra );
453 return {};
454}
455
456}
float area(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns the area of given face
Definition MRMeshMath.h:168
MRMESH_API float circumcircleDiameter(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns circumcircle diameter of given mesh triangle
MRMESH_API float dihedralAngle(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
Vector3< T > circumcircleCenter(const Vector3< T > &a, const Vector3< T > &b)
Computes the center of the the triangle's 0AB circumcircle.
Definition MRTriMath.h:48
MRMESH_API float dihedralAngleCos(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
std::optional< T > findTriExitPos(const Vector3< T > &b, const Vector3< T > &c, const Vector3< T > &grad)
Definition MRTriMath.h:402
std::optional< T > quadrangleOtherDiagonal(T a, T b, T c, T a1, T b1)
Definition MRTriMath.h:314
Triangle3< T > makeDegenerate(const Triangle3< T > &t)
make degenerate triangle (all 3 points on a line) that maximally resembles the input one and has the ...
Definition MRTriMath.h:214
T tanSqOfHalfAngle(T a, T b, T c)
Definition MRTriMath.h:334
Triangle3< T > triangleWithNormal(const Triangle3< T > &t, const Vector3< T > &n)
Definition MRTriMath.h:239
Vector3f dirDblArea(const MeshTopology &topology, const VertCoords &points, FaceId f)
computes directed double area for a triangular face from its vertices
Definition MRMeshMath.h:153
MRMESH_API float triangleAspectRatio(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns aspect ratio of given mesh triangle equal to the ratio of the circum-radius to twice its in-r...
std::optional< Vector3< T > > gradientInTri(const Vector3< T > &b, const Vector3< T > &c, T vb, T vc)
Definition MRTriMath.h:378
std::optional< Plane3< T > > tangentPlaneToSpheres(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c, T ra, T rb, T rc)
Definition MRTriMath.h:449
MRMESH_API float circumcircleDiameterSq(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns squared circumcircle diameter of given mesh triangle
bool circumballCenters(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c, T radius, Vector3< T > &centerPos, Vector3< T > &centerNeg)
ball's center from the negative side of triangle
Definition MRTriMath.h:74
std::optional< Vector3< T > > tangentPlaneNormalToSpheres(const Vector3< T > &b, const Vector3< T > &c, T rb, T rc)
Definition MRTriMath.h:431
Vector3 normalized() const MR_REQUIRES_IF_SUPPORTED(std
Definition MRVector3.h:77
float dblArea(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns twice the area of given face
Definition MRMeshMath.h:162
T minTriangleAngle(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c)
Definition MRTriMath.h:107
std::optional< Vector2< T > > posFromTriEdgeLengths(T a, T b, T c)
Definition MRTriMath.h:292
float cotan(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId ue)
Definition MRMeshMath.h:315
T lengthSq() const
Definition MRVector3.h:68
T dblAreaSq(const Vector3< T > &p, const Vector3< T > &q, const Vector3< T > &r)
computes the square of double area of given triangle
Definition MRTriMath.h:172
length
Definition MRObjectDimensionsEnum.h:17
MRMESH_API float dihedralAngleSin(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
T minTriangleAngleSin(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c)
Definition MRTriMath.h:95
only for bindings generation
Definition MRCameraOrientationPlugin.h:8
Definition MRPlane3.h:21
Definition MRVector2.h:29
Definition MRVector3.h:33