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 "MRPch/MRBindingMacros.h"
7#include <algorithm>
8#include <cassert>
9#include <cmath>
10#include <limits>
11#include <optional>
12
13namespace MR
14{
17
18
21template <typename T>
22[[nodiscard]] T circumcircleDiameterSq( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
23{
24 const auto ab = ( b - a ).lengthSq();
25 const auto ca = ( a - c ).lengthSq();
26 const auto bc = ( c - b ).lengthSq();
27 if ( ab <= 0 )
28 return ca;
29 if ( ca <= 0 )
30 return bc;
31 if ( bc <= 0 )
32 return ab;
33 const auto f = cross( b - a, c - a ).lengthSq();
34 if ( f <= 0 )
35 return std::numeric_limits<T>::infinity();
36 return ab * ca * bc / f;
37}
38
39MR_BIND_TEMPLATE( float circumcircleDiameterSq( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
40MR_BIND_TEMPLATE( double circumcircleDiameterSq( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
41
44template <typename T>
45[[nodiscard]] inline T circumcircleDiameter( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
46{
47 return std::sqrt( circumcircleDiameterSq( a, b, c ) );
48}
49
50MR_BIND_TEMPLATE( float circumcircleDiameter( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
51MR_BIND_TEMPLATE( double circumcircleDiameter( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
52
56template <typename T>
57[[nodiscard]] T mincircleDiameterSq( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
58{
59 const auto ab = ( b - a ).lengthSq();
60 const auto ca = ( a - c ).lengthSq();
61 const auto bc = ( c - b ).lengthSq();
62 if ( ca >= bc + ab )
63 return ca;
64 if ( bc >= ab + ca )
65 return bc;
66 if ( ab >= ca + bc )
67 return ab;
68 const auto f = cross( b - a, c - a ).lengthSq();
69 assert( f > 0 );
70 return ab * ca * bc / f;
71}
72
73MR_BIND_TEMPLATE( float mincircleDiameterSq( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
74MR_BIND_TEMPLATE( double mincircleDiameterSq( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
75
77template <typename T>
78[[nodiscard]] Vector3<T> circumcircleCenter( const Vector3<T> & a, const Vector3<T> & b )
79{
80 const auto xabSq = cross( a, b ).lengthSq();
81 const auto aa = a.lengthSq();
82 const auto bb = b.lengthSq();
83 if ( xabSq <= 0 )
84 {
85 if ( aa <= 0 )
86 return b / T(2);
88 return a / T(2);
89 }
90 const auto ab = dot( a, b );
91 return ( bb * ( aa - ab ) * a + aa * ( bb - ab ) * b ) / ( 2 * xabSq );
92}
93
94MR_BIND_TEMPLATE( Vector3f circumcircleCenter( const Vector3f & a, const Vector3f & b ) );
95MR_BIND_TEMPLATE( Vector3d circumcircleCenter( const Vector3d & a, const Vector3d & b ) );
96
98template <typename T>
99[[nodiscard]] inline Vector3<T> circumcircleCenter( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
100{
101 return circumcircleCenter( a - c, b - c ) + c;
102}
103
104MR_BIND_TEMPLATE( Vector3f circumcircleCenter( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
105MR_BIND_TEMPLATE( Vector3d circumcircleCenter( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
106
109template <typename T>
110[[nodiscard]] bool circumballCenters( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T radius,
111 Vector3<T> & centerPos,
112 Vector3<T> & centerNeg )
113{
114 const auto rr = sqr( radius );
115 const auto circRadSq = circumcircleDiameterSq( a, b, c ) / T( 4 );
116 if ( rr < circRadSq )
117 return false;
118
119 const auto x = std::sqrt( rr - circRadSq );
120 const auto xn = x * normal( a, b, c );
121 const auto circCenter = circumcircleCenter( a, b, c );
122 centerPos = circCenter + xn;
123 centerNeg = circCenter - xn;
124
125 return true;
126}
127
128MR_BIND_TEMPLATE( bool circumballCenters( const Vector3f & a, const Vector3f & b, const Vector3f & c, float radius, Vector3f & centerPos, Vector3f & centerNeg ) );
129MR_BIND_TEMPLATE( bool circumballCenters( const Vector3d & a, const Vector3d & b, const Vector3d & c, double radius, Vector3d & centerPos, Vector3d & centerNeg ) );
130
133template <typename T>
134[[nodiscard]] T minTriangleAngleSin( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
135{
136 const auto ab = ( b - a ).length();
137 const auto ca = ( a - c ).length();
138 const auto bc = ( c - b ).length();
139 if ( ab <= 0 || ca <= 0 || bc <= 0 )
140 return 0;
141 const auto f = cross( b - a, c - a ).length();
142 return f * std::min( { ab, ca, bc } ) / ( ab * ca * bc );
143}
144
145MR_BIND_TEMPLATE( float minTriangleAngleSin( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
146MR_BIND_TEMPLATE( double minTriangleAngleSin( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
147
148template <typename T>
149[[nodiscard]] T minTriangleAngle( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
150{
151 return std::asin( minTriangleAngleSin( a, b, c ) );
152}
153
154MR_BIND_TEMPLATE( float minTriangleAngle( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
155MR_BIND_TEMPLATE( double minTriangleAngle( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
156
159template<typename T>
160[[nodiscard]] T triangleAspectRatio( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c )
161{
162 const auto bc = ( c - b ).length();
163 const auto ca = ( a - c ).length();
164 const auto ab = ( b - a ).length();
165 auto halfPerimeter = ( bc + ca + ab ) / 2;
166 auto den = 8 * ( halfPerimeter - bc ) * ( halfPerimeter - ca ) * ( halfPerimeter - ab );
167 if ( den <= 0 )
168 return std::numeric_limits<T>::max();
169
170 return bc * ca * ab / den;
171}
172
173MR_BIND_TEMPLATE( float triangleAspectRatio( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
174MR_BIND_TEMPLATE( double triangleAspectRatio( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
175
177template<typename T>
178[[nodiscard]] inline Vector3<T> dirDblArea( const Triangle3<T> & t )
179{
180 return cross( t[1] - t[0], t[2] - t[0] );
181}
182
183MR_BIND_TEMPLATE( Vector3f dirDblArea( const Triangle3f & t ) );
184MR_BIND_TEMPLATE( Vector3d dirDblArea( const Triangle3d & t ) );
185
187template<typename T>
188[[nodiscard]] inline Vector3<T> dirDblArea( const Vector3<T> & q, const Vector3<T> & r )
189{
190 return cross( q, r );
191}
192
193MR_BIND_TEMPLATE( Vector3f dirDblArea( const Vector3f & a, const Vector3f & b ) );
194MR_BIND_TEMPLATE( Vector3d dirDblArea( const Vector3d & a, const Vector3d & b ) );
195
197template<typename T>
198[[nodiscard]] inline Vector3<T> dirDblArea( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
199{
200 return cross( q - p, r - p );
201}
202
203MR_BIND_TEMPLATE( Vector3f dirDblArea( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
204MR_BIND_TEMPLATE( Vector3d dirDblArea( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
205
207template<typename T>
208[[nodiscard]] inline Vector3<T> normal( const Vector3<T> & q, const Vector3<T> & r )
209{
210 return dirDblArea( q, r ).normalized();
211}
212
213MR_BIND_TEMPLATE( Vector3f normal( const Vector3f & a, const Vector3f & b ) );
214MR_BIND_TEMPLATE( Vector3d normal( const Vector3d & a, const Vector3d & b ) );
215
217template<typename T>
218[[nodiscard]] inline Vector3<T> normal( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
219{
220 return normal( q - p, r - p );
221}
222
223MR_BIND_TEMPLATE( Vector3f normal( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
224MR_BIND_TEMPLATE( Vector3d normal( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
225
227template<typename T>
228[[nodiscard]] inline Vector3<T> normal( const Triangle3<T> & t )
229{
230 return normal( t[1] - t[0], t[2] - t[0] );
231}
232
233MR_BIND_TEMPLATE( Vector3f normal( const Triangle3f & t ) );
234MR_BIND_TEMPLATE( Vector3d normal( const Triangle3d & t ) );
235
237template<typename T>
238[[nodiscard]] inline T dblAreaSq( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
239{
240 return dirDblArea( p, q, r ).lengthSq();
241}
242
243MR_BIND_TEMPLATE( float dblAreaSq( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
244MR_BIND_TEMPLATE( double dblAreaSq( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
245
247template<typename T>
248[[nodiscard]] inline T dblArea( const Triangle3<T> & t )
249{
250 return dirDblArea( t ).length();
251}
252
253MR_BIND_TEMPLATE( float dblArea( const Triangle3f & t ) );
254MR_BIND_TEMPLATE( double dblArea( const Triangle3d & t ) );
255
257template<typename T>
258[[nodiscard]] inline T dblArea( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
259{
260 return dirDblArea( p, q, r ).length();
261}
262
263MR_BIND_TEMPLATE( float dblArea( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
264MR_BIND_TEMPLATE( double dblArea( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
265
267template<typename T>
268[[nodiscard]] inline T area( const Vector3<T> & p, const Vector3<T> & q, const Vector3<T> & r )
269{
270 return dblArea( p, q, r ) / 2;
271}
272
273MR_BIND_TEMPLATE( float area( const Vector3f & a, const Vector3f & b, const Vector3f & c ) );
274MR_BIND_TEMPLATE( double area( const Vector3d & a, const Vector3d & b, const Vector3d & c ) );
275
277template<typename T>
278[[nodiscard]] inline T dblArea( const Vector2<T> & p, const Vector2<T> & q, const Vector2<T> & r )
279{
280 return std::abs( cross( q - p, r - p ) );
281}
282
283MR_BIND_TEMPLATE( float dblArea( const Vector2f & a, const Vector2f & b, const Vector2f & c ) );
284MR_BIND_TEMPLATE( double dblArea( const Vector2d & a, const Vector2d & b, const Vector2d & c ) );
285
287template<typename T>
288[[nodiscard]] inline T area( const Vector2<T> & p, const Vector2<T> & q, const Vector2<T> & r )
289{
290 return dblArea( p, q, r ) / 2;
291}
292
293MR_BIND_TEMPLATE( float area( const Vector2f & a, const Vector2f & b, const Vector2f & c ) );
294MR_BIND_TEMPLATE( double area( const Vector2d & a, const Vector2d & b, const Vector2d & c ) );
295
297template <typename T>
298[[nodiscard]] Triangle3<T> makeDegenerate( const Triangle3<T> & t )
299{
300 const auto c = ( t[0] + t[1] + t[2] ) / T(3);
301 int longest = 0;
302 T longestSq = 0;
303 for ( int i = 0; i < 3; ++i )
304 {
305 const auto sq = ( t[i] - c ).lengthSq();
306 if ( longestSq >= sq )
307 continue;
308 longest = i;
309 longestSq = sq;
310 }
311 const auto d = ( t[longest] - c ).normalized();
312
314 Triangle3<T> res;
315 for ( int i = 0; i < 3; ++i )
316 res[i] = c + d * dot( d, t[i] - c );
317 return res;
318}
319
322
325template <typename T>
326[[nodiscard]] Triangle3<T> triangleWithNormal( const Triangle3<T> & t, const Vector3<T> & n )
327{
328 const auto c = ( t[0] + t[1] + t[2] ) / T(3);
329 Triangle3<T> res;
330 for ( int i = 0; i < 3; ++i )
331 res[i] = t[i] - n * dot( n, t[i] - c );
332
333 if ( dot( n, dirDblArea( res ) ) < 0 )
334 res = makeDegenerate( res );
335 return res;
336}
337
338MR_BIND_TEMPLATE( Triangle3f triangleWithNormal( const Triangle3f & t, const Vector3f & n ) );
339MR_BIND_TEMPLATE( Triangle3d triangleWithNormal( const Triangle3d & t, const Vector3d & n ) );
340
346template <typename T>
347[[nodiscard]] T dihedralAngleSin( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm, const Vector3<T>& edgeVec )
348{
349 auto edgeDir = edgeVec.normalized();
350 return dot( edgeDir, cross( leftNorm, rightNorm ) );
351}
352
353MR_BIND_TEMPLATE( float dihedralAngleSin( const Vector3f&, const Vector3f&, const Vector3f& ) );
354MR_BIND_TEMPLATE( double dihedralAngleSin( const Vector3d&, const Vector3d&, const Vector3d& ) );
355
361template <typename T>
362[[nodiscard]] T dihedralAngleCos( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm )
363{
364 return dot( leftNorm, rightNorm );
365}
366
367MR_BIND_TEMPLATE( float dihedralAngleCos( const Vector3f&, const Vector3f& ) );
368MR_BIND_TEMPLATE( double dihedralAngleCos( const Vector3d&, const Vector3d& ) );
369
376template <typename T>
377[[nodiscard]] T dihedralAngle( const Vector3<T>& leftNorm, const Vector3<T>& rightNorm, const Vector3<T>& edgeVec )
378{
379 auto sin = dihedralAngleSin( leftNorm, rightNorm, edgeVec );
380 auto cos = dihedralAngleCos( leftNorm, rightNorm );
381 return std::atan2( sin, cos );
382}
383
384MR_BIND_TEMPLATE( float dihedralAngle( const Vector3f&, const Vector3f&, const Vector3f& ) );
385MR_BIND_TEMPLATE( double dihedralAngle( const Vector3d&, const Vector3d&, const Vector3d& ) );
386
390template <typename T>
391[[nodiscard]] std::optional<Vector2<T>> posFromTriEdgeLengths( T a, T b, T c )
392{
393 if ( c == 0 )
394 {
395 if ( a == b )
396 return Vector2<T>{ a, 0 };
397 else
398 return {};
399 }
400 const auto aa = sqr( a );
401 const auto y = ( aa - sqr( b ) + sqr( c ) ) / ( 2 * c );
402 const auto yy = sqr( y );
403 if ( yy > aa )
404 return {};
405 const auto x = std::sqrt( aa - yy );
406 return Vector2<T>{ x, y };
407}
408
409MR_BIND_TEMPLATE( std::optional<Vector2f> posFromTriEdgeLengths( float, float, float ) );
410MR_BIND_TEMPLATE( std::optional<Vector2d> posFromTriEdgeLengths( double, double, double ) );
411
415template <typename T>
416[[nodiscard]] std::optional<T> quadrangleOtherDiagonal( T a, T b, T c, T a1, T b1 )
417{
418 const auto p = posFromTriEdgeLengths( a, b, c );
419 if ( !p )
420 return {};
421 auto p1 = posFromTriEdgeLengths( a1, b1, c );
422 if ( !p1 )
423 return {};
424 p1->x = -p1->x;
426 auto y = ( p->x * p1->y - p1->x * p->y ) / ( p->x - p1->x );
427 if ( y < 0 || y > c )
428 return {};
429 return ( *p - *p1 ).length();
430}
431
432MR_BIND_TEMPLATE( std::optional<float > quadrangleOtherDiagonal( float, float, float, float, float ) );
433MR_BIND_TEMPLATE( std::optional<double> quadrangleOtherDiagonal( double, double, double, double, double ) );
434
438template <typename T>
439[[nodiscard]] inline T tanSqOfHalfAngle( T a, T b, T c )
440{
441 const T den = ( a + b + c ) * ( b + c - a );
442 if ( den <= 0 )
443 return std::numeric_limits<T>::infinity();
444 const T num = ( a + c - b ) * ( a + b - c );
445 if ( num <= 0 )
446 return 0;
447 return num / den;
448}
449
450MR_BIND_TEMPLATE( float tanSqOfHalfAngle( float, float, float ) );
451MR_BIND_TEMPLATE( double tanSqOfHalfAngle( double, double, double ) );
452
455template <typename T>
456[[nodiscard]] inline T cotan( const Triangle3<T> & t, T absMaxVal = std::numeric_limits<T>::max() )
457{
458 auto a = t[0] - t[2];
459 auto b = t[1] - t[2];
460 auto nom = dot( a, b );
461 auto den = cross( a, b ).length();
462 if ( fabs( nom ) >= absMaxVal * den )
463 return absMaxVal * sgn( nom );
464 return nom / den;
465}
466
467MR_BIND_TEMPLATE( float cotan( const Triangle3f&, float ) );
468MR_BIND_TEMPLATE( double cotan( const Triangle3d&, double ) );
469
473template <typename T>
474[[nodiscard]] inline T cotan( T a, T b, T c )
475{
476 const T den = ( a + b + c ) * ( b + c - a );
477 if ( den <= 0 )
478 return -std::numeric_limits<T>::infinity();
479 const T num = ( a + c - b ) * ( a + b - c );
480 if ( num <= 0 )
481 return std::numeric_limits<T>::infinity();
482 const auto tanSq = num / den;
483 return ( 1 - tanSq ) / ( 2 * std::sqrt( tanSq ) );
484}
485
486MR_BIND_TEMPLATE( float cotan( float, float, float ) );
487MR_BIND_TEMPLATE( double cotan( double, double, double ) );
488
491template <typename T>
492[[nodiscard]] std::optional<Vector3<T>> gradientInTri( const Vector3<T> & b, const Vector3<T> & c, T vb, T vc )
493{
494 const auto bb = dot( b, b );
495 const auto bc = dot( b, c );
496 const auto cc = dot( c, c );
497 const auto det = bb * cc - bc * bc;
498 if ( det <= 0 )
499 return {};
500 const auto kb = ( 1 / det ) * ( cc * vb - bc * vc );
501 const auto kc = ( 1 / det ) * (-bc * vb + bb * vc );
502 return kb * b + kc * c;
503}
504
505MR_BIND_TEMPLATE( std::optional<Vector3f> gradientInTri( const Vector3f&, const Vector3f&, float, float ) );
506MR_BIND_TEMPLATE( std::optional<Vector3d> gradientInTri( const Vector3d&, const Vector3d&, double, double ) );
507
510template <typename T>
511[[nodiscard]] std::optional<Vector3<T>> gradientInTri( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T va, T vb, T vc )
512{
513 return gradientInTri( b - a, c - a, vb - va, vc - va );
514}
515
516MR_BIND_TEMPLATE( std::optional<Vector3f> gradientInTri( const Vector3f&, const Vector3f&, const Vector3f&, float, float, float ) );
517MR_BIND_TEMPLATE( std::optional<Vector3d> gradientInTri( const Vector3d&, const Vector3d&, const Vector3d&, double, double, double ) );
518
521template <typename T>
522[[nodiscard]] std::optional<T> findTriExitPos( const Vector3<T> & b, const Vector3<T> & c, const Vector3<T> & grad )
523{
524 const auto gradSq = grad.lengthSq();
525 if ( gradSq <= 0 )
526 return {};
527 const auto d = c - b;
529 const auto gort = d - ( dot( d, grad ) / gradSq ) * grad;
530 const auto god = dot( gort, d );
531 if ( god <= 0 )
532 return {};
533 const auto gob = -dot( gort, b );
534 if ( gob <= 0 || gob >= god )
535 return {};
536 const auto a = gob / god;
537 assert( a < std::numeric_limits<T>::max() );
538 const auto ip = a * c + ( 1 - a ) * b;
539 if ( dot( grad, ip ) >= 0 )
540 return {};
541 return a;
542}
543
544MR_BIND_TEMPLATE( std::optional<float > findTriExitPos( const Vector3f&, const Vector3f&, const Vector3f& ) );
545MR_BIND_TEMPLATE( std::optional<double> findTriExitPos( const Vector3d&, const Vector3d&, const Vector3d& ) );
546
553template <typename T>
554[[nodiscard]] std::optional<Vector3<T>> tangentPlaneNormalToSpheres( const Vector3<T> & b, const Vector3<T> & c, T rb, T rc )
555{
556 auto grad = gradientInTri( b, c, rb, rc );
557 if ( !grad.has_value() )
558 return {};
559 auto gradSq = grad->lengthSq();
560 if ( gradSq >= 1 )
561 return {};
562 return sqrt( 1 - gradSq ) * normal( b, c ) - *grad;
563}
564
565MR_BIND_TEMPLATE( std::optional<Vector3f> tangentPlaneNormalToSpheres( const Vector3f&, const Vector3f&, float, float ) );
566MR_BIND_TEMPLATE( std::optional<Vector3d> tangentPlaneNormalToSpheres( const Vector3d&, const Vector3d&, double, double ) );
567
574template <typename T>
575[[nodiscard]] std::optional<Plane3<T>> tangentPlaneToSpheres( const Vector3<T> & a, const Vector3<T> & b, const Vector3<T> & c, T ra, T rb, T rc )
576{
577 if ( auto n = tangentPlaneNormalToSpheres( b - a, c - a, rb - ra, rc - ra ) )
578 return Plane3<T>( *n, dot( *n, a ) + ra );
579 return {};
580}
581
582MR_BIND_TEMPLATE( std::optional<Plane3f> tangentPlaneToSpheres( const Vector3f&, const Vector3f&, const Vector3f&, float, float, float ) );
583MR_BIND_TEMPLATE( std::optional<Plane3d> tangentPlaneToSpheres( const Vector3d&, const Vector3d&, const Vector3d&, double, double, double ) );
584
585}
Plane3
Definition MRMeshFwd.h:394
float area(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns the area of given face
Definition MRMeshMath.h:168
constexpr T sqr(T x) noexcept
squared value
Definition MRMeshFwd.h:768
float dihedralAngleSin(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
auto dot(const Matrix2< T > &a, const Matrix2< T > &b) -> decltype(dot(a.x, b.x))
double-dot product: x = a : b
Definition MRMatrix2.h:142
Vector3< T > circumcircleCenter(const Vector3< T > &a, const Vector3< T > &b)
Computes the center of the the triangle's 0AB circumcircle.
Definition MRTriMath.h:78
std::optional< T > findTriExitPos(const Vector3< T > &b, const Vector3< T > &c, const Vector3< T > &grad)
Definition MRTriMath.h:522
std::optional< T > quadrangleOtherDiagonal(T a, T b, T c, T a1, T b1)
Definition MRTriMath.h:416
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:298
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...
T tanSqOfHalfAngle(T a, T b, T c)
Definition MRTriMath.h:439
Triangle3< T > triangleWithNormal(const Triangle3< T > &t, const Vector3< T > &n)
Definition MRTriMath.h:326
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
float dihedralAngleCos(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
constexpr int sgn(T x) noexcept
sign of given value in { -1, 0, 1 }
Definition MRMeshFwd.h:772
std::optional< Vector3< T > > gradientInTri(const Vector3< T > &b, const Vector3< T > &c, T vb, T vc)
Definition MRTriMath.h:492
Triangle3< double > Triangle3d
Definition MRMeshFwd.h:471
float circumcircleDiameter(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns circumcircle diameter of given mesh triangle
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:575
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:110
std::optional< Vector3< T > > tangentPlaneNormalToSpheres(const Vector3< T > &b, const Vector3< T > &c, T rb, T rc)
Definition MRTriMath.h:554
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
MR_BIND_TEMPLATE(std::pair< Vector3f, TriPointf > closestPointInTriangle(const Vector3f &p, const Vector3f &a, const Vector3f &b, const Vector3f &c))
T minTriangleAngle(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c)
Definition MRTriMath.h:149
Triangle3< float > Triangle3f
Definition MRMeshFwd.h:470
std::optional< Vector2< T > > posFromTriEdgeLengths(T a, T b, T c)
Definition MRTriMath.h:391
float circumcircleDiameterSq(const MeshTopology &topology, const VertCoords &points, FaceId f)
returns squared circumcircle diameter of given mesh triangle
float cotan(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId ue)
Definition MRMeshMath.h:315
T lengthSq() const
Definition MRVector3.h:68
float dihedralAngle(const MeshTopology &topology, const VertCoords &points, UndirectedEdgeId e)
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:238
length
Definition MRObjectDimensionsEnum.h:17
@ normal
Definition MRUnits.h:45
T mincircleDiameterSq(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c)
Definition MRTriMath.h:57
T minTriangleAngleSin(const Vector3< T > &a, const Vector3< T > &b, const Vector3< T > &c)
Definition MRTriMath.h:134
only for bindings generation
Definition MRCameraOrientationPlugin.h:8
Definition MRVector2.h:29
T cross(const Vector2< T > &a, const Vector2< T > &b)
cross product
Definition MRVector2.h:160
Definition MRVector3.h:33