VTK  9.4.20251007
vtkMath.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright 2011 Sandia Corporation
3// SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
23
24#ifndef vtkMath_h
25#define vtkMath_h
26
27#include "vtkCommonCoreModule.h" // For export macro
28#include "vtkMathPrivate.hxx" // For Matrix meta-class helpers
29#include "vtkMatrixUtilities.h" // For Matrix wrapping / mapping
30#include "vtkObject.h"
31#include "vtkSmartPointer.h" // For vtkSmartPointer.
32#include "vtkTypeTraits.h" // For type traits
33
34#include "vtkMathConfigure.h" // For <cmath> and VTK_HAS_ISNAN etc.
35
36#include <algorithm> // for std::clamp
37#include <cassert> // assert() in inline implementations.
38#include <type_traits> // for type_traits
39
40#ifndef DBL_MIN
41#define VTK_DBL_MIN 2.2250738585072014e-308
42#else // DBL_MIN
43#define VTK_DBL_MIN DBL_MIN
44#endif // DBL_MIN
45
46#ifndef DBL_EPSILON
47#define VTK_DBL_EPSILON 2.2204460492503131e-16
48#else // DBL_EPSILON
49#define VTK_DBL_EPSILON DBL_EPSILON
50#endif // DBL_EPSILON
51
52#ifndef VTK_DBL_EPSILON
53#ifndef DBL_EPSILON
54#define VTK_DBL_EPSILON 2.2204460492503131e-16
55#else // DBL_EPSILON
56#define VTK_DBL_EPSILON DBL_EPSILON
57#endif // DBL_EPSILON
58#endif // VTK_DBL_EPSILON
59
60VTK_ABI_NAMESPACE_BEGIN
61class vtkDataArray;
62class vtkPoints;
63class vtkMathInternal;
66VTK_ABI_NAMESPACE_END
67
68namespace vtk_detail
69{
70VTK_ABI_NAMESPACE_BEGIN
71// forward declaration
72template <typename OutT>
73void RoundDoubleToIntegralIfNecessary(double val, OutT* ret);
74VTK_ABI_NAMESPACE_END
75} // end namespace vtk_detail
76
77VTK_ABI_NAMESPACE_BEGIN
78class VTKCOMMONCORE_EXPORT vtkMath : public vtkObject
79{
80public:
81 static vtkMath* New();
82 vtkTypeMacro(vtkMath, vtkObject);
83 void PrintSelf(ostream& os, vtkIndent indent) override;
84
85private:
86 template <class VectorT, class = void>
87 struct VectorImplementsSize : std::false_type
88 {
89 };
90
91 template <class VectorT>
92 struct VectorImplementsSize<VectorT, decltype((void)std::declval<VectorT>().size(), void())>
93 : std::true_type
94 {
95 };
96
100 template <class VectorT>
101 using EnableIfVectorImplementsSize =
102 typename std::enable_if<VectorImplementsSize<VectorT>::value>::type;
103
104public:
113 static constexpr int DYNAMIC_VECTOR_SIZE() { return 0; }
114
118 static constexpr double Pi() { return 3.141592653589793; }
119
121
124 static float RadiansFromDegrees(float degrees);
125 static double RadiansFromDegrees(double degrees);
127
129
132 static float DegreesFromRadians(float radians);
133 static double DegreesFromRadians(double radians);
135
139#if 1
140 static int Round(float f) { return static_cast<int>(f + (f >= 0.0 ? 0.5 : -0.5)); }
141 static int Round(double f) { return static_cast<int>(f + (f >= 0.0 ? 0.5 : -0.5)); }
142#endif
143
148 template <typename OutT>
149 static void RoundDoubleToIntegralIfNecessary(double val, OutT* ret)
150 {
151 // Can't specialize template methods in a template class, so we move the
152 // implementations to a external namespace.
154 }
155
161 static int Floor(double x);
162
168 static int Ceil(double x);
169
175 static int CeilLog2(vtkTypeUInt64 x);
176
181 template <class T>
182 static T Min(const T& a, const T& b);
183
188 template <class T>
189 static T Max(const T& a, const T& b);
190
194 static bool IsPowerOfTwo(vtkTypeUInt64 x);
195
201 static int NearestPowerOfTwo(int x);
202
207 static vtkTypeInt64 Factorial(int N);
208
214 static vtkTypeInt64 Binomial(int m, int n);
215
227 static int* BeginCombination(int m, int n);
228
239 static int NextCombination(int m, int n, int* combination);
240
244 static void FreeCombination(int* combination);
245
261 static void RandomSeed(int s);
262
274 static int GetSeed();
275
289 static double Random();
290
303 static double Random(double min, double max);
304
317 static double Gaussian();
318
331 static double Gaussian(double mean, double std);
332
337 template <class VectorT1, class VectorT2>
338 static void Assign(const VectorT1& a, VectorT2&& b)
339 {
340 b[0] = a[0];
341 b[1] = a[1];
342 b[2] = a[2];
343 }
344
348 static void Assign(const double a[3], double b[3]) { vtkMath::Assign<>(a, b); }
349
353 static void Add(const float a[3], const float b[3], float c[3])
354 {
355 for (int i = 0; i < 3; ++i)
356 {
357 c[i] = a[i] + b[i];
358 }
359 }
360
364 static void Add(const double a[3], const double b[3], double c[3])
365 {
366 for (int i = 0; i < 3; ++i)
367 {
368 c[i] = a[i] + b[i];
369 }
370 }
371
377 template <class VectorT1, class VectorT2, class VectorT3>
378 static void Add(VectorT1&& a, VectorT2&& b, VectorT3& c)
379 {
380 for (int i = 0; i < 3; ++i)
381 {
382 c[i] = a[i] + b[i];
383 }
384 }
385
389 static void Subtract(const float a[3], const float b[3], float c[3])
390 {
391 for (int i = 0; i < 3; ++i)
392 {
393 c[i] = a[i] - b[i];
394 }
395 }
396
400 static void Subtract(const double a[3], const double b[3], double c[3])
401 {
402 for (int i = 0; i < 3; ++i)
403 {
404 c[i] = a[i] - b[i];
405 }
406 }
407
413 template <class VectorT1, class VectorT2, class VectorT3>
414 static void Subtract(const VectorT1& a, const VectorT2& b, VectorT3&& c)
415 {
416 c[0] = a[0] - b[0];
417 c[1] = a[1] - b[1];
418 c[2] = a[2] - b[2];
419 }
420
425 static void MultiplyScalar(float a[3], float s)
426 {
427 for (int i = 0; i < 3; ++i)
428 {
429 a[i] *= s;
430 }
431 }
432
437 static void MultiplyScalar2D(float a[2], float s)
438 {
439 for (int i = 0; i < 2; ++i)
440 {
441 a[i] *= s;
442 }
443 }
444
449 static void MultiplyScalar(double a[3], double s)
450 {
451 for (int i = 0; i < 3; ++i)
452 {
453 a[i] *= s;
454 }
455 }
456
461 static void MultiplyScalar2D(double a[2], double s)
462 {
463 for (int i = 0; i < 2; ++i)
464 {
465 a[i] *= s;
466 }
467 }
468
472 static float Dot(const float a[3], const float b[3])
473 {
474 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
475 }
476
480 static double Dot(const double a[3], const double b[3])
481 {
482 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
483 }
484
500 template <typename ReturnTypeT = double, typename TupleRangeT1, typename TupleRangeT2,
501 typename EnableT = typename std::conditional<!std::is_pointer<TupleRangeT1>::value &&
502 !std::is_array<TupleRangeT1>::value,
503 TupleRangeT1, TupleRangeT2>::type::value_type>
504 static ReturnTypeT Dot(const TupleRangeT1& a, const TupleRangeT2& b)
505 {
506 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
507 }
508
512 static void Outer(const float a[3], const float b[3], float c[3][3])
513 {
514 for (int i = 0; i < 3; ++i)
515 {
516 for (int j = 0; j < 3; ++j)
517 {
518 c[i][j] = a[i] * b[j];
519 }
520 }
521 }
522
526 static void Outer(const double a[3], const double b[3], double c[3][3])
527 {
528 for (int i = 0; i < 3; ++i)
529 {
530 for (int j = 0; j < 3; ++j)
531 {
532 c[i][j] = a[i] * b[j];
533 }
534 }
535 }
536
542 template <class VectorT1, class VectorT2, class VectorT3>
543 static void Cross(VectorT1&& a, VectorT2&& b, VectorT3& c);
544
549 static void Cross(const float a[3], const float b[3], float c[3]);
550
555 static void Cross(const double a[3], const double b[3], double c[3]);
556
558
561 static float Norm(const float* x, int n);
562 static double Norm(const double* x, int n);
564
568 static float Norm(const float v[3]) { return std::sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]); }
569
573 static double Norm(const double v[3])
574 {
575 return std::sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
576 }
577
587 template <typename ReturnTypeT = double, typename TupleRangeT>
588 static ReturnTypeT SquaredNorm(const TupleRangeT& v)
589 {
590 return v[0] * v[0] + v[1] * v[1] + v[2] * v[2];
591 }
592
597 static float Normalize(float v[3]);
598
603 static double Normalize(double v[3]);
604
606
613 static void Perpendiculars(const double v1[3], double v2[3], double v3[3], double theta);
614 static void Perpendiculars(const float v1[3], float v2[3], float v3[3], double theta);
616
618
623 static bool ProjectVector(const float a[3], const float b[3], float projection[3]);
624 static bool ProjectVector(const double a[3], const double b[3], double projection[3]);
626
628
634 static bool ProjectVector2D(const float a[2], const float b[2], float projection[2]);
635 static bool ProjectVector2D(const double a[2], const double b[2], double projection[2]);
637
653 template <typename ReturnTypeT = double, typename TupleRangeT1, typename TupleRangeT2,
654 typename EnableT = typename std::conditional<!std::is_pointer<TupleRangeT1>::value &&
655 !std::is_array<TupleRangeT1>::value,
656 TupleRangeT1, TupleRangeT2>::type::value_type>
657 static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1& p1, const TupleRangeT2& p2);
658
663 static float Distance2BetweenPoints(const float p1[3], const float p2[3]);
664
669 static double Distance2BetweenPoints(const double p1[3], const double p2[3]);
670
674 static double AngleBetweenVectors(const double v1[3], const double v2[3]);
675
680 const double v1[3], const double v2[3], const double vn[3]);
681
686 static double GaussianAmplitude(double variance, double distanceFromMean);
687
692 static double GaussianAmplitude(double mean, double variance, double position);
693
699 static double GaussianWeight(double variance, double distanceFromMean);
700
706 static double GaussianWeight(double mean, double variance, double position);
707
711 static float Dot2D(const float x[2], const float y[2]) { return x[0] * y[0] + x[1] * y[1]; }
712
716 static double Dot2D(const double x[2], const double y[2]) { return x[0] * y[0] + x[1] * y[1]; }
717
721 static void Outer2D(const float x[2], const float y[2], float A[2][2])
722 {
723 for (int i = 0; i < 2; ++i)
724 {
725 for (int j = 0; j < 2; ++j)
726 {
727 A[i][j] = x[i] * y[j];
728 }
729 }
730 }
731
735 static void Outer2D(const double x[2], const double y[2], double A[2][2])
736 {
737 for (int i = 0; i < 2; ++i)
738 {
739 for (int j = 0; j < 2; ++j)
740 {
741 A[i][j] = x[i] * y[j];
742 }
743 }
744 }
745
750 static float Norm2D(const float x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
751
756 static double Norm2D(const double x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
757
762 static float Normalize2D(float v[2]);
763
768 static double Normalize2D(double v[2]);
769
773 static float Determinant2x2(const float c1[2], const float c2[2])
774 {
775 return c1[0] * c2[1] - c2[0] * c1[1];
776 }
777
779
782 static double Determinant2x2(double a, double b, double c, double d) { return a * d - b * c; }
783 static double Determinant2x2(const double c1[2], const double c2[2])
784 {
785 return c1[0] * c2[1] - c2[0] * c1[1];
786 }
787
788
790
793 static void LUFactor3x3(float A[3][3], int index[3]);
794 static void LUFactor3x3(double A[3][3], int index[3]);
796
798
801 static void LUSolve3x3(const float A[3][3], const int index[3], float x[3]);
802 static void LUSolve3x3(const double A[3][3], const int index[3], double x[3]);
804
806
810 static void LinearSolve3x3(const float A[3][3], const float x[3], float y[3]);
811 static void LinearSolve3x3(const double A[3][3], const double x[3], double y[3]);
813
815
818 static void Multiply3x3(const float A[3][3], const float v[3], float u[3]);
819 static void Multiply3x3(const double A[3][3], const double v[3], double u[3]);
821
823
826 static void Multiply3x3(const float A[3][3], const float B[3][3], float C[3][3]);
827 static void Multiply3x3(const double A[3][3], const double B[3][3], double C[3][3]);
829
853 template <int RowsT, int MidDimT, int ColsT,
854 class LayoutT1 = vtkMatrixUtilities::Layout::Identity,
855 class LayoutT2 = vtkMatrixUtilities::Layout::Identity, class MatrixT1, class MatrixT2,
856 class MatrixT3>
857 static void MultiplyMatrix(MatrixT1&& M1, MatrixT2&& M2, MatrixT3&& M3)
858 {
859 vtkMathPrivate::MultiplyMatrix<RowsT, MidDimT, ColsT, LayoutT1, LayoutT2>::Compute(
860 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2), std::forward<MatrixT3>(M3));
861 }
862
883 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
884 class MatrixT, class VectorT1, class VectorT2>
885 static void MultiplyMatrixWithVector(MatrixT&& M, VectorT1&& X, VectorT2&& Y)
886 {
887 vtkMathPrivate::MultiplyMatrix<RowsT, ColsT, 1, LayoutT>::Compute(
888 std::forward<MatrixT>(M), std::forward<VectorT1>(X), std::forward<VectorT2>(Y));
889 }
890
896 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
897 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
898 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
899 {
900 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
901 vtkMatrixUtilities::Layout::Identity,
902 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x),
903 std::forward<VectorT2>(y));
904 }
905
912 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
913 class = typename std::enable_if<SizeT == DYNAMIC_VECTOR_SIZE()>::type,
914 class = EnableIfVectorImplementsSize<VectorT1>>
915 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
916 {
917 ScalarT dot = 0.0;
918 using SizeType = decltype(std::declval<VectorT1>().size());
919 for (SizeType dim = 0; dim < x.size(); ++dim)
920 {
921 dot += x[dim] * y[dim];
922 }
923 return dot;
924 }
925
933 template <int SizeT, class VectorT>
935 VectorT&& x)
936 {
938 return vtkMath::Dot<Scalar, SizeT>(std::forward<VectorT>(x), std::forward<VectorT>(x));
939 }
940
957 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT>
959 MatrixT&& M)
960 {
961 return vtkMathPrivate::Determinant<SizeT, LayoutT>::Compute(std::forward<MatrixT>(M));
962 }
963
979 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT1,
980 class MatrixT2>
981 static void InvertMatrix(MatrixT1&& M1, MatrixT2&& M2)
982 {
983 vtkMathPrivate::InvertMatrix<SizeT, LayoutT>::Compute(
984 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2));
985 }
986
1000 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1001 class MatrixT, class VectorT1, class VectorT2>
1002 static void LinearSolve(MatrixT&& M, VectorT1&& x, VectorT2&& y)
1003 {
1004 vtkMathPrivate::LinearSolve<RowsT, ColsT, LayoutT>::Compute(
1005 std::forward<MatrixT>(M), std::forward<VectorT1>(x), std::forward<VectorT2>(y));
1006 }
1007
1022 template <class ScalarT, int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1023 class VectorT1, class MatrixT, class VectorT2,
1024 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
1025 static ScalarT Dot(VectorT1&& x, MatrixT&& M, VectorT2&& y)
1026 {
1027 ScalarT tmp[SizeT];
1028 vtkMathPrivate::MultiplyMatrix<SizeT, SizeT, 1, LayoutT>::Compute(
1029 std::forward<MatrixT>(M), std::forward<VectorT2>(y), tmp);
1030 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
1031 vtkMatrixUtilities::Layout::Identity,
1032 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x), tmp);
1033 }
1034
1040 static void MultiplyMatrix(const double* const* A, const double* const* B, unsigned int rowA,
1041 unsigned int colA, unsigned int rowB, unsigned int colB, double** C);
1042
1044
1048 static void Transpose3x3(const float A[3][3], float AT[3][3]);
1049 static void Transpose3x3(const double A[3][3], double AT[3][3]);
1051
1053
1057 static void Invert3x3(const float A[3][3], float AI[3][3]);
1058 static void Invert3x3(const double A[3][3], double AI[3][3]);
1060
1062
1065 static void Identity3x3(float A[3][3]);
1066 static void Identity3x3(double A[3][3]);
1068
1070
1073 static double Determinant3x3(const float A[3][3]);
1074 static double Determinant3x3(const double A[3][3]);
1076
1080 static float Determinant3x3(const float c1[3], const float c2[3], const float c3[3]);
1081
1085 static double Determinant3x3(const double c1[3], const double c2[3], const double c3[3]);
1086
1093 static double Determinant3x3(double a1, double a2, double a3, double b1, double b2, double b3,
1094 double c1, double c2, double c3);
1095
1097
1104 static void QuaternionToMatrix3x3(const float quat[4], float A[3][3]);
1105 static void QuaternionToMatrix3x3(const double quat[4], double A[3][3]);
1106 template <class QuaternionT, class MatrixT,
1107 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1108 static void QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A);
1110
1112
1121 static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4]);
1122 static void Matrix3x3ToQuaternion(const double A[3][3], double quat[4]);
1123 template <class MatrixT, class QuaternionT,
1124 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1125 static void Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q);
1127
1129
1135 static void MultiplyQuaternion(const float q1[4], const float q2[4], float q[4]);
1136 static void MultiplyQuaternion(const double q1[4], const double q2[4], double q[4]);
1138
1140
1144 static void RotateVectorByNormalizedQuaternion(const float v[3], const float q[4], float r[3]);
1145 static void RotateVectorByNormalizedQuaternion(const double v[3], const double q[4], double r[3]);
1147
1149
1153 static void RotateVectorByWXYZ(const float v[3], const float q[4], float r[3]);
1154 static void RotateVectorByWXYZ(const double v[3], const double q[4], double r[3]);
1156
1158
1163 static void Orthogonalize3x3(const float A[3][3], float B[3][3]);
1164 static void Orthogonalize3x3(const double A[3][3], double B[3][3]);
1166
1168
1174 static void Diagonalize3x3(const float A[3][3], float w[3], float V[3][3]);
1175 static void Diagonalize3x3(const double A[3][3], double w[3], double V[3][3]);
1177
1179
1189 const float A[3][3], float U[3][3], float w[3], float VT[3][3]);
1191 const double A[3][3], double U[3][3], double w[3], double VT[3][3]);
1193
1202 double a00, double a01, double a10, double a11, double b0, double b1, double& x0, double& x1);
1203
1212 static vtkTypeBool SolveLinearSystem(double** A, double* x, int size);
1213
1220 static vtkTypeBool InvertMatrix(double** A, double** AI, int size);
1221
1228 double** A, double** AI, int size, int* tmp1Size, double* tmp2Size);
1229
1252 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size);
1253
1259 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size, double* tmpSize);
1260
1269 static void LUSolveLinearSystem(double** A, int* index, double* x, int size);
1270
1279 static double EstimateMatrixCondition(const double* const* A, int size);
1280
1282
1290 static vtkTypeBool Jacobi(float** a, float* w, float** v);
1291 static vtkTypeBool Jacobi(double** a, double* w, double** v);
1293
1295
1304 static vtkTypeBool JacobiN(float** a, int n, float* w, float** v);
1305 static vtkTypeBool JacobiN(double** a, int n, double* w, double** v);
1307
1322 int numberOfSamples, double** xt, int xOrder, double** mt);
1323
1338 static vtkTypeBool SolveLeastSquares(int numberOfSamples, double** xt, int xOrder, double** yt,
1339 int yOrder, double** mt, int checkHomogeneous = 1);
1340
1342
1349 static void RGBToHSV(const float rgb[3], float hsv[3])
1350 {
1351 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1352 }
1353 static void RGBToHSV(float r, float g, float b, float* h, float* s, float* v);
1354 static void RGBToHSV(const double rgb[3], double hsv[3])
1355 {
1356 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1357 }
1358 static void RGBToHSV(double r, double g, double b, double* h, double* s, double* v);
1360
1362
1369 static void HSVToRGB(const float hsv[3], float rgb[3])
1370 {
1371 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1372 }
1373 static void HSVToRGB(float h, float s, float v, float* r, float* g, float* b);
1374 static void HSVToRGB(const double hsv[3], double rgb[3])
1375 {
1376 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1377 }
1378 static void HSVToRGB(double h, double s, double v, double* r, double* g, double* b);
1380
1382
1385 static void LabToXYZ(const double lab[3], double xyz[3])
1386 {
1387 LabToXYZ(lab[0], lab[1], lab[2], xyz + 0, xyz + 1, xyz + 2);
1388 }
1389 static void LabToXYZ(double L, double a, double b, double* x, double* y, double* z);
1391
1393
1396 static void XYZToLab(const double xyz[3], double lab[3])
1397 {
1398 XYZToLab(xyz[0], xyz[1], xyz[2], lab + 0, lab + 1, lab + 2);
1399 }
1400 static void XYZToLab(double x, double y, double z, double* L, double* a, double* b);
1402
1404
1407 static void XYZToRGB(const double xyz[3], double rgb[3])
1408 {
1409 XYZToRGB(xyz[0], xyz[1], xyz[2], rgb + 0, rgb + 1, rgb + 2);
1410 }
1411 static void XYZToRGB(double x, double y, double z, double* r, double* g, double* b);
1413
1415
1418 static void RGBToXYZ(const double rgb[3], double xyz[3])
1419 {
1420 RGBToXYZ(rgb[0], rgb[1], rgb[2], xyz + 0, xyz + 1, xyz + 2);
1421 }
1422 static void RGBToXYZ(double r, double g, double b, double* x, double* y, double* z);
1424
1426
1432 static void RGBToLab(const double rgb[3], double lab[3])
1433 {
1434 RGBToLab(rgb[0], rgb[1], rgb[2], lab + 0, lab + 1, lab + 2);
1435 }
1436 static void RGBToLab(double red, double green, double blue, double* L, double* a, double* b);
1438
1440
1443 static void LabToRGB(const double lab[3], double rgb[3])
1444 {
1445 LabToRGB(lab[0], lab[1], lab[2], rgb + 0, rgb + 1, rgb + 2);
1446 }
1447 static void LabToRGB(double L, double a, double b, double* red, double* green, double* blue);
1449
1451
1454 static void UninitializeBounds(double bounds[6])
1455 {
1456 bounds[0] = 1.0;
1457 bounds[1] = -1.0;
1458 bounds[2] = 1.0;
1459 bounds[3] = -1.0;
1460 bounds[4] = 1.0;
1461 bounds[5] = -1.0;
1462 }
1463
1464
1466
1469 static vtkTypeBool AreBoundsInitialized(const double bounds[6])
1470 {
1471 if (bounds[1] - bounds[0] < 0.0)
1472 {
1473 return 0;
1474 }
1475 return 1;
1476 }
1477
1478
1483 template <class T>
1484 static T ClampValue(const T& value, const T& min, const T& max);
1485
1487
1491 static void ClampValue(double* value, const double range[2]);
1492 static void ClampValue(double value, const double range[2], double* clamped_value);
1493 static void ClampValues(double* values, int nb_values, const double range[2]);
1494 static void ClampValues(
1495 const double* values, int nb_values, const double range[2], double* clamped_values);
1497
1504 static double ClampAndNormalizeValue(double value, const double range[2]);
1505
1510 template <class T1, class T2>
1511 static void TensorFromSymmetricTensor(const T1 symmTensor[6], T2 tensor[9]);
1512
1518 template <class T>
1519 static void TensorFromSymmetricTensor(T tensor[9]);
1520
1530 double range_min, double range_max, double scale = 1.0, double shift = 0.0);
1531
1540 static vtkTypeBool GetAdjustedScalarRange(vtkDataArray* array, int comp, double range[2]);
1541
1546 static vtkTypeBool ExtentIsWithinOtherExtent(const int extent1[6], const int extent2[6]);
1547
1554 const double bounds1[6], const double bounds2[6], const double delta[3]);
1555
1562 const double point[3], const double bounds[6], const double delta[3]);
1563
1574 const double bounds[6], const double normal[3], const double point[3]);
1575
1585 static double Solve3PointCircle(
1586 const double p1[3], const double p2[3], const double p3[3], double center[3]);
1587
1591 static double Inf();
1592
1596 static double NegInf();
1597
1601 static double Nan();
1602
1606 static vtkTypeBool IsInf(double x);
1607
1611 static vtkTypeBool IsNan(double x);
1612
1617 static bool IsFinite(double x);
1618
1623 static int QuadraticRoot(double a, double b, double c, double min, double max, double* u);
1624
1630 static vtkIdType ComputeGCD(vtkIdType m, vtkIdType n) { return (n ? ComputeGCD(n, m % n) : m); }
1631
1636 {
1640 };
1641
1664 template <class Iter1, class Iter2, class Iter3>
1665 static void Convolve1D(Iter1 beginSample, Iter1 endSample, Iter2 beginKernel, Iter2 endKernel,
1666 Iter3 beginOut, Iter3 endOut, ConvolutionMode mode = ConvolutionMode::FULL)
1667 {
1668 int sampleSize = std::distance(beginSample, endSample);
1669 int kernelSize = std::distance(beginKernel, endKernel);
1670 int outSize = std::distance(beginOut, endOut);
1671
1672 if (sampleSize <= 0 || kernelSize <= 0 || outSize <= 0)
1673 {
1674 return;
1675 }
1676
1677 int begin = 0;
1678 int end = outSize;
1679
1680 switch (mode)
1681 {
1683 begin = static_cast<int>(std::ceil((std::min)(sampleSize, kernelSize) / 2.0)) - 1;
1684 end = begin + (std::max)(sampleSize, kernelSize);
1685 break;
1687 begin = (std::min)(sampleSize, kernelSize) - 1;
1688 end = begin + std::abs(sampleSize - kernelSize) + 1;
1689 break;
1691 default:
1692 break;
1693 }
1694
1695 for (int i = begin; i < end; i++)
1696 {
1697 Iter3 out = beginOut + i - begin;
1698 *out = 0;
1699 for (int j = (std::max)(i - sampleSize + 1, 0); j <= (std::min)(i, kernelSize - 1); j++)
1700 {
1701 *out += *(beginSample + (i - j)) * *(beginKernel + j);
1702 }
1703 }
1704 }
1705
1710 static void GetPointAlongLine(double result[3], double p1[3], double p2[3], const double offset)
1711 {
1712 double directionVector[3] = { p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2] };
1713 vtkMath::Normalize(directionVector);
1714 result[0] = p2[0] + (offset * directionVector[0]);
1715 result[1] = p2[1] + (offset * directionVector[1]);
1716 result[2] = p2[2] + (offset * directionVector[2]);
1717 }
1718
1719protected:
1720 vtkMath() = default;
1721 ~vtkMath() override = default;
1722
1724
1725private:
1726 vtkMath(const vtkMath&) = delete;
1727 void operator=(const vtkMath&) = delete;
1728};
1729
1730//----------------------------------------------------------------------------
1731inline float vtkMath::RadiansFromDegrees(float x)
1732{
1733 return x * 0.017453292f;
1734}
1735
1736//----------------------------------------------------------------------------
1737inline double vtkMath::RadiansFromDegrees(double x)
1738{
1739 return x * 0.017453292519943295;
1740}
1741
1742//----------------------------------------------------------------------------
1743inline float vtkMath::DegreesFromRadians(float x)
1744{
1745 return x * 57.2957795131f;
1746}
1747
1748//----------------------------------------------------------------------------
1749inline double vtkMath::DegreesFromRadians(double x)
1750{
1751 return x * 57.29577951308232;
1752}
1753
1754//----------------------------------------------------------------------------
1755inline bool vtkMath::IsPowerOfTwo(vtkTypeUInt64 x)
1756{
1757 return ((x != 0) & ((x & (x - 1)) == 0));
1758}
1759
1760//----------------------------------------------------------------------------
1761// Credit goes to Peter Hart and William Lewis on comp.lang.python 1997
1763{
1764 unsigned int z = static_cast<unsigned int>(((x > 0) ? x - 1 : 0));
1765 z |= z >> 1;
1766 z |= z >> 2;
1767 z |= z >> 4;
1768 z |= z >> 8;
1769 z |= z >> 16;
1770 return static_cast<int>(z + 1);
1771}
1772
1773//----------------------------------------------------------------------------
1774// Modify the trunc() operation provided by static_cast<int>() to get floor(),
1775// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1776inline int vtkMath::Floor(double x)
1777{
1778 int i = static_cast<int>(x);
1779 return i - (i > x);
1780}
1781
1782//----------------------------------------------------------------------------
1783// Modify the trunc() operation provided by static_cast<int>() to get ceil(),
1784// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1785inline int vtkMath::Ceil(double x)
1786{
1787 int i = static_cast<int>(x);
1788 return i + (i < x);
1789}
1790
1791//----------------------------------------------------------------------------
1792template <class T>
1793inline T vtkMath::Min(const T& a, const T& b)
1794{
1795 return (b <= a ? b : a);
1796}
1797
1798//----------------------------------------------------------------------------
1799template <class T>
1800inline T vtkMath::Max(const T& a, const T& b)
1801{
1802 return (b > a ? b : a);
1803}
1804
1805//----------------------------------------------------------------------------
1806inline float vtkMath::Normalize(float v[3])
1807{
1808 float den = vtkMath::Norm(v);
1809 if (den != 0.0)
1810 {
1811 for (int i = 0; i < 3; ++i)
1812 {
1813 v[i] /= den;
1814 }
1815 }
1816 return den;
1817}
1818
1819//----------------------------------------------------------------------------
1820inline double vtkMath::Normalize(double v[3])
1821{
1822 double den = vtkMath::Norm(v);
1823 if (den != 0.0)
1824 {
1825 for (int i = 0; i < 3; ++i)
1826 {
1827 v[i] /= den;
1828 }
1829 }
1830 return den;
1831}
1832
1833//----------------------------------------------------------------------------
1834inline float vtkMath::Normalize2D(float v[2])
1835{
1836 float den = vtkMath::Norm2D(v);
1837 if (den != 0.0)
1838 {
1839 for (int i = 0; i < 2; ++i)
1840 {
1841 v[i] /= den;
1842 }
1843 }
1844 return den;
1845}
1846
1847//----------------------------------------------------------------------------
1848inline double vtkMath::Normalize2D(double v[2])
1849{
1850 double den = vtkMath::Norm2D(v);
1851 if (den != 0.0)
1852 {
1853 for (int i = 0; i < 2; ++i)
1854 {
1855 v[i] /= den;
1856 }
1857 }
1858 return den;
1859}
1860
1861//----------------------------------------------------------------------------
1862inline float vtkMath::Determinant3x3(const float c1[3], const float c2[3], const float c3[3])
1863{
1864 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
1865 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
1866}
1867
1868//----------------------------------------------------------------------------
1869inline double vtkMath::Determinant3x3(const double c1[3], const double c2[3], const double c3[3])
1870{
1871 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
1872 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
1873}
1874
1875//----------------------------------------------------------------------------
1877 double a1, double a2, double a3, double b1, double b2, double b3, double c1, double c2, double c3)
1878{
1879 return (a1 * vtkMath::Determinant2x2(b2, b3, c2, c3) -
1880 b1 * vtkMath::Determinant2x2(a2, a3, c2, c3) + c1 * vtkMath::Determinant2x2(a2, a3, b2, b3));
1881}
1882
1883//----------------------------------------------------------------------------
1884inline float vtkMath::Distance2BetweenPoints(const float p1[3], const float p2[3])
1885{
1886 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
1887 (p1[2] - p2[2]) * (p1[2] - p2[2]));
1888}
1889
1890//----------------------------------------------------------------------------
1891inline double vtkMath::Distance2BetweenPoints(const double p1[3], const double p2[3])
1892{
1893 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
1894 (p1[2] - p2[2]) * (p1[2] - p2[2]));
1895}
1896
1897//----------------------------------------------------------------------------
1898template <typename ReturnTypeT, typename TupleRangeT1, typename TupleRangeT2, typename EnableT>
1899inline ReturnTypeT vtkMath::Distance2BetweenPoints(const TupleRangeT1& p1, const TupleRangeT2& p2)
1900{
1901 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
1902 (p1[2] - p2[2]) * (p1[2] - p2[2]));
1903}
1904
1905//----------------------------------------------------------------------------
1906template <class VectorT1, class VectorT2, class VectorT3>
1907void vtkMath::Cross(VectorT1&& a, VectorT2&& b, VectorT3& c)
1908{
1910 ValueType Cx = a[1] * b[2] - a[2] * b[1];
1911 ValueType Cy = a[2] * b[0] - a[0] * b[2];
1912 ValueType Cz = a[0] * b[1] - a[1] * b[0];
1913 c[0] = Cx;
1914 c[1] = Cy;
1915 c[2] = Cz;
1916}
1917
1918//----------------------------------------------------------------------------
1919// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
1920inline void vtkMath::Cross(const float a[3], const float b[3], float c[3])
1921{
1922 float Cx = a[1] * b[2] - a[2] * b[1];
1923 float Cy = a[2] * b[0] - a[0] * b[2];
1924 float Cz = a[0] * b[1] - a[1] * b[0];
1925 c[0] = Cx;
1926 c[1] = Cy;
1927 c[2] = Cz;
1928}
1929
1930//----------------------------------------------------------------------------
1931// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
1932inline void vtkMath::Cross(const double a[3], const double b[3], double c[3])
1933{
1934 double Cx = a[1] * b[2] - a[2] * b[1];
1935 double Cy = a[2] * b[0] - a[0] * b[2];
1936 double Cz = a[0] * b[1] - a[1] * b[0];
1937 c[0] = Cx;
1938 c[1] = Cy;
1939 c[2] = Cz;
1940}
1941
1942//----------------------------------------------------------------------------
1943template <class T>
1944inline double vtkDeterminant3x3(const T A[3][3])
1945{
1946 return A[0][0] * A[1][1] * A[2][2] + A[1][0] * A[2][1] * A[0][2] + A[2][0] * A[0][1] * A[1][2] -
1947 A[0][0] * A[2][1] * A[1][2] - A[1][0] * A[0][1] * A[2][2] - A[2][0] * A[1][1] * A[0][2];
1948}
1949
1950//----------------------------------------------------------------------------
1951inline double vtkMath::Determinant3x3(const float A[3][3])
1952{
1953 return vtkDeterminant3x3(A);
1954}
1955
1956//----------------------------------------------------------------------------
1957inline double vtkMath::Determinant3x3(const double A[3][3])
1958{
1959 return vtkDeterminant3x3(A);
1960}
1961
1962//----------------------------------------------------------------------------
1963template <class T>
1964inline T vtkMath::ClampValue(const T& value, const T& min, const T& max)
1965{
1966 assert("pre: valid_range" && min <= max);
1967
1968#if __cplusplus >= 201703L
1969 return std::clamp(value, min, max);
1970#else
1971 // compilers are good at optimizing the ternary operator,
1972 // use '<' since it is preferred by STL for custom types
1973 T v = (min < value ? value : min);
1974 return (v < max ? v : max);
1975#endif
1976}
1977
1978//----------------------------------------------------------------------------
1979inline void vtkMath::ClampValue(double* value, const double range[2])
1980{
1981 if (value && range)
1982 {
1983 assert("pre: valid_range" && range[0] <= range[1]);
1984
1985 *value = vtkMath::ClampValue(*value, range[0], range[1]);
1986 }
1987}
1988
1989//----------------------------------------------------------------------------
1990inline void vtkMath::ClampValue(double value, const double range[2], double* clamped_value)
1991{
1992 if (range && clamped_value)
1993 {
1994 assert("pre: valid_range" && range[0] <= range[1]);
1995
1996 *clamped_value = vtkMath::ClampValue(value, range[0], range[1]);
1997 }
1998}
1999
2000// ---------------------------------------------------------------------------
2001inline double vtkMath::ClampAndNormalizeValue(double value, const double range[2])
2002{
2003 assert("pre: valid_range" && range[0] <= range[1]);
2004
2005 double result;
2006 if (range[0] == range[1])
2007 {
2008 result = 0.0;
2009 }
2010 else
2011 {
2012 // clamp
2013 result = vtkMath::ClampValue(value, range[0], range[1]);
2014
2015 // normalize
2016 result = (result - range[0]) / (range[1] - range[0]);
2017 }
2018
2019 assert("post: valid_result" && result >= 0.0 && result <= 1.0);
2020
2021 return result;
2022}
2023
2024//-----------------------------------------------------------------------------
2025template <class T1, class T2>
2026inline void vtkMath::TensorFromSymmetricTensor(const T1 symmTensor[9], T2 tensor[9])
2027{
2028 for (int i = 0; i < 3; ++i)
2029 {
2030 tensor[4 * i] = symmTensor[i];
2031 }
2032 tensor[1] = tensor[3] = symmTensor[3];
2033 tensor[2] = tensor[6] = symmTensor[5];
2034 tensor[5] = tensor[7] = symmTensor[4];
2035}
2036
2037//-----------------------------------------------------------------------------
2038template <class T>
2040{
2041 tensor[6] = tensor[5]; // XZ
2042 tensor[7] = tensor[4]; // YZ
2043 tensor[8] = tensor[2]; // ZZ
2044 tensor[4] = tensor[1]; // YY
2045 tensor[5] = tensor[7]; // YZ
2046 tensor[2] = tensor[6]; // XZ
2047 tensor[1] = tensor[3]; // XY
2048}
2049VTK_ABI_NAMESPACE_END
2050
2051namespace
2052{
2053template <class QuaternionT, class MatrixT>
2054inline void vtkQuaternionToMatrix3x3(QuaternionT&& quat, MatrixT&& A)
2055{
2057
2058 Scalar ww = quat[0] * quat[0];
2059 Scalar wx = quat[0] * quat[1];
2060 Scalar wy = quat[0] * quat[2];
2061 Scalar wz = quat[0] * quat[3];
2062
2063 Scalar xx = quat[1] * quat[1];
2064 Scalar yy = quat[2] * quat[2];
2065 Scalar zz = quat[3] * quat[3];
2066
2067 Scalar xy = quat[1] * quat[2];
2068 Scalar xz = quat[1] * quat[3];
2069 Scalar yz = quat[2] * quat[3];
2070
2071 Scalar rr = xx + yy + zz;
2072 // normalization factor, just in case quaternion was not normalized
2073 Scalar f = 1 / (ww + rr);
2074 Scalar s = (ww - rr) * f;
2075 f *= 2;
2076
2078
2079 Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) = xx * f + s;
2080 Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) = (xy + wz) * f;
2081 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A)) = (xz - wy) * f;
2082
2083 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A)) = (xy - wz) * f;
2084 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) = yy * f + s;
2085 Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) = (yz + wx) * f;
2086
2087 Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) = (xz + wy) * f;
2088 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A)) = (yz - wx) * f;
2089 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A)) = zz * f + s;
2090}
2091} // anonymous namespace
2092
2093VTK_ABI_NAMESPACE_BEGIN
2094//------------------------------------------------------------------------------
2095inline void vtkMath::QuaternionToMatrix3x3(const float quat[4], float A[3][3])
2096{
2097 vtkQuaternionToMatrix3x3(quat, A);
2098}
2099
2100//------------------------------------------------------------------------------
2101inline void vtkMath::QuaternionToMatrix3x3(const double quat[4], double A[3][3])
2102{
2103 vtkQuaternionToMatrix3x3(quat, A);
2104}
2105
2106//-----------------------------------------------------------------------------
2107template <class QuaternionT, class MatrixT, class EnableT>
2108inline void vtkMath::QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A)
2109{
2110 vtkQuaternionToMatrix3x3(std::forward<QuaternionT>(q), std::forward<MatrixT>(A));
2111}
2112VTK_ABI_NAMESPACE_END
2113
2114namespace
2115{
2116//------------------------------------------------------------------------------
2117// The solution is based on
2118// Berthold K. P. Horn (1987),
2119// "Closed-form solution of absolute orientation using unit quaternions,"
2120// Journal of the Optical Society of America A, 4:629-642
2121template <class MatrixT, class QuaternionT>
2122inline void vtkMatrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& quat)
2123{
2125
2126 Scalar N[4][4];
2127
2129
2130 // on-diagonal elements
2131 N[0][0] = Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) +
2132 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) +
2133 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2134 N[1][1] = Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) -
2135 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) -
2136 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2137 N[2][2] = -Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) +
2138 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) -
2139 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2140 N[3][3] = -Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) -
2141 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) +
2142 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2143
2144 // off-diagonal elements
2145 N[0][1] = N[1][0] = Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) -
2146 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A));
2147 N[0][2] = N[2][0] = Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) -
2148 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A));
2149 N[0][3] = N[3][0] = Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) -
2150 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A));
2151
2152 N[1][2] = N[2][1] = Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) +
2153 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A));
2154 N[1][3] = N[3][1] = Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) +
2155 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A));
2156 N[2][3] = N[3][2] = Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) +
2157 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A));
2158
2159 Scalar eigenvectors[4][4], eigenvalues[4];
2160
2161 // convert into format that JacobiN can use,
2162 // then use Jacobi to find eigenvalues and eigenvectors
2163 Scalar *NTemp[4], *eigenvectorsTemp[4];
2164 for (int i = 0; i < 4; ++i)
2165 {
2166 NTemp[i] = N[i];
2167 eigenvectorsTemp[i] = eigenvectors[i];
2168 }
2169 vtkMath::JacobiN(NTemp, 4, eigenvalues, eigenvectorsTemp);
2170
2171 // the first eigenvector is the one we want
2172 quat[0] = eigenvectors[0][0];
2173 quat[1] = eigenvectors[1][0];
2174 quat[2] = eigenvectors[2][0];
2175 quat[3] = eigenvectors[3][0];
2176}
2177} // anonymous namespace
2178
2179VTK_ABI_NAMESPACE_BEGIN
2180//------------------------------------------------------------------------------
2181inline void vtkMath::Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
2182{
2183 vtkMatrix3x3ToQuaternion(A, quat);
2184}
2185
2186//------------------------------------------------------------------------------
2187inline void vtkMath::Matrix3x3ToQuaternion(const double A[3][3], double quat[4])
2188{
2189 vtkMatrix3x3ToQuaternion(A, quat);
2190}
2191
2192//-----------------------------------------------------------------------------
2193template <class MatrixT, class QuaternionT, class EnableT>
2194inline void vtkMath::Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q)
2195{
2196 vtkMatrix3x3ToQuaternion(std::forward<MatrixT>(A), std::forward<QuaternionT>(q));
2197}
2198VTK_ABI_NAMESPACE_END
2199
2200namespace vtk_detail
2201{
2202VTK_ABI_NAMESPACE_BEGIN
2203// Can't specialize templates inside a template class, so we move the impl here.
2204template <typename OutT>
2205void RoundDoubleToIntegralIfNecessary(double val, OutT* ret)
2206{ // OutT is integral -- clamp and round
2207 if (!vtkMath::IsNan(val))
2208 {
2209 double min = static_cast<double>(vtkTypeTraits<OutT>::Min());
2210 double max = static_cast<double>(vtkTypeTraits<OutT>::Max());
2211 val = vtkMath::ClampValue(val, min, max);
2212 *ret = static_cast<OutT>((val >= 0.0) ? (val + 0.5) : (val - 0.5));
2213 }
2214 else
2215 *ret = 0;
2216}
2217template <>
2218inline void RoundDoubleToIntegralIfNecessary(double val, double* retVal)
2219{ // OutT is double: passthrough
2220 *retVal = val;
2221}
2222template <>
2223inline void RoundDoubleToIntegralIfNecessary(double val, float* retVal)
2224{ // OutT is float -- just clamp (as doubles, then the cast to float is well-defined.)
2225 if (!vtkMath::IsNan(val))
2226 {
2227 double min = static_cast<double>(vtkTypeTraits<float>::Min());
2228 double max = static_cast<double>(vtkTypeTraits<float>::Max());
2229 val = vtkMath::ClampValue(val, min, max);
2230 }
2231
2232 *retVal = static_cast<float>(val);
2233}
2234VTK_ABI_NAMESPACE_END
2235} // end namespace vtk_detail
2236
2237VTK_ABI_NAMESPACE_BEGIN
2238//-----------------------------------------------------------------------------
2239#if defined(VTK_HAS_ISINF) || defined(VTK_HAS_STD_ISINF)
2240#define VTK_MATH_ISINF_IS_INLINE
2241inline vtkTypeBool vtkMath::IsInf(double x)
2242{
2243#if defined(VTK_HAS_STD_ISINF)
2244 return std::isinf(x);
2245#else
2246 return (isinf(x) != 0); // Force conversion to bool
2247#endif
2248}
2249#endif
2250
2251//-----------------------------------------------------------------------------
2252#if defined(VTK_HAS_ISNAN) || defined(VTK_HAS_STD_ISNAN)
2253#define VTK_MATH_ISNAN_IS_INLINE
2254inline vtkTypeBool vtkMath::IsNan(double x)
2255{
2256#if defined(VTK_HAS_STD_ISNAN)
2257 return std::isnan(x);
2258#else
2259 return (isnan(x) != 0); // Force conversion to bool
2260#endif
2261}
2262#endif
2263
2264//-----------------------------------------------------------------------------
2265#if defined(VTK_HAS_ISFINITE) || defined(VTK_HAS_STD_ISFINITE) || defined(VTK_HAS_FINITE)
2266#define VTK_MATH_ISFINITE_IS_INLINE
2267inline bool vtkMath::IsFinite(double x)
2268{
2269#if defined(VTK_HAS_STD_ISFINITE)
2270 return std::isfinite(x);
2271#elif defined(VTK_HAS_ISFINITE)
2272 return (isfinite(x) != 0); // Force conversion to bool
2273#else
2274 return (finite(x) != 0); // Force conversion to bool
2275#endif
2276}
2277#endif
2278
2279VTK_ABI_NAMESPACE_END
2280#endif
RealT mt
Definition PyrC2Basis.h:39
RealT ww
Definition PyrI2Basis.h:13
Gaussian sequence of pseudo random numbers implemented with the Box-Mueller transform.
a simple class to control print indentation
Definition vtkIndent.h:29
static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1 &p1, const TupleRangeT2 &p2)
Compute distance squared between two points p1 and p2.
Definition vtkMath.h:1899
static void Multiply3x3(const float A[3][3], const float B[3][3], float C[3][3])
Multiply one 3x3 matrix by another according to C = AB.
static double Dot(const double a[3], const double b[3])
Dot product of two 3-vectors (double version).
Definition vtkMath.h:480
static int GetScalarTypeFittingRange(double range_min, double range_max, double scale=1.0, double shift=0.0)
Return the scalar type that is most likely to have enough precision to store a given range of data on...
static void RGBToXYZ(double r, double g, double b, double *x, double *y, double *z)
Convert color from the RGB system to CIE XYZ.
static void Multiply3x3(const double A[3][3], const double B[3][3], double C[3][3])
Multiply one 3x3 matrix by another according to C = AB.
static double Norm(const double *x, int n)
Compute the norm of n-vector.
static int Round(float f)
Rounds a float to the nearest integer.
Definition vtkMath.h:140
static vtkIdType ComputeGCD(vtkIdType m, vtkIdType n)
Compute the greatest common divisor (GCD) of two positive integers m and n.
Definition vtkMath.h:1630
static double Norm2D(const double x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:756
static double GaussianAmplitude(double variance, double distanceFromMean)
Compute the amplitude of a Gaussian function with mean=0 and specified variance.
static void XYZToRGB(double x, double y, double z, double *r, double *g, double *b)
Convert color from the CIE XYZ system to RGB.
static void GetPointAlongLine(double result[3], double p1[3], double p2[3], const double offset)
Get the coordinates of a point along a line defined by p1 and p2, at a specified offset relative to p...
Definition vtkMath.h:1710
static void Subtract(const float a[3], const float b[3], float c[3])
Subtraction of two 3-vectors (float version).
Definition vtkMath.h:389
static void LUSolve3x3(const double A[3][3], const int index[3], double x[3])
LU back substitution for a 3x3 matrix.
static vtkTypeBool SolveHomogeneousLeastSquares(int numberOfSamples, double **xt, int xOrder, double **mt)
Solves for the least squares best fit matrix for the homogeneous equation X'M' = 0'.
static void Outer2D(const float x[2], const float y[2], float A[2][2])
Outer product of two 2-vectors (float version).
Definition vtkMath.h:721
static bool ProjectVector(const double a[3], const double b[3], double projection[3])
Compute the projection of vector a on vector b and return it in projection[3].
static vtkSmartPointer< vtkMathInternal > Internal
Definition vtkMath.h:1723
static float Norm(const float *x, int n)
Compute the norm of n-vector.
static vtkTypeBool ExtentIsWithinOtherExtent(const int extent1[6], const int extent2[6])
Return true if first 3D extent is within second 3D extent Extent is x-min, x-max, y-min,...
static double GaussianAmplitude(double mean, double variance, double position)
Compute the amplitude of a Gaussian function with specified mean and variance.
static void Add(const double a[3], const double b[3], double c[3])
Addition of two 3-vectors (double version).
Definition vtkMath.h:364
static void RGBToHSV(float r, float g, float b, float *h, float *s, float *v)
Convert color in RGB format (Red, Green, Blue) to HSV format (Hue, Saturation, Value).
static float Norm(const float v[3])
Compute the norm of 3-vector (float version).
Definition vtkMath.h:568
static ReturnTypeT Dot(const TupleRangeT1 &a, const TupleRangeT2 &b)
Compute dot product between two points p1 and p2.
Definition vtkMath.h:504
static vtkTypeBool Jacobi(double **a, double *w, double **v)
Jacobi iteration for the solution of eigenvectors/eigenvalues of a 3x3 real symmetric matrix.
static ScalarT Dot(VectorT1 &&x, VectorT2 &&y)
Computes the dot product between 2 vectors x and y.
Definition vtkMath.h:915
static void XYZToLab(const double xyz[3], double lab[3])
Convert Color from the CIE XYZ system to CIE-L*ab.
Definition vtkMath.h:1396
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static vtkTypeInt64 Factorial(int N)
Compute N factorial, N!
static vtkTypeInt64 Binomial(int m, int n)
The number of combinations of n objects from a pool of m objects (m>n).
static double Random()
Generate pseudo-random numbers distributed according to the uniform distribution between 0....
static void Identity3x3(float A[3][3])
Set A to the identity matrix.
static void SingularValueDecomposition3x3(const float A[3][3], float U[3][3], float w[3], float VT[3][3])
Perform singular value decomposition on a 3x3 matrix.
static double Nan()
Special IEEE-754 number used to represent Not-A-Number (Nan).
static void Perpendiculars(const float v1[3], float v2[3], float v3[3], double theta)
Given a unit vector v1, find two unit vectors v2 and v3 such that v1 cross v2 = v3 (i....
static double Gaussian(double mean, double std)
Generate pseudo-random numbers distributed according to the Gaussian distribution with mean mean and ...
static bool IsFinite(double x)
Test if a number has finite value i.e.
static void LUSolveLinearSystem(double **A, int *index, double *x, int size)
Solve linear equations Ax = b using LU decomposition A = LU where L is lower triangular matrix and U ...
static double EstimateMatrixCondition(const double *const *A, int size)
Estimate the condition number of a LU factored matrix.
static void LUFactor3x3(float A[3][3], int index[3])
LU Factorization of a 3x3 matrix.
static void LinearSolve(MatrixT &&M, VectorT1 &&x, VectorT2 &&y)
This method solves linear systems M * x = y.
Definition vtkMath.h:1002
static void FreeCombination(int *combination)
Free the "iterator" array created by vtkMath::BeginCombination.
static double Random(double min, double max)
Generate pseudo-random numbers distributed according to the uniform distribution between min and max.
static void TensorFromSymmetricTensor(const T1 symmTensor[6], T2 tensor[9])
Convert a 6-Component symmetric tensor into a 9-Component tensor, no allocation performed.
static void LabToXYZ(const double lab[3], double xyz[3])
Convert color from the CIE-L*ab system to CIE XYZ.
Definition vtkMath.h:1385
static double Solve3PointCircle(const double p1[3], const double p2[3], const double p3[3], double center[3])
In Euclidean space, there is a unique circle passing through any given three non-collinear points P1,...
static vtkTypeBool PointIsWithinBounds(const double point[3], const double bounds[6], const double delta[3])
Return true if point is within the given 3D bounds Bounds is x-min, x-max, y-min, y-max,...
static float Dot(const float a[3], const float b[3])
Dot product of two 3-vectors (float version).
Definition vtkMath.h:472
static void Diagonalize3x3(const float A[3][3], float w[3], float V[3][3])
Diagonalize a symmetric 3x3 matrix and return the eigenvalues in w and the eigenvectors in the column...
static void LabToXYZ(double L, double a, double b, double *x, double *y, double *z)
Convert color from the CIE-L*ab system to CIE XYZ.
static vtkTypeBool GetAdjustedScalarRange(vtkDataArray *array, int comp, double range[2])
Get a vtkDataArray's scalar range for a given component.
static bool ProjectVector(const float a[3], const float b[3], float projection[3])
Compute the projection of vector a on vector b and return it in projection[3].
static void MultiplyScalar2D(float a[2], float s)
Multiplies a 2-vector by a scalar (float version).
Definition vtkMath.h:437
static void HSVToRGB(const float hsv[3], float rgb[3])
Convert color in HSV format (Hue, Saturation, Value) to RGB format (Red, Green, Blue).
Definition vtkMath.h:1369
static void Assign(const double a[3], double b[3])
Assign values to a 3-vector (double version).
Definition vtkMath.h:348
static double Determinant2x2(const double c1[2], const double c2[2])
Calculate the determinant of a 2x2 matrix: | a b | | c d |.
Definition vtkMath.h:783
static T Max(const T &a, const T &b)
Returns the maximum of the two arguments provided.
Definition vtkMath.h:1800
static void Outer2D(const double x[2], const double y[2], double A[2][2])
Outer product of two 2-vectors (double version).
Definition vtkMath.h:735
static void RandomSeed(int s)
Initialize seed value.
static double NegInf()
Special IEEE-754 number used to represent negative infinity.
static void MultiplyScalar2D(double a[2], double s)
Multiplies a 2-vector by a scalar (double version).
Definition vtkMath.h:461
static void LabToRGB(double L, double a, double b, double *red, double *green, double *blue)
Convert color from the CIE-L*ab system to RGB.
static double Gaussian()
Generate pseudo-random numbers distributed according to the standard normal distribution.
static int Ceil(double x)
Rounds a double to the nearest integer not less than itself.
Definition vtkMath.h:1785
static void HSVToRGB(const double hsv[3], double rgb[3])
Convert color in HSV format (Hue, Saturation, Value) to RGB format (Red, Green, Blue).
Definition vtkMath.h:1374
~vtkMath() override=default
static double Inf()
Special IEEE-754 number used to represent positive infinity.
static vtkTypeBool Jacobi(float **a, float *w, float **v)
Jacobi iteration for the solution of eigenvectors/eigenvalues of a 3x3 real symmetric matrix.
static int PlaneIntersectsAABB(const double bounds[6], const double normal[3], const double point[3])
Implements Plane / Axis-Aligned Bounding-Box intersection as described in Graphics Gems IV,...
static ScalarT Dot(VectorT1 &&x, VectorT2 &&y)
Computes the dot product between 2 vectors x and y.
Definition vtkMath.h:898
static void RGBToXYZ(const double rgb[3], double xyz[3])
Convert color from the RGB system to CIE XYZ.
Definition vtkMath.h:1418
static void QuaternionToMatrix3x3(const float quat[4], float A[3][3])
Convert a quaternion to a 3x3 rotation matrix.
Definition vtkMath.h:2095
static int NearestPowerOfTwo(int x)
Compute the nearest power of two that is not less than x.
Definition vtkMath.h:1762
static void HSVToRGB(double h, double s, double v, double *r, double *g, double *b)
Convert color in HSV format (Hue, Saturation, Value) to RGB format (Red, Green, Blue).
static void SingularValueDecomposition3x3(const double A[3][3], double U[3][3], double w[3], double VT[3][3])
Perform singular value decomposition on a 3x3 matrix.
static double SignedAngleBetweenVectors(const double v1[3], const double v2[3], const double vn[3])
Compute signed angle in radians between two vectors with regard to a third orthogonal vector.
static ScalarT Dot(VectorT1 &&x, MatrixT &&M, VectorT2 &&y)
Computes the dot product x^T M y, where x and y are vectors and M is a metric matrix.
Definition vtkMath.h:1025
static float Normalize2D(float v[2])
Normalize (in place) a 2-vector.
Definition vtkMath.h:1834
static void Invert3x3(const double A[3][3], double AI[3][3])
Invert a 3x3 matrix.
static void HSVToRGB(float h, float s, float v, float *r, float *g, float *b)
Convert color in HSV format (Hue, Saturation, Value) to RGB format (Red, Green, Blue).
static constexpr int DYNAMIC_VECTOR_SIZE()
When this value is passed to a select templated functions in vtkMath, the computation can be performe...
Definition vtkMath.h:113
static void MultiplyQuaternion(const double q1[4], const double q2[4], double q[4])
Multiply two quaternions.
static void Multiply3x3(const double A[3][3], const double v[3], double u[3])
Multiply a vector by a 3x3 matrix.
static void Outer(const double a[3], const double b[3], double c[3][3])
Outer product of two 3-vectors (double version).
Definition vtkMath.h:526
static vtkTypeBool InvertMatrix(double **A, double **AI, int size, int *tmp1Size, double *tmp2Size)
Thread safe version of InvertMatrix method.
static vtkTypeBool InvertMatrix(double **A, double **AI, int size)
Invert input square matrix A into matrix AI.
static void LUSolve3x3(const float A[3][3], const int index[3], float x[3])
LU back substitution for a 3x3 matrix.
static int GetSeed()
Return the current seed used by the random number generator.
static void Assign(const VectorT1 &a, VectorT2 &&b)
Assign values to a 3-vector (templated version).
Definition vtkMath.h:338
static float RadiansFromDegrees(float degrees)
Convert degrees into radians.
Definition vtkMath.h:1731
static void Convolve1D(Iter1 beginSample, Iter1 endSample, Iter2 beginKernel, Iter2 endKernel, Iter3 beginOut, Iter3 endOut, ConvolutionMode mode=ConvolutionMode::FULL)
Compute the convolution of a sampled 1D signal by a given kernel.
Definition vtkMath.h:1665
static void RotateVectorByWXYZ(const double v[3], const double q[4], double r[3])
rotate a vector by WXYZ using // https://en.wikipedia.org/wiki/Rodrigues%27_rotation_formula
static void Add(const float a[3], const float b[3], float c[3])
Addition of two 3-vectors (float version).
Definition vtkMath.h:353
static int CeilLog2(vtkTypeUInt64 x)
Gives the exponent of the lowest power of two not less than x.
static vtkTypeBool AreBoundsInitialized(const double bounds[6])
Are the bounds initialized?
Definition vtkMath.h:1469
static bool ProjectVector2D(const double a[2], const double b[2], double projection[2])
Compute the projection of 2D vector a on 2D vector b and returns the result in projection[2].
static vtkTypeBool JacobiN(float **a, int n, float *w, float **v)
JacobiN iteration for the solution of eigenvectors/eigenvalues of a nxn real symmetric matrix.
static int NextCombination(int m, int n, int *combination)
Given m, n, and a valid combination of n integers in the range [0,m[, this function alters the intege...
static constexpr double Pi()
A mathematical constant.
Definition vtkMath.h:118
static void Multiply3x3(const float A[3][3], const float v[3], float u[3])
Multiply a vector by a 3x3 matrix.
static void Subtract(const double a[3], const double b[3], double c[3])
Subtraction of two 3-vectors (double version).
Definition vtkMath.h:400
static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
Convert a 3x3 matrix into a quaternion.
Definition vtkMath.h:2181
static vtkMath * New()
static void Orthogonalize3x3(const double A[3][3], double B[3][3])
Orthogonalize a 3x3 matrix and put the result in B.
static void XYZToRGB(const double xyz[3], double rgb[3])
Convert color from the CIE XYZ system to RGB.
Definition vtkMath.h:1407
static double ClampAndNormalizeValue(double value, const double range[2])
Clamp a value against a range and then normalize it between 0 and 1.
Definition vtkMath.h:2001
static void MultiplyScalar(double a[3], double s)
Multiplies a 3-vector by a scalar (double version).
Definition vtkMath.h:449
static double Dot2D(const double x[2], const double y[2])
Dot product of two 2-vectors.
Definition vtkMath.h:716
static void LinearSolve3x3(const float A[3][3], const float x[3], float y[3])
Solve Ay = x for y and place the result in y.
static vtkTypeBool IsNan(double x)
Test if a number is equal to the special floating point value Not-A-Number (Nan).
static void Diagonalize3x3(const double A[3][3], double w[3], double V[3][3])
Diagonalize a symmetric 3x3 matrix and return the eigenvalues in w and the eigenvectors in the column...
static void RGBToLab(const double rgb[3], double lab[3])
Convert color from the RGB system to CIE-L*ab.
Definition vtkMath.h:1432
static int Floor(double x)
Rounds a double to the nearest integer not greater than itself.
Definition vtkMath.h:1776
static void RotateVectorByNormalizedQuaternion(const double v[3], const double q[4], double r[3])
rotate a vector by a normalized quaternion using // https://en.wikipedia.org/wiki/Rodrigues%27_rotati...
static void Subtract(const VectorT1 &a, const VectorT2 &b, VectorT3 &&c)
Subtraction of two 3-vectors (templated version).
Definition vtkMath.h:414
static vtkTypeBool BoundsIsWithinOtherBounds(const double bounds1[6], const double bounds2[6], const double delta[3])
Return true if first 3D bounds is within the second 3D bounds Bounds is x-min, x-max,...
static double Determinant2x2(double a, double b, double c, double d)
Calculate the determinant of a 2x2 matrix: | a b | | c d |.
Definition vtkMath.h:782
static void RGBToHSV(const double rgb[3], double hsv[3])
Convert color in RGB format (Red, Green, Blue) to HSV format (Hue, Saturation, Value).
Definition vtkMath.h:1354
static vtkTypeBool JacobiN(double **a, int n, double *w, double **v)
JacobiN iteration for the solution of eigenvectors/eigenvalues of a nxn real symmetric matrix.
static double AngleBetweenVectors(const double v1[3], const double v2[3])
Compute angle in radians between two vectors.
static void MultiplyMatrix(const double *const *A, const double *const *B, unsigned int rowA, unsigned int colA, unsigned int rowB, unsigned int colB, double **C)
General matrix multiplication.
static float DegreesFromRadians(float radians)
Convert radians into degrees.
Definition vtkMath.h:1743
static float Determinant2x2(const float c1[2], const float c2[2])
Compute determinant of 2x2 matrix.
Definition vtkMath.h:773
static int Round(double f)
Definition vtkMath.h:141
static vtkTypeBool IsInf(double x)
Test if a number is equal to the special floating point value infinity.
static double GaussianWeight(double mean, double variance, double position)
Compute the amplitude of an unnormalized Gaussian function with specified mean and variance.
static void UninitializeBounds(double bounds[6])
Set the bounds to an uninitialized state.
Definition vtkMath.h:1454
vtkMath()=default
static void RGBToHSV(double r, double g, double b, double *h, double *s, double *v)
Convert color in RGB format (Red, Green, Blue) to HSV format (Hue, Saturation, Value).
static void Outer(const float a[3], const float b[3], float c[3][3])
Outer product of two 3-vectors (float version).
Definition vtkMath.h:512
static int * BeginCombination(int m, int n)
Start iterating over "m choose n" objects.
static double Norm(const double v[3])
Compute the norm of 3-vector (double version).
Definition vtkMath.h:573
static void RoundDoubleToIntegralIfNecessary(double val, OutT *ret)
Round a double to type OutT if OutT is integral, otherwise simply clamp the value to the output range...
Definition vtkMath.h:149
static void RotateVectorByWXYZ(const float v[3], const float q[4], float r[3])
rotate a vector by WXYZ using // https://en.wikipedia.org/wiki/Rodrigues%27_rotation_formula
static bool IsPowerOfTwo(vtkTypeUInt64 x)
Returns true if integer is a power of two.
Definition vtkMath.h:1755
static void Invert3x3(const float A[3][3], float AI[3][3])
Invert a 3x3 matrix.
static float Normalize(float v[3])
Normalize (in place) a 3-vector.
Definition vtkMath.h:1806
static void Transpose3x3(const double A[3][3], double AT[3][3])
Transpose a 3x3 matrix.
static ReturnTypeT SquaredNorm(const TupleRangeT &v)
Compute the squared norm of a 3-vector.
Definition vtkMath.h:588
static double Determinant3x3(const float A[3][3])
Return the determinant of a 3x3 matrix.
Definition vtkMath.h:1951
static float Dot2D(const float x[2], const float y[2])
Dot product of two 2-vectors.
Definition vtkMath.h:711
ConvolutionMode
Support the convolution operations.
Definition vtkMath.h:1636
static void RotateVectorByNormalizedQuaternion(const float v[3], const float q[4], float r[3])
rotate a vector by a normalized quaternion using // https://en.wikipedia.org/wiki/Rodrigues%27_rotati...
static void RGBToHSV(const float rgb[3], float hsv[3])
Convert color in RGB format (Red, Green, Blue) to HSV format (Hue, Saturation, Value).
Definition vtkMath.h:1349
static void Add(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Addition of two 3-vectors (double version).
Definition vtkMath.h:378
static void Orthogonalize3x3(const float A[3][3], float B[3][3])
Orthogonalize a 3x3 matrix and put the result in B.
static bool ProjectVector2D(const float a[2], const float b[2], float projection[2])
Compute the projection of 2D vector a on 2D vector b and returns the result in projection[2].
static vtkTypeBool SolveLinearSystemGEPP2x2(double a00, double a01, double a10, double a11, double b0, double b1, double &x0, double &x1)
Solve linear equation Ax = b using Gaussian Elimination with Partial Pivoting for a 2x2 system.
static vtkMatrixUtilities::ScalarTypeExtractor< MatrixT >::value_type Determinant(MatrixT &&M)
Computes the determinant of input square SizeT x SizeT matrix M.
Definition vtkMath.h:958
static vtkTypeBool SolveLinearSystem(double **A, double *x, int size)
Solve linear equations Ax = b using Crout's method.
static void LabToRGB(const double lab[3], double rgb[3])
Convert color from the CIE-L*ab system to RGB.
Definition vtkMath.h:1443
static float Norm2D(const float x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:750
static vtkTypeBool LUFactorLinearSystem(double **A, int *index, int size, double *tmpSize)
Thread safe version of LUFactorLinearSystem method.
static void LinearSolve3x3(const double A[3][3], const double x[3], double y[3])
Solve Ay = x for y and place the result in y.
static void XYZToLab(double x, double y, double z, double *L, double *a, double *b)
Convert Color from the CIE XYZ system to CIE-L*ab.
static void MultiplyScalar(float a[3], float s)
Multiplies a 3-vector by a scalar (float version).
Definition vtkMath.h:425
static T Min(const T &a, const T &b)
Returns the minimum of the two arguments provided.
Definition vtkMath.h:1793
static void InvertMatrix(MatrixT1 &&M1, MatrixT2 &&M2)
Computes the inverse of input matrix M1 into M2.
Definition vtkMath.h:981
static void Cross(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Cross product of two 3-vectors.
Definition vtkMath.h:1907
static void MultiplyMatrix(MatrixT1 &&M1, MatrixT2 &&M2, MatrixT3 &&M3)
Multiply matrices such that M3 = M1 x M2.
Definition vtkMath.h:857
static void Perpendiculars(const double v1[3], double v2[3], double v3[3], double theta)
Given a unit vector v1, find two unit vectors v2 and v3 such that v1 cross v2 = v3 (i....
static T ClampValue(const T &value, const T &min, const T &max)
Clamp some value against a range, return the result.
Definition vtkMath.h:1964
static vtkTypeBool SolveLeastSquares(int numberOfSamples, double **xt, int xOrder, double **yt, int yOrder, double **mt, int checkHomogeneous=1)
Solves for the least squares best fit matrix for the equation X'M' = Y'.
static void Identity3x3(double A[3][3])
Set A to the identity matrix.
static void LUFactor3x3(double A[3][3], int index[3])
LU Factorization of a 3x3 matrix.
static vtkTypeBool LUFactorLinearSystem(double **A, int *index, int size)
Factor linear equations Ax = b using LU decomposition into the form A = LU where L is a unit lower tr...
static void RGBToLab(double red, double green, double blue, double *L, double *a, double *b)
Convert color from the RGB system to CIE-L*ab.
static void MultiplyQuaternion(const float q1[4], const float q2[4], float q[4])
Multiply two quaternions.
static double GaussianWeight(double variance, double distanceFromMean)
Compute the amplitude of an unnormalized Gaussian function with mean=0 and specified variance.
static void ClampValues(const double *values, int nb_values, const double range[2], double *clamped_values)
Clamp some values against a range The method without 'clamped_values' will perform in-place clamping.
static void Transpose3x3(const float A[3][3], float AT[3][3])
Transpose a 3x3 matrix.
static vtkMatrixUtilities::ScalarTypeExtractor< VectorT >::value_type SquaredNorm(VectorT &&x)
Computes the dot product between 2 vectors x and y.
Definition vtkMath.h:934
static void ClampValues(double *values, int nb_values, const double range[2])
Clamp some values against a range The method without 'clamped_values' will perform in-place clamping.
static int QuadraticRoot(double a, double b, double c, double min, double max, double *u)
find roots of ax^2+bx+c=0 in the interval min,max.
static void MultiplyMatrixWithVector(MatrixT &&M, VectorT1 &&X, VectorT2 &&Y)
Multiply matrix M with vector Y such that Y = M x X.
Definition vtkMath.h:885
Park and Miller Sequence of pseudo random numbers.
represent and manipulate 3D points
Definition vtkPoints.h:30
Computes the portion of a dataset which is inside a selection.
Hold a reference to a vtkObjectBase instance.
void RoundDoubleToIntegralIfNecessary(double val, OutT *ret)
Definition vtkMath.h:2205
typename detail::ScalarTypeExtractor< std::is_array< DerefContainer >::value||std::is_pointer< DerefContainer >::value, ContainerT >::value_type value_type
value_type is the underlying arithmetic type held in ContainerT
Template defining traits of native types used by VTK.
int vtkTypeBool
Definition vtkABI.h:64
#define vtkDataArray
double vtkDeterminant3x3(const T A[3][3])
Definition vtkMath.h:1944
int vtkIdType
Definition vtkType.h:315
#define max(a, b)