VTK  9.4.20241223
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
133#ifndef vtkMath_h
134#define vtkMath_h
135
136#include "vtkCommonCoreModule.h" // For export macro
137#include "vtkMathPrivate.hxx" // For Matrix meta-class helpers
138#include "vtkMatrixUtilities.h" // For Matrix wrapping / mapping
139#include "vtkObject.h"
140#include "vtkSmartPointer.h" // For vtkSmartPointer.
141#include "vtkTypeTraits.h" // For type traits
142
143#include "vtkMathConfigure.h" // For <cmath> and VTK_HAS_ISNAN etc.
144
145#include <algorithm> // for std::clamp
146#include <cassert> // assert() in inline implementations.
147#include <type_traits> // for type_traits
148
149#ifndef DBL_MIN
150#define VTK_DBL_MIN 2.2250738585072014e-308
151#else // DBL_MIN
152#define VTK_DBL_MIN DBL_MIN
153#endif // DBL_MIN
154
155#ifndef DBL_EPSILON
156#define VTK_DBL_EPSILON 2.2204460492503131e-16
157#else // DBL_EPSILON
158#define VTK_DBL_EPSILON DBL_EPSILON
159#endif // DBL_EPSILON
160
161#ifndef VTK_DBL_EPSILON
162#ifndef DBL_EPSILON
163#define VTK_DBL_EPSILON 2.2204460492503131e-16
164#else // DBL_EPSILON
165#define VTK_DBL_EPSILON DBL_EPSILON
166#endif // DBL_EPSILON
167#endif // VTK_DBL_EPSILON
168
169VTK_ABI_NAMESPACE_BEGIN
170class vtkDataArray;
171class vtkPoints;
172class vtkMathInternal;
175VTK_ABI_NAMESPACE_END
176
177namespace vtk_detail
178{
179VTK_ABI_NAMESPACE_BEGIN
180// forward declaration
181template <typename OutT>
182void RoundDoubleToIntegralIfNecessary(double val, OutT* ret);
183VTK_ABI_NAMESPACE_END
184} // end namespace vtk_detail
185
186VTK_ABI_NAMESPACE_BEGIN
187class VTKCOMMONCORE_EXPORT vtkMath : public vtkObject
188{
189public:
190 static vtkMath* New();
191 vtkTypeMacro(vtkMath, vtkObject);
192 void PrintSelf(ostream& os, vtkIndent indent) override;
193
194private:
195 template <class VectorT, class = void>
196 struct VectorImplementsSize : std::false_type
197 {
198 };
199
200 template <class VectorT>
201 struct VectorImplementsSize<VectorT, decltype((void)std::declval<VectorT>().size(), void())>
202 : std::true_type
203 {
204 };
205
209 template <class VectorT>
210 using EnableIfVectorImplementsSize =
211 typename std::enable_if<VectorImplementsSize<VectorT>::value>::type;
212
213public:
222 static constexpr int DYNAMIC_VECTOR_SIZE() { return 0; }
223
227 static constexpr double Pi() { return 3.141592653589793; }
228
230
233 static float RadiansFromDegrees(float degrees);
234 static double RadiansFromDegrees(double degrees);
236
238
241 static float DegreesFromRadians(float radians);
242 static double DegreesFromRadians(double radians);
244
248#if 1
249 static int Round(float f) { return static_cast<int>(f + (f >= 0.0 ? 0.5 : -0.5)); }
250 static int Round(double f) { return static_cast<int>(f + (f >= 0.0 ? 0.5 : -0.5)); }
251#endif
252
257 template <typename OutT>
258 static void RoundDoubleToIntegralIfNecessary(double val, OutT* ret)
259 {
260 // Can't specialize template methods in a template class, so we move the
261 // implementations to a external namespace.
263 }
264
270 static int Floor(double x);
271
277 static int Ceil(double x);
278
284 static int CeilLog2(vtkTypeUInt64 x);
285
290 template <class T>
291 static T Min(const T& a, const T& b);
292
297 template <class T>
298 static T Max(const T& a, const T& b);
299
303 static bool IsPowerOfTwo(vtkTypeUInt64 x);
304
310 static int NearestPowerOfTwo(int x);
311
316 static vtkTypeInt64 Factorial(int N);
317
323 static vtkTypeInt64 Binomial(int m, int n);
324
336 static int* BeginCombination(int m, int n);
337
348 static int NextCombination(int m, int n, int* combination);
349
353 static void FreeCombination(int* combination);
354
370 static void RandomSeed(int s);
371
383 static int GetSeed();
384
398 static double Random();
399
412 static double Random(double min, double max);
413
426 static double Gaussian();
427
440 static double Gaussian(double mean, double std);
441
446 template <class VectorT1, class VectorT2>
447 static void Assign(const VectorT1& a, VectorT2&& b)
448 {
449 b[0] = a[0];
450 b[1] = a[1];
451 b[2] = a[2];
452 }
453
457 static void Assign(const double a[3], double b[3]) { vtkMath::Assign<>(a, b); }
458
462 static void Add(const float a[3], const float b[3], float c[3])
463 {
464 for (int i = 0; i < 3; ++i)
465 {
466 c[i] = a[i] + b[i];
467 }
468 }
469
473 static void Add(const double a[3], const double b[3], double c[3])
474 {
475 for (int i = 0; i < 3; ++i)
476 {
477 c[i] = a[i] + b[i];
478 }
479 }
480
486 template <class VectorT1, class VectorT2, class VectorT3>
487 static void Add(VectorT1&& a, VectorT2&& b, VectorT3& c)
488 {
489 for (int i = 0; i < 3; ++i)
490 {
491 c[i] = a[i] + b[i];
492 }
493 }
494
498 static void Subtract(const float a[3], const float b[3], float c[3])
499 {
500 for (int i = 0; i < 3; ++i)
501 {
502 c[i] = a[i] - b[i];
503 }
504 }
505
509 static void Subtract(const double a[3], const double b[3], double c[3])
510 {
511 for (int i = 0; i < 3; ++i)
512 {
513 c[i] = a[i] - b[i];
514 }
515 }
516
522 template <class VectorT1, class VectorT2, class VectorT3>
523 static void Subtract(const VectorT1& a, const VectorT2& b, VectorT3&& c)
524 {
525 c[0] = a[0] - b[0];
526 c[1] = a[1] - b[1];
527 c[2] = a[2] - b[2];
528 }
529
534 static void MultiplyScalar(float a[3], float s)
535 {
536 for (int i = 0; i < 3; ++i)
537 {
538 a[i] *= s;
539 }
540 }
541
546 static void MultiplyScalar2D(float a[2], float s)
547 {
548 for (int i = 0; i < 2; ++i)
549 {
550 a[i] *= s;
551 }
552 }
553
558 static void MultiplyScalar(double a[3], double s)
559 {
560 for (int i = 0; i < 3; ++i)
561 {
562 a[i] *= s;
563 }
564 }
565
570 static void MultiplyScalar2D(double a[2], double s)
571 {
572 for (int i = 0; i < 2; ++i)
573 {
574 a[i] *= s;
575 }
576 }
577
581 static float Dot(const float a[3], const float b[3])
582 {
583 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
584 }
585
589 static double Dot(const double a[3], const double b[3])
590 {
591 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
592 }
593
609 template <typename ReturnTypeT = double, typename TupleRangeT1, typename TupleRangeT2,
610 typename EnableT = typename std::conditional<!std::is_pointer<TupleRangeT1>::value &&
611 !std::is_array<TupleRangeT1>::value,
612 TupleRangeT1, TupleRangeT2>::type::value_type>
613 static ReturnTypeT Dot(const TupleRangeT1& a, const TupleRangeT2& b)
614 {
615 return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
616 }
617
621 static void Outer(const float a[3], const float b[3], float c[3][3])
622 {
623 for (int i = 0; i < 3; ++i)
624 {
625 for (int j = 0; j < 3; ++j)
626 {
627 c[i][j] = a[i] * b[j];
628 }
629 }
630 }
631
635 static void Outer(const double a[3], const double b[3], double c[3][3])
636 {
637 for (int i = 0; i < 3; ++i)
638 {
639 for (int j = 0; j < 3; ++j)
640 {
641 c[i][j] = a[i] * b[j];
642 }
643 }
644 }
645
651 template <class VectorT1, class VectorT2, class VectorT3>
652 static void Cross(VectorT1&& a, VectorT2&& b, VectorT3& c);
653
658 static void Cross(const float a[3], const float b[3], float c[3]);
659
664 static void Cross(const double a[3], const double b[3], double c[3]);
665
667
670 static float Norm(const float* x, int n);
671 static double Norm(const double* x, int n);
673
677 static float Norm(const float v[3]) { return std::sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]); }
678
682 static double Norm(const double v[3])
683 {
684 return std::sqrt(v[0] * v[0] + v[1] * v[1] + v[2] * v[2]);
685 }
686
696 template <typename ReturnTypeT = double, typename TupleRangeT>
697 static ReturnTypeT SquaredNorm(const TupleRangeT& v)
698 {
699 return v[0] * v[0] + v[1] * v[1] + v[2] * v[2];
700 }
701
706 static float Normalize(float v[3]);
707
712 static double Normalize(double v[3]);
713
715
722 static void Perpendiculars(const double v1[3], double v2[3], double v3[3], double theta);
723 static void Perpendiculars(const float v1[3], float v2[3], float v3[3], double theta);
725
727
732 static bool ProjectVector(const float a[3], const float b[3], float projection[3]);
733 static bool ProjectVector(const double a[3], const double b[3], double projection[3]);
735
737
743 static bool ProjectVector2D(const float a[2], const float b[2], float projection[2]);
744 static bool ProjectVector2D(const double a[2], const double b[2], double projection[2]);
746
762 template <typename ReturnTypeT = double, typename TupleRangeT1, typename TupleRangeT2,
763 typename EnableT = typename std::conditional<!std::is_pointer<TupleRangeT1>::value &&
764 !std::is_array<TupleRangeT1>::value,
765 TupleRangeT1, TupleRangeT2>::type::value_type>
766 static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1& p1, const TupleRangeT2& p2);
767
772 static float Distance2BetweenPoints(const float p1[3], const float p2[3]);
773
778 static double Distance2BetweenPoints(const double p1[3], const double p2[3]);
779
783 static double AngleBetweenVectors(const double v1[3], const double v2[3]);
784
789 const double v1[3], const double v2[3], const double vn[3]);
790
795 static double GaussianAmplitude(double variance, double distanceFromMean);
796
801 static double GaussianAmplitude(double mean, double variance, double position);
802
808 static double GaussianWeight(double variance, double distanceFromMean);
809
815 static double GaussianWeight(double mean, double variance, double position);
816
820 static float Dot2D(const float x[2], const float y[2]) { return x[0] * y[0] + x[1] * y[1]; }
821
825 static double Dot2D(const double x[2], const double y[2]) { return x[0] * y[0] + x[1] * y[1]; }
826
830 static void Outer2D(const float x[2], const float y[2], float A[2][2])
831 {
832 for (int i = 0; i < 2; ++i)
833 {
834 for (int j = 0; j < 2; ++j)
835 {
836 A[i][j] = x[i] * y[j];
837 }
838 }
839 }
840
844 static void Outer2D(const double x[2], const double y[2], double A[2][2])
845 {
846 for (int i = 0; i < 2; ++i)
847 {
848 for (int j = 0; j < 2; ++j)
849 {
850 A[i][j] = x[i] * y[j];
851 }
852 }
853 }
854
859 static float Norm2D(const float x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
860
865 static double Norm2D(const double x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
866
871 static float Normalize2D(float v[2]);
872
877 static double Normalize2D(double v[2]);
878
882 static float Determinant2x2(const float c1[2], const float c2[2])
883 {
884 return c1[0] * c2[1] - c2[0] * c1[1];
885 }
886
888
891 static double Determinant2x2(double a, double b, double c, double d) { return a * d - b * c; }
892 static double Determinant2x2(const double c1[2], const double c2[2])
893 {
894 return c1[0] * c2[1] - c2[0] * c1[1];
895 }
897
899
902 static void LUFactor3x3(float A[3][3], int index[3]);
903 static void LUFactor3x3(double A[3][3], int index[3]);
905
907
910 static void LUSolve3x3(const float A[3][3], const int index[3], float x[3]);
911 static void LUSolve3x3(const double A[3][3], const int index[3], double x[3]);
913
915
919 static void LinearSolve3x3(const float A[3][3], const float x[3], float y[3]);
920 static void LinearSolve3x3(const double A[3][3], const double x[3], double y[3]);
922
924
927 static void Multiply3x3(const float A[3][3], const float v[3], float u[3]);
928 static void Multiply3x3(const double A[3][3], const double v[3], double u[3]);
930
932
935 static void Multiply3x3(const float A[3][3], const float B[3][3], float C[3][3]);
936 static void Multiply3x3(const double A[3][3], const double B[3][3], double C[3][3]);
938
962 template <int RowsT, int MidDimT, int ColsT,
963 class LayoutT1 = vtkMatrixUtilities::Layout::Identity,
964 class LayoutT2 = vtkMatrixUtilities::Layout::Identity, class MatrixT1, class MatrixT2,
965 class MatrixT3>
966 static void MultiplyMatrix(MatrixT1&& M1, MatrixT2&& M2, MatrixT3&& M3)
967 {
968 vtkMathPrivate::MultiplyMatrix<RowsT, MidDimT, ColsT, LayoutT1, LayoutT2>::Compute(
969 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2), std::forward<MatrixT3>(M3));
970 }
971
992 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
993 class MatrixT, class VectorT1, class VectorT2>
994 static void MultiplyMatrixWithVector(MatrixT&& M, VectorT1&& X, VectorT2&& Y)
995 {
996 vtkMathPrivate::MultiplyMatrix<RowsT, ColsT, 1, LayoutT>::Compute(
997 std::forward<MatrixT>(M), std::forward<VectorT1>(X), std::forward<VectorT2>(Y));
998 }
999
1005 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
1006 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
1007 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
1008 {
1009 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
1010 vtkMatrixUtilities::Layout::Identity,
1011 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x),
1012 std::forward<VectorT2>(y));
1013 }
1014
1021 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
1022 class = typename std::enable_if<SizeT == DYNAMIC_VECTOR_SIZE()>::type,
1023 class = EnableIfVectorImplementsSize<VectorT1>>
1024 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
1025 {
1026 ScalarT dot = 0.0;
1027 using SizeType = decltype(std::declval<VectorT1>().size());
1028 for (SizeType dim = 0; dim < x.size(); ++dim)
1029 {
1030 dot += x[dim] * y[dim];
1031 }
1032 return dot;
1033 }
1034
1042 template <int SizeT, class VectorT>
1044 VectorT&& x)
1045 {
1047 return vtkMath::Dot<Scalar, SizeT>(std::forward<VectorT>(x), std::forward<VectorT>(x));
1048 }
1049
1066 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT>
1068 MatrixT&& M)
1069 {
1070 return vtkMathPrivate::Determinant<SizeT, LayoutT>::Compute(std::forward<MatrixT>(M));
1071 }
1072
1088 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT1,
1089 class MatrixT2>
1090 static void InvertMatrix(MatrixT1&& M1, MatrixT2&& M2)
1091 {
1092 vtkMathPrivate::InvertMatrix<SizeT, LayoutT>::Compute(
1093 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2));
1094 }
1095
1109 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1110 class MatrixT, class VectorT1, class VectorT2>
1111 static void LinearSolve(MatrixT&& M, VectorT1&& x, VectorT2&& y)
1112 {
1113 vtkMathPrivate::LinearSolve<RowsT, ColsT, LayoutT>::Compute(
1114 std::forward<MatrixT>(M), std::forward<VectorT1>(x), std::forward<VectorT2>(y));
1115 }
1116
1131 template <class ScalarT, int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1132 class VectorT1, class MatrixT, class VectorT2,
1133 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
1134 static ScalarT Dot(VectorT1&& x, MatrixT&& M, VectorT2&& y)
1135 {
1136 ScalarT tmp[SizeT];
1137 vtkMathPrivate::MultiplyMatrix<SizeT, SizeT, 1, LayoutT>::Compute(
1138 std::forward<MatrixT>(M), std::forward<VectorT2>(y), tmp);
1139 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
1140 vtkMatrixUtilities::Layout::Identity,
1141 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x), tmp);
1142 }
1143
1149 static void MultiplyMatrix(const double* const* A, const double* const* B, unsigned int rowA,
1150 unsigned int colA, unsigned int rowB, unsigned int colB, double** C);
1151
1153
1157 static void Transpose3x3(const float A[3][3], float AT[3][3]);
1158 static void Transpose3x3(const double A[3][3], double AT[3][3]);
1160
1162
1166 static void Invert3x3(const float A[3][3], float AI[3][3]);
1167 static void Invert3x3(const double A[3][3], double AI[3][3]);
1169
1171
1174 static void Identity3x3(float A[3][3]);
1175 static void Identity3x3(double A[3][3]);
1177
1179
1182 static double Determinant3x3(const float A[3][3]);
1183 static double Determinant3x3(const double A[3][3]);
1185
1189 static float Determinant3x3(const float c1[3], const float c2[3], const float c3[3]);
1190
1194 static double Determinant3x3(const double c1[3], const double c2[3], const double c3[3]);
1195
1202 static double Determinant3x3(double a1, double a2, double a3, double b1, double b2, double b3,
1203 double c1, double c2, double c3);
1204
1206
1213 static void QuaternionToMatrix3x3(const float quat[4], float A[3][3]);
1214 static void QuaternionToMatrix3x3(const double quat[4], double A[3][3]);
1215 template <class QuaternionT, class MatrixT,
1216 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1217 static void QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A);
1219
1221
1230 static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4]);
1231 static void Matrix3x3ToQuaternion(const double A[3][3], double quat[4]);
1232 template <class MatrixT, class QuaternionT,
1233 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1234 static void Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q);
1236
1238
1244 static void MultiplyQuaternion(const float q1[4], const float q2[4], float q[4]);
1245 static void MultiplyQuaternion(const double q1[4], const double q2[4], double q[4]);
1247
1249
1253 static void RotateVectorByNormalizedQuaternion(const float v[3], const float q[4], float r[3]);
1254 static void RotateVectorByNormalizedQuaternion(const double v[3], const double q[4], double r[3]);
1256
1258
1262 static void RotateVectorByWXYZ(const float v[3], const float q[4], float r[3]);
1263 static void RotateVectorByWXYZ(const double v[3], const double q[4], double r[3]);
1265
1267
1272 static void Orthogonalize3x3(const float A[3][3], float B[3][3]);
1273 static void Orthogonalize3x3(const double A[3][3], double B[3][3]);
1275
1277
1283 static void Diagonalize3x3(const float A[3][3], float w[3], float V[3][3]);
1284 static void Diagonalize3x3(const double A[3][3], double w[3], double V[3][3]);
1286
1288
1298 const float A[3][3], float U[3][3], float w[3], float VT[3][3]);
1300 const double A[3][3], double U[3][3], double w[3], double VT[3][3]);
1302
1311 double a00, double a01, double a10, double a11, double b0, double b1, double& x0, double& x1);
1312
1321 static vtkTypeBool SolveLinearSystem(double** A, double* x, int size);
1322
1329 static vtkTypeBool InvertMatrix(double** A, double** AI, int size);
1330
1337 double** A, double** AI, int size, int* tmp1Size, double* tmp2Size);
1338
1361 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size);
1362
1368 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size, double* tmpSize);
1369
1378 static void LUSolveLinearSystem(double** A, int* index, double* x, int size);
1379
1388 static double EstimateMatrixCondition(const double* const* A, int size);
1389
1391
1399 static vtkTypeBool Jacobi(float** a, float* w, float** v);
1400 static vtkTypeBool Jacobi(double** a, double* w, double** v);
1402
1404
1413 static vtkTypeBool JacobiN(float** a, int n, float* w, float** v);
1414 static vtkTypeBool JacobiN(double** a, int n, double* w, double** v);
1416
1431 int numberOfSamples, double** xt, int xOrder, double** mt);
1432
1447 static vtkTypeBool SolveLeastSquares(int numberOfSamples, double** xt, int xOrder, double** yt,
1448 int yOrder, double** mt, int checkHomogeneous = 1);
1449
1451
1458 static void RGBToHSV(const float rgb[3], float hsv[3])
1459 {
1460 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1461 }
1462 static void RGBToHSV(float r, float g, float b, float* h, float* s, float* v);
1463 static void RGBToHSV(const double rgb[3], double hsv[3])
1464 {
1465 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1466 }
1467 static void RGBToHSV(double r, double g, double b, double* h, double* s, double* v);
1469
1471
1478 static void HSVToRGB(const float hsv[3], float rgb[3])
1479 {
1480 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1481 }
1482 static void HSVToRGB(float h, float s, float v, float* r, float* g, float* b);
1483 static void HSVToRGB(const double hsv[3], double rgb[3])
1484 {
1485 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1486 }
1487 static void HSVToRGB(double h, double s, double v, double* r, double* g, double* b);
1489
1491
1494 static void LabToXYZ(const double lab[3], double xyz[3])
1495 {
1496 LabToXYZ(lab[0], lab[1], lab[2], xyz + 0, xyz + 1, xyz + 2);
1497 }
1498 static void LabToXYZ(double L, double a, double b, double* x, double* y, double* z);
1500
1502
1505 static void XYZToLab(const double xyz[3], double lab[3])
1506 {
1507 XYZToLab(xyz[0], xyz[1], xyz[2], lab + 0, lab + 1, lab + 2);
1508 }
1509 static void XYZToLab(double x, double y, double z, double* L, double* a, double* b);
1511
1513
1516 static void XYZToRGB(const double xyz[3], double rgb[3])
1517 {
1518 XYZToRGB(xyz[0], xyz[1], xyz[2], rgb + 0, rgb + 1, rgb + 2);
1519 }
1520 static void XYZToRGB(double x, double y, double z, double* r, double* g, double* b);
1522
1524
1527 static void RGBToXYZ(const double rgb[3], double xyz[3])
1528 {
1529 RGBToXYZ(rgb[0], rgb[1], rgb[2], xyz + 0, xyz + 1, xyz + 2);
1530 }
1531 static void RGBToXYZ(double r, double g, double b, double* x, double* y, double* z);
1533
1535
1541 static void RGBToLab(const double rgb[3], double lab[3])
1542 {
1543 RGBToLab(rgb[0], rgb[1], rgb[2], lab + 0, lab + 1, lab + 2);
1544 }
1545 static void RGBToLab(double red, double green, double blue, double* L, double* a, double* b);
1547
1549
1552 static void LabToRGB(const double lab[3], double rgb[3])
1553 {
1554 LabToRGB(lab[0], lab[1], lab[2], rgb + 0, rgb + 1, rgb + 2);
1555 }
1556 static void LabToRGB(double L, double a, double b, double* red, double* green, double* blue);
1558
1560
1563 static void UninitializeBounds(double bounds[6])
1564 {
1565 bounds[0] = 1.0;
1566 bounds[1] = -1.0;
1567 bounds[2] = 1.0;
1568 bounds[3] = -1.0;
1569 bounds[4] = 1.0;
1570 bounds[5] = -1.0;
1571 }
1573
1575
1578 static vtkTypeBool AreBoundsInitialized(const double bounds[6])
1579 {
1580 if (bounds[1] - bounds[0] < 0.0)
1581 {
1582 return 0;
1583 }
1584 return 1;
1585 }
1587
1592 template <class T>
1593 static T ClampValue(const T& value, const T& min, const T& max);
1594
1596
1600 static void ClampValue(double* value, const double range[2]);
1601 static void ClampValue(double value, const double range[2], double* clamped_value);
1602 static void ClampValues(double* values, int nb_values, const double range[2]);
1603 static void ClampValues(
1604 const double* values, int nb_values, const double range[2], double* clamped_values);
1606
1613 static double ClampAndNormalizeValue(double value, const double range[2]);
1614
1619 template <class T1, class T2>
1620 static void TensorFromSymmetricTensor(const T1 symmTensor[6], T2 tensor[9]);
1621
1627 template <class T>
1628 static void TensorFromSymmetricTensor(T tensor[9]);
1629
1639 double range_min, double range_max, double scale = 1.0, double shift = 0.0);
1640
1649 static vtkTypeBool GetAdjustedScalarRange(vtkDataArray* array, int comp, double range[2]);
1650
1655 static vtkTypeBool ExtentIsWithinOtherExtent(const int extent1[6], const int extent2[6]);
1656
1663 const double bounds1[6], const double bounds2[6], const double delta[3]);
1664
1671 const double point[3], const double bounds[6], const double delta[3]);
1672
1683 const double bounds[6], const double normal[3], const double point[3]);
1684
1694 static double Solve3PointCircle(
1695 const double p1[3], const double p2[3], const double p3[3], double center[3]);
1696
1700 static double Inf();
1701
1705 static double NegInf();
1706
1710 static double Nan();
1711
1715 static vtkTypeBool IsInf(double x);
1716
1720 static vtkTypeBool IsNan(double x);
1721
1726 static bool IsFinite(double x);
1727
1732 static int QuadraticRoot(double a, double b, double c, double min, double max, double* u);
1733
1739 static vtkIdType ComputeGCD(vtkIdType m, vtkIdType n) { return (n ? ComputeGCD(n, m % n) : m); }
1740
1745 {
1746 FULL,
1747 SAME,
1748 VALID
1749 };
1750
1773 template <class Iter1, class Iter2, class Iter3>
1774 static void Convolve1D(Iter1 beginSample, Iter1 endSample, Iter2 beginKernel, Iter2 endKernel,
1775 Iter3 beginOut, Iter3 endOut, ConvolutionMode mode = ConvolutionMode::FULL)
1776 {
1777 int sampleSize = std::distance(beginSample, endSample);
1778 int kernelSize = std::distance(beginKernel, endKernel);
1779 int outSize = std::distance(beginOut, endOut);
1780
1781 if (sampleSize <= 0 || kernelSize <= 0 || outSize <= 0)
1782 {
1783 return;
1784 }
1785
1786 int begin = 0;
1787 int end = outSize;
1788
1789 switch (mode)
1790 {
1791 case ConvolutionMode::SAME:
1792 begin = static_cast<int>(std::ceil((std::min)(sampleSize, kernelSize) / 2.0)) - 1;
1793 end = begin + (std::max)(sampleSize, kernelSize);
1794 break;
1795 case ConvolutionMode::VALID:
1796 begin = (std::min)(sampleSize, kernelSize) - 1;
1797 end = begin + std::abs(sampleSize - kernelSize) + 1;
1798 break;
1799 case ConvolutionMode::FULL:
1800 default:
1801 break;
1802 }
1803
1804 for (int i = begin; i < end; i++)
1805 {
1806 Iter3 out = beginOut + i - begin;
1807 *out = 0;
1808 for (int j = (std::max)(i - sampleSize + 1, 0); j <= (std::min)(i, kernelSize - 1); j++)
1809 {
1810 *out += *(beginSample + (i - j)) * *(beginKernel + j);
1811 }
1812 }
1813 }
1814
1819 static void GetPointAlongLine(double result[3], double p1[3], double p2[3], const double offset)
1820 {
1821 double directionVector[3] = { p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2] };
1822 vtkMath::Normalize(directionVector);
1823 result[0] = p2[0] + (offset * directionVector[0]);
1824 result[1] = p2[1] + (offset * directionVector[1]);
1825 result[2] = p2[2] + (offset * directionVector[2]);
1826 }
1827
1828protected:
1829 vtkMath() = default;
1830 ~vtkMath() override = default;
1831
1833
1834private:
1835 vtkMath(const vtkMath&) = delete;
1836 void operator=(const vtkMath&) = delete;
1837};
1838
1839//----------------------------------------------------------------------------
1840inline float vtkMath::RadiansFromDegrees(float x)
1841{
1842 return x * 0.017453292f;
1843}
1844
1845//----------------------------------------------------------------------------
1846inline double vtkMath::RadiansFromDegrees(double x)
1847{
1848 return x * 0.017453292519943295;
1849}
1850
1851//----------------------------------------------------------------------------
1852inline float vtkMath::DegreesFromRadians(float x)
1853{
1854 return x * 57.2957795131f;
1855}
1856
1857//----------------------------------------------------------------------------
1858inline double vtkMath::DegreesFromRadians(double x)
1859{
1860 return x * 57.29577951308232;
1861}
1862
1863//----------------------------------------------------------------------------
1864inline bool vtkMath::IsPowerOfTwo(vtkTypeUInt64 x)
1865{
1866 return ((x != 0) & ((x & (x - 1)) == 0));
1867}
1868
1869//----------------------------------------------------------------------------
1870// Credit goes to Peter Hart and William Lewis on comp.lang.python 1997
1872{
1873 unsigned int z = static_cast<unsigned int>(((x > 0) ? x - 1 : 0));
1874 z |= z >> 1;
1875 z |= z >> 2;
1876 z |= z >> 4;
1877 z |= z >> 8;
1878 z |= z >> 16;
1879 return static_cast<int>(z + 1);
1880}
1881
1882//----------------------------------------------------------------------------
1883// Modify the trunc() operation provided by static_cast<int>() to get floor(),
1884// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1885inline int vtkMath::Floor(double x)
1886{
1887 int i = static_cast<int>(x);
1888 return i - (i > x);
1889}
1890
1891//----------------------------------------------------------------------------
1892// Modify the trunc() operation provided by static_cast<int>() to get ceil(),
1893// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1894inline int vtkMath::Ceil(double x)
1895{
1896 int i = static_cast<int>(x);
1897 return i + (i < x);
1898}
1899
1900//----------------------------------------------------------------------------
1901template <class T>
1902inline T vtkMath::Min(const T& a, const T& b)
1903{
1904 return (b <= a ? b : a);
1905}
1906
1907//----------------------------------------------------------------------------
1908template <class T>
1909inline T vtkMath::Max(const T& a, const T& b)
1910{
1911 return (b > a ? b : a);
1912}
1913
1914//----------------------------------------------------------------------------
1915inline float vtkMath::Normalize(float v[3])
1916{
1917 float den = vtkMath::Norm(v);
1918 if (den != 0.0)
1919 {
1920 for (int i = 0; i < 3; ++i)
1921 {
1922 v[i] /= den;
1923 }
1924 }
1925 return den;
1926}
1927
1928//----------------------------------------------------------------------------
1929inline double vtkMath::Normalize(double v[3])
1930{
1931 double den = vtkMath::Norm(v);
1932 if (den != 0.0)
1933 {
1934 for (int i = 0; i < 3; ++i)
1935 {
1936 v[i] /= den;
1937 }
1938 }
1939 return den;
1940}
1941
1942//----------------------------------------------------------------------------
1943inline float vtkMath::Normalize2D(float v[2])
1944{
1945 float den = vtkMath::Norm2D(v);
1946 if (den != 0.0)
1947 {
1948 for (int i = 0; i < 2; ++i)
1949 {
1950 v[i] /= den;
1951 }
1952 }
1953 return den;
1954}
1955
1956//----------------------------------------------------------------------------
1957inline double vtkMath::Normalize2D(double v[2])
1958{
1959 double den = vtkMath::Norm2D(v);
1960 if (den != 0.0)
1961 {
1962 for (int i = 0; i < 2; ++i)
1963 {
1964 v[i] /= den;
1965 }
1966 }
1967 return den;
1968}
1969
1970//----------------------------------------------------------------------------
1971inline float vtkMath::Determinant3x3(const float c1[3], const float c2[3], const float c3[3])
1972{
1973 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
1974 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
1975}
1976
1977//----------------------------------------------------------------------------
1978inline double vtkMath::Determinant3x3(const double c1[3], const double c2[3], const double c3[3])
1979{
1980 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
1981 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
1982}
1983
1984//----------------------------------------------------------------------------
1986 double a1, double a2, double a3, double b1, double b2, double b3, double c1, double c2, double c3)
1987{
1988 return (a1 * vtkMath::Determinant2x2(b2, b3, c2, c3) -
1989 b1 * vtkMath::Determinant2x2(a2, a3, c2, c3) + c1 * vtkMath::Determinant2x2(a2, a3, b2, b3));
1990}
1991
1992//----------------------------------------------------------------------------
1993inline float vtkMath::Distance2BetweenPoints(const float p1[3], const float p2[3])
1994{
1995 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
1996 (p1[2] - p2[2]) * (p1[2] - p2[2]));
1997}
1998
1999//----------------------------------------------------------------------------
2000inline double vtkMath::Distance2BetweenPoints(const double p1[3], const double p2[3])
2001{
2002 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
2003 (p1[2] - p2[2]) * (p1[2] - p2[2]));
2004}
2005
2006//----------------------------------------------------------------------------
2007template <typename ReturnTypeT, typename TupleRangeT1, typename TupleRangeT2, typename EnableT>
2008inline ReturnTypeT vtkMath::Distance2BetweenPoints(const TupleRangeT1& p1, const TupleRangeT2& p2)
2009{
2010 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
2011 (p1[2] - p2[2]) * (p1[2] - p2[2]));
2012}
2013
2014//----------------------------------------------------------------------------
2015template <class VectorT1, class VectorT2, class VectorT3>
2016void vtkMath::Cross(VectorT1&& a, VectorT2&& b, VectorT3& c)
2017{
2019 ValueType Cx = a[1] * b[2] - a[2] * b[1];
2020 ValueType Cy = a[2] * b[0] - a[0] * b[2];
2021 ValueType Cz = a[0] * b[1] - a[1] * b[0];
2022 c[0] = Cx;
2023 c[1] = Cy;
2024 c[2] = Cz;
2025}
2026
2027//----------------------------------------------------------------------------
2028// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
2029inline void vtkMath::Cross(const float a[3], const float b[3], float c[3])
2030{
2031 float Cx = a[1] * b[2] - a[2] * b[1];
2032 float Cy = a[2] * b[0] - a[0] * b[2];
2033 float Cz = a[0] * b[1] - a[1] * b[0];
2034 c[0] = Cx;
2035 c[1] = Cy;
2036 c[2] = Cz;
2037}
2038
2039//----------------------------------------------------------------------------
2040// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
2041inline void vtkMath::Cross(const double a[3], const double b[3], double c[3])
2042{
2043 double Cx = a[1] * b[2] - a[2] * b[1];
2044 double Cy = a[2] * b[0] - a[0] * b[2];
2045 double Cz = a[0] * b[1] - a[1] * b[0];
2046 c[0] = Cx;
2047 c[1] = Cy;
2048 c[2] = Cz;
2049}
2050
2051//----------------------------------------------------------------------------
2052template <class T>
2053inline double vtkDeterminant3x3(const T A[3][3])
2054{
2055 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] -
2056 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];
2057}
2058
2059//----------------------------------------------------------------------------
2060inline double vtkMath::Determinant3x3(const float A[3][3])
2061{
2062 return vtkDeterminant3x3(A);
2063}
2064
2065//----------------------------------------------------------------------------
2066inline double vtkMath::Determinant3x3(const double A[3][3])
2067{
2068 return vtkDeterminant3x3(A);
2069}
2070
2071//----------------------------------------------------------------------------
2072template <class T>
2073inline T vtkMath::ClampValue(const T& value, const T& min, const T& max)
2074{
2075 assert("pre: valid_range" && min <= max);
2076
2077#if __cplusplus >= 201703L
2078 return std::clamp(value, min, max);
2079#else
2080 // compilers are good at optimizing the ternary operator,
2081 // use '<' since it is preferred by STL for custom types
2082 T v = (min < value ? value : min);
2083 return (v < max ? v : max);
2084#endif
2085}
2086
2087//----------------------------------------------------------------------------
2088inline void vtkMath::ClampValue(double* value, const double range[2])
2089{
2090 if (value && range)
2091 {
2092 assert("pre: valid_range" && range[0] <= range[1]);
2093
2094 *value = vtkMath::ClampValue(*value, range[0], range[1]);
2095 }
2096}
2097
2098//----------------------------------------------------------------------------
2099inline void vtkMath::ClampValue(double value, const double range[2], double* clamped_value)
2100{
2101 if (range && clamped_value)
2102 {
2103 assert("pre: valid_range" && range[0] <= range[1]);
2104
2105 *clamped_value = vtkMath::ClampValue(value, range[0], range[1]);
2106 }
2107}
2108
2109// ---------------------------------------------------------------------------
2110inline double vtkMath::ClampAndNormalizeValue(double value, const double range[2])
2111{
2112 assert("pre: valid_range" && range[0] <= range[1]);
2113
2114 double result;
2115 if (range[0] == range[1])
2116 {
2117 result = 0.0;
2118 }
2119 else
2120 {
2121 // clamp
2122 result = vtkMath::ClampValue(value, range[0], range[1]);
2123
2124 // normalize
2125 result = (result - range[0]) / (range[1] - range[0]);
2126 }
2127
2128 assert("post: valid_result" && result >= 0.0 && result <= 1.0);
2129
2130 return result;
2131}
2132
2133//-----------------------------------------------------------------------------
2134template <class T1, class T2>
2135inline void vtkMath::TensorFromSymmetricTensor(const T1 symmTensor[9], T2 tensor[9])
2136{
2137 for (int i = 0; i < 3; ++i)
2138 {
2139 tensor[4 * i] = symmTensor[i];
2140 }
2141 tensor[1] = tensor[3] = symmTensor[3];
2142 tensor[2] = tensor[6] = symmTensor[5];
2143 tensor[5] = tensor[7] = symmTensor[4];
2144}
2145
2146//-----------------------------------------------------------------------------
2147template <class T>
2149{
2150 tensor[6] = tensor[5]; // XZ
2151 tensor[7] = tensor[4]; // YZ
2152 tensor[8] = tensor[2]; // ZZ
2153 tensor[4] = tensor[1]; // YY
2154 tensor[5] = tensor[7]; // YZ
2155 tensor[2] = tensor[6]; // XZ
2156 tensor[1] = tensor[3]; // XY
2157}
2158VTK_ABI_NAMESPACE_END
2159
2160namespace
2161{
2162template <class QuaternionT, class MatrixT>
2163inline void vtkQuaternionToMatrix3x3(QuaternionT&& quat, MatrixT&& A)
2164{
2166
2167 Scalar ww = quat[0] * quat[0];
2168 Scalar wx = quat[0] * quat[1];
2169 Scalar wy = quat[0] * quat[2];
2170 Scalar wz = quat[0] * quat[3];
2171
2172 Scalar xx = quat[1] * quat[1];
2173 Scalar yy = quat[2] * quat[2];
2174 Scalar zz = quat[3] * quat[3];
2175
2176 Scalar xy = quat[1] * quat[2];
2177 Scalar xz = quat[1] * quat[3];
2178 Scalar yz = quat[2] * quat[3];
2179
2180 Scalar rr = xx + yy + zz;
2181 // normalization factor, just in case quaternion was not normalized
2182 Scalar f = 1 / (ww + rr);
2183 Scalar s = (ww - rr) * f;
2184 f *= 2;
2185
2187
2188 Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) = xx * f + s;
2189 Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) = (xy + wz) * f;
2190 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A)) = (xz - wy) * f;
2191
2192 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A)) = (xy - wz) * f;
2193 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) = yy * f + s;
2194 Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) = (yz + wx) * f;
2195
2196 Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) = (xz + wy) * f;
2197 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A)) = (yz - wx) * f;
2198 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A)) = zz * f + s;
2199}
2200} // anonymous namespace
2201
2202VTK_ABI_NAMESPACE_BEGIN
2203//------------------------------------------------------------------------------
2204inline void vtkMath::QuaternionToMatrix3x3(const float quat[4], float A[3][3])
2205{
2206 vtkQuaternionToMatrix3x3(quat, A);
2207}
2208
2209//------------------------------------------------------------------------------
2210inline void vtkMath::QuaternionToMatrix3x3(const double quat[4], double A[3][3])
2211{
2212 vtkQuaternionToMatrix3x3(quat, A);
2213}
2214
2215//-----------------------------------------------------------------------------
2216template <class QuaternionT, class MatrixT, class EnableT>
2217inline void vtkMath::QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A)
2218{
2219 vtkQuaternionToMatrix3x3(std::forward<QuaternionT>(q), std::forward<MatrixT>(A));
2220}
2221VTK_ABI_NAMESPACE_END
2222
2223namespace
2224{
2225//------------------------------------------------------------------------------
2226// The solution is based on
2227// Berthold K. P. Horn (1987),
2228// "Closed-form solution of absolute orientation using unit quaternions,"
2229// Journal of the Optical Society of America A, 4:629-642
2230template <class MatrixT, class QuaternionT>
2231inline void vtkMatrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& quat)
2232{
2234
2235 Scalar N[4][4];
2236
2238
2239 // on-diagonal elements
2240 N[0][0] = Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) +
2241 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) +
2242 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2243 N[1][1] = Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) -
2244 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) -
2245 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2246 N[2][2] = -Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) +
2247 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) -
2248 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2249 N[3][3] = -Wrapper::template Get<0, 0>(std::forward<MatrixT>(A)) -
2250 Wrapper::template Get<1, 1>(std::forward<MatrixT>(A)) +
2251 Wrapper::template Get<2, 2>(std::forward<MatrixT>(A));
2252
2253 // off-diagonal elements
2254 N[0][1] = N[1][0] = Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) -
2255 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A));
2256 N[0][2] = N[2][0] = Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) -
2257 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A));
2258 N[0][3] = N[3][0] = Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) -
2259 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A));
2260
2261 N[1][2] = N[2][1] = Wrapper::template Get<1, 0>(std::forward<MatrixT>(A)) +
2262 Wrapper::template Get<0, 1>(std::forward<MatrixT>(A));
2263 N[1][3] = N[3][1] = Wrapper::template Get<0, 2>(std::forward<MatrixT>(A)) +
2264 Wrapper::template Get<2, 0>(std::forward<MatrixT>(A));
2265 N[2][3] = N[3][2] = Wrapper::template Get<2, 1>(std::forward<MatrixT>(A)) +
2266 Wrapper::template Get<1, 2>(std::forward<MatrixT>(A));
2267
2268 Scalar eigenvectors[4][4], eigenvalues[4];
2269
2270 // convert into format that JacobiN can use,
2271 // then use Jacobi to find eigenvalues and eigenvectors
2272 Scalar *NTemp[4], *eigenvectorsTemp[4];
2273 for (int i = 0; i < 4; ++i)
2274 {
2275 NTemp[i] = N[i];
2276 eigenvectorsTemp[i] = eigenvectors[i];
2277 }
2278 vtkMath::JacobiN(NTemp, 4, eigenvalues, eigenvectorsTemp);
2279
2280 // the first eigenvector is the one we want
2281 quat[0] = eigenvectors[0][0];
2282 quat[1] = eigenvectors[1][0];
2283 quat[2] = eigenvectors[2][0];
2284 quat[3] = eigenvectors[3][0];
2285}
2286} // anonymous namespace
2287
2288VTK_ABI_NAMESPACE_BEGIN
2289//------------------------------------------------------------------------------
2290inline void vtkMath::Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
2291{
2292 vtkMatrix3x3ToQuaternion(A, quat);
2293}
2294
2295//------------------------------------------------------------------------------
2296inline void vtkMath::Matrix3x3ToQuaternion(const double A[3][3], double quat[4])
2297{
2298 vtkMatrix3x3ToQuaternion(A, quat);
2299}
2300
2301//-----------------------------------------------------------------------------
2302template <class MatrixT, class QuaternionT, class EnableT>
2303inline void vtkMath::Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q)
2304{
2305 vtkMatrix3x3ToQuaternion(std::forward<MatrixT>(A), std::forward<QuaternionT>(q));
2306}
2307VTK_ABI_NAMESPACE_END
2308
2309namespace vtk_detail
2310{
2311VTK_ABI_NAMESPACE_BEGIN
2312// Can't specialize templates inside a template class, so we move the impl here.
2313template <typename OutT>
2314void RoundDoubleToIntegralIfNecessary(double val, OutT* ret)
2315{ // OutT is integral -- clamp and round
2316 if (!vtkMath::IsNan(val))
2317 {
2318 double min = static_cast<double>(vtkTypeTraits<OutT>::Min());
2319 double max = static_cast<double>(vtkTypeTraits<OutT>::Max());
2320 val = vtkMath::ClampValue(val, min, max);
2321 *ret = static_cast<OutT>((val >= 0.0) ? (val + 0.5) : (val - 0.5));
2322 }
2323 else
2324 *ret = 0;
2325}
2326template <>
2327inline void RoundDoubleToIntegralIfNecessary(double val, double* retVal)
2328{ // OutT is double: passthrough
2329 *retVal = val;
2330}
2331template <>
2332inline void RoundDoubleToIntegralIfNecessary(double val, float* retVal)
2333{ // OutT is float -- just clamp (as doubles, then the cast to float is well-defined.)
2334 if (!vtkMath::IsNan(val))
2335 {
2336 double min = static_cast<double>(vtkTypeTraits<float>::Min());
2337 double max = static_cast<double>(vtkTypeTraits<float>::Max());
2338 val = vtkMath::ClampValue(val, min, max);
2339 }
2340
2341 *retVal = static_cast<float>(val);
2342}
2343VTK_ABI_NAMESPACE_END
2344} // end namespace vtk_detail
2345
2346VTK_ABI_NAMESPACE_BEGIN
2347//-----------------------------------------------------------------------------
2348#if defined(VTK_HAS_ISINF) || defined(VTK_HAS_STD_ISINF)
2349#define VTK_MATH_ISINF_IS_INLINE
2350inline vtkTypeBool vtkMath::IsInf(double x)
2351{
2352#if defined(VTK_HAS_STD_ISINF)
2353 return std::isinf(x);
2354#else
2355 return (isinf(x) != 0); // Force conversion to bool
2356#endif
2357}
2358#endif
2359
2360//-----------------------------------------------------------------------------
2361#if defined(VTK_HAS_ISNAN) || defined(VTK_HAS_STD_ISNAN)
2362#define VTK_MATH_ISNAN_IS_INLINE
2363inline vtkTypeBool vtkMath::IsNan(double x)
2364{
2365#if defined(VTK_HAS_STD_ISNAN)
2366 return std::isnan(x);
2367#else
2368 return (isnan(x) != 0); // Force conversion to bool
2369#endif
2370}
2371#endif
2372
2373//-----------------------------------------------------------------------------
2374#if defined(VTK_HAS_ISFINITE) || defined(VTK_HAS_STD_ISFINITE) || defined(VTK_HAS_FINITE)
2375#define VTK_MATH_ISFINITE_IS_INLINE
2376inline bool vtkMath::IsFinite(double x)
2377{
2378#if defined(VTK_HAS_STD_ISFINITE)
2379 return std::isfinite(x);
2380#elif defined(VTK_HAS_ISFINITE)
2381 return (isfinite(x) != 0); // Force conversion to bool
2382#else
2383 return (finite(x) != 0); // Force conversion to bool
2384#endif
2385}
2386#endif
2387
2388VTK_ABI_NAMESPACE_END
2389#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.
abstract superclass for arrays of numeric data
a simple class to control print indentation
Definition vtkIndent.h:108
performs common math operations
Definition vtkMath.h:188
static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1 &p1, const TupleRangeT2 &p2)
Compute distance squared between two points p1 and p2.
Definition vtkMath.h:2008
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:589
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:249
static vtkIdType ComputeGCD(vtkIdType m, vtkIdType n)
Compute the greatest common divisor (GCD) of two positive integers m and n.
Definition vtkMath.h:1739
static double Norm2D(const double x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:865
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:1819
static void Subtract(const float a[3], const float b[3], float c[3])
Subtraction of two 3-vectors (float version).
Definition vtkMath.h:498
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:830
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:1832
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:473
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:677
static ReturnTypeT Dot(const TupleRangeT1 &a, const TupleRangeT2 &b)
Compute dot product between two points p1 and p2.
Definition vtkMath.h:613
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:1024
static void XYZToLab(const double xyz[3], double lab[3])
Convert Color from the CIE XYZ system to CIE-L*ab.
Definition vtkMath.h:1505
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! = N*(N-1) * (N-2)...*3*2*1.
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:1111
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:1494
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:581
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:546
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:1478
static void Assign(const double a[3], double b[3])
Assign values to a 3-vector (double version).
Definition vtkMath.h:457
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:892
static T Max(const T &a, const T &b)
Returns the maximum of the two arguments provided.
Definition vtkMath.h:1909
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:844
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:570
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:1894
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:1483
~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:1007
static void RGBToXYZ(const double rgb[3], double xyz[3])
Convert color from the RGB system to CIE XYZ.
Definition vtkMath.h:1527
static void QuaternionToMatrix3x3(const float quat[4], float A[3][3])
Convert a quaternion to a 3x3 rotation matrix.
Definition vtkMath.h:2204
static int NearestPowerOfTwo(int x)
Compute the nearest power of two that is not less than x.
Definition vtkMath.h:1871
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:1134
static float Normalize2D(float v[2])
Normalize (in place) a 2-vector.
Definition vtkMath.h:1943
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:222
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:635
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:447
static float RadiansFromDegrees(float degrees)
Convert degrees into radians.
Definition vtkMath.h:1840
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:1774
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:462
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:1578
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:227
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:509
static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
Convert a 3x3 matrix into a quaternion.
Definition vtkMath.h:2290
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:1516
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:2110
static void MultiplyScalar(double a[3], double s)
Multiplies a 3-vector by a scalar (double version).
Definition vtkMath.h:558
static double Dot2D(const double x[2], const double y[2])
Dot product of two 2-vectors.
Definition vtkMath.h:825
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:1541
static int Floor(double x)
Rounds a double to the nearest integer not greater than itself.
Definition vtkMath.h:1885
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:523
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:891
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:1463
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:1852
static float Determinant2x2(const float c1[2], const float c2[2])
Compute determinant of 2x2 matrix.
Definition vtkMath.h:882
static int Round(double f)
Definition vtkMath.h:250
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:1563
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:621
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:682
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:258
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:1864
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:1915
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:697
static double Determinant3x3(const float A[3][3])
Return the determinant of a 3x3 matrix.
Definition vtkMath.h:2060
static float Dot2D(const float x[2], const float y[2])
Dot product of two 2-vectors.
Definition vtkMath.h:820
ConvolutionMode
Support the convolution operations.
Definition vtkMath.h:1745
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:1458
static void Add(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Addition of two 3-vectors (double version).
Definition vtkMath.h:487
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:1067
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:1552
static float Norm2D(const float x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:859
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:534
static T Min(const T &a, const T &b)
Returns the minimum of the two arguments provided.
Definition vtkMath.h:1902
static void InvertMatrix(MatrixT1 &&M1, MatrixT2 &&M2)
Computes the inverse of input matrix M1 into M2.
Definition vtkMath.h:1090
static void Cross(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Cross product of two 3-vectors.
Definition vtkMath.h:2016
static void MultiplyMatrix(MatrixT1 &&M1, MatrixT2 &&M2, MatrixT3 &&M3)
Multiply matrices such that M3 = M1 x M2.
Definition vtkMath.h:966
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:2073
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:1043
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:994
Park and Miller Sequence of pseudo random numbers.
abstract base class for most VTK objects
Definition vtkObject.h:162
represent and manipulate 3D points
Definition vtkPoints.h:139
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:2314
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
double vtkDeterminant3x3(const T A[3][3])
Definition vtkMath.h:2053
int vtkIdType
Definition vtkType.h:315
#define max(a, b)