VTK  9.6.20260222
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
132
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 inline float Normalize(float v[3]);
707
712 static inline 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
784 static double Distance2BetweenPoints2D(const double p1[2], const double p2[2]);
785
789 static double AngleBetweenVectors(const double v1[3], const double v2[3]);
790
795 const double v1[3], const double v2[3], const double vn[3]);
796
801 static double GaussianAmplitude(double variance, double distanceFromMean);
802
807 static double GaussianAmplitude(double mean, double variance, double position);
808
814 static double GaussianWeight(double variance, double distanceFromMean);
815
821 static double GaussianWeight(double mean, double variance, double position);
822
826 static float Dot2D(const float x[2], const float y[2]) { return x[0] * y[0] + x[1] * y[1]; }
827
831 static double Dot2D(const double x[2], const double y[2]) { return x[0] * y[0] + x[1] * y[1]; }
832
836 static void Outer2D(const float x[2], const float y[2], float A[2][2])
837 {
838 for (int i = 0; i < 2; ++i)
839 {
840 for (int j = 0; j < 2; ++j)
841 {
842 A[i][j] = x[i] * y[j];
843 }
844 }
845 }
846
850 static void Outer2D(const double x[2], const double y[2], double A[2][2])
851 {
852 for (int i = 0; i < 2; ++i)
853 {
854 for (int j = 0; j < 2; ++j)
855 {
856 A[i][j] = x[i] * y[j];
857 }
858 }
859 }
860
865 static float Norm2D(const float x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
866
871 static double Norm2D(const double x[2]) { return std::sqrt(x[0] * x[0] + x[1] * x[1]); }
872
877 static float Normalize2D(float v[2]);
878
883 static double Normalize2D(double v[2]);
884
888 static float Determinant2x2(const float c1[2], const float c2[2])
889 {
890 return c1[0] * c2[1] - c2[0] * c1[1];
891 }
892
894
897 static double Determinant2x2(double a, double b, double c, double d) { return a * d - b * c; }
898 static double Determinant2x2(const double c1[2], const double c2[2])
899 {
900 return c1[0] * c2[1] - c2[0] * c1[1];
901 }
902
903
905
908 static void LUFactor3x3(float A[3][3], int index[3]);
909 static void LUFactor3x3(double A[3][3], int index[3]);
911
913
916 static void LUSolve3x3(const float A[3][3], const int index[3], float x[3]);
917 static void LUSolve3x3(const double A[3][3], const int index[3], double x[3]);
919
921
925 static void LinearSolve3x3(const float A[3][3], const float x[3], float y[3]);
926 static void LinearSolve3x3(const double A[3][3], const double x[3], double y[3]);
928
930
933 static void Multiply3x3(const float A[3][3], const float v[3], float u[3]);
934 static void Multiply3x3(const double A[3][3], const double v[3], double u[3]);
936
938
941 static void Multiply3x3(const float A[3][3], const float B[3][3], float C[3][3]);
942 static void Multiply3x3(const double A[3][3], const double B[3][3], double C[3][3]);
944
968 template <int RowsT, int MidDimT, int ColsT,
969 class LayoutT1 = vtkMatrixUtilities::Layout::Identity,
970 class LayoutT2 = vtkMatrixUtilities::Layout::Identity, class MatrixT1, class MatrixT2,
971 class MatrixT3>
972 static void MultiplyMatrix(MatrixT1&& M1, MatrixT2&& M2, MatrixT3&& M3)
973 {
974 vtkMathPrivate::MultiplyMatrix<RowsT, MidDimT, ColsT, LayoutT1, LayoutT2>::Compute(
975 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2), std::forward<MatrixT3>(M3));
976 }
977
998 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
999 class MatrixT, class VectorT1, class VectorT2>
1000 static void MultiplyMatrixWithVector(MatrixT&& M, VectorT1&& X, VectorT2&& Y)
1001 {
1002 vtkMathPrivate::MultiplyMatrix<RowsT, ColsT, 1, LayoutT>::Compute(
1003 std::forward<MatrixT>(M), std::forward<VectorT1>(X), std::forward<VectorT2>(Y));
1004 }
1005
1011 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
1012 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
1013 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
1014 {
1015 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
1016 vtkMatrixUtilities::Layout::Identity,
1017 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x),
1018 std::forward<VectorT2>(y));
1019 }
1020
1027 template <class ScalarT, int SizeT, class VectorT1, class VectorT2,
1028 class = typename std::enable_if<SizeT == DYNAMIC_VECTOR_SIZE()>::type,
1029 class = EnableIfVectorImplementsSize<VectorT1>>
1030 static ScalarT Dot(VectorT1&& x, VectorT2&& y)
1031 {
1032 ScalarT dot = 0.0;
1033 using SizeType = decltype(std::declval<VectorT1>().size());
1034 for (SizeType dim = 0; dim < x.size(); ++dim)
1035 {
1036 dot += x[dim] * y[dim];
1037 }
1038 return dot;
1039 }
1040
1048 template <int SizeT, class VectorT>
1050 VectorT&& x)
1051 {
1053 return vtkMath::Dot<Scalar, SizeT>(std::forward<VectorT>(x), std::forward<VectorT>(x));
1054 }
1055
1072 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT>
1074 MatrixT&& M)
1075 {
1076 return vtkMathPrivate::Determinant<SizeT, LayoutT>::Compute(std::forward<MatrixT>(M));
1077 }
1078
1094 template <int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity, class MatrixT1,
1095 class MatrixT2>
1096 static void InvertMatrix(MatrixT1&& M1, MatrixT2&& M2)
1097 {
1098 vtkMathPrivate::InvertMatrix<SizeT, LayoutT>::Compute(
1099 std::forward<MatrixT1>(M1), std::forward<MatrixT2>(M2));
1100 }
1101
1115 template <int RowsT, int ColsT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1116 class MatrixT, class VectorT1, class VectorT2>
1117 static void LinearSolve(MatrixT&& M, VectorT1&& x, VectorT2&& y)
1118 {
1119 vtkMathPrivate::LinearSolve<RowsT, ColsT, LayoutT>::Compute(
1120 std::forward<MatrixT>(M), std::forward<VectorT1>(x), std::forward<VectorT2>(y));
1121 }
1122
1137 template <class ScalarT, int SizeT, class LayoutT = vtkMatrixUtilities::Layout::Identity,
1138 class VectorT1, class MatrixT, class VectorT2,
1139 class = typename std::enable_if<SizeT != DYNAMIC_VECTOR_SIZE()>::type>
1140 static ScalarT Dot(VectorT1&& x, MatrixT&& M, VectorT2&& y)
1141 {
1142 ScalarT tmp[SizeT];
1143 vtkMathPrivate::MultiplyMatrix<SizeT, SizeT, 1, LayoutT>::Compute(
1144 std::forward<MatrixT>(M), std::forward<VectorT2>(y), tmp);
1145 return vtkMathPrivate::ContractRowWithCol<ScalarT, 1, SizeT, 1, 0, 0,
1146 vtkMatrixUtilities::Layout::Identity,
1147 vtkMatrixUtilities::Layout::Transpose>::Compute(std::forward<VectorT1>(x), tmp);
1148 }
1149
1155 static void MultiplyMatrix(const double* const* A, const double* const* B, unsigned int rowA,
1156 unsigned int colA, unsigned int rowB, unsigned int colB, double** C);
1157
1159
1163 static void Transpose3x3(const float A[3][3], float AT[3][3]);
1164 static void Transpose3x3(const double A[3][3], double AT[3][3]);
1166
1168
1172 static void Invert3x3(const float A[3][3], float AI[3][3]);
1173 static void Invert3x3(const double A[3][3], double AI[3][3]);
1175
1177
1180 static void Identity3x3(float A[3][3]);
1181 static void Identity3x3(double A[3][3]);
1183
1185
1188 static double Determinant3x3(const float A[3][3]);
1189 static double Determinant3x3(const double A[3][3]);
1191
1195 static float Determinant3x3(const float c1[3], const float c2[3], const float c3[3]);
1196
1200 static double Determinant3x3(const double c1[3], const double c2[3], const double c3[3]);
1201
1208 static double Determinant3x3(double a1, double a2, double a3, double b1, double b2, double b3,
1209 double c1, double c2, double c3);
1210
1212
1219 static void QuaternionToMatrix3x3(const float quat[4], float A[3][3]);
1220 static void QuaternionToMatrix3x3(const double quat[4], double A[3][3]);
1221 template <class QuaternionT, class MatrixT,
1222 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1223 static void QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A);
1225
1227
1236 static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4]);
1237 static void Matrix3x3ToQuaternion(const double A[3][3], double quat[4]);
1238 template <class MatrixT, class QuaternionT,
1239 class EnableT = typename std::enable_if<!vtkMatrixUtilities::MatrixIs2DArray<MatrixT>()>::type>
1240 static void Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q);
1242
1244
1250 static void MultiplyQuaternion(const float q1[4], const float q2[4], float q[4]);
1251 static void MultiplyQuaternion(const double q1[4], const double q2[4], double q[4]);
1253
1255
1259 static void RotateVectorByNormalizedQuaternion(const float v[3], const float q[4], float r[3]);
1260 static void RotateVectorByNormalizedQuaternion(const double v[3], const double q[4], double r[3]);
1262
1264
1268 static void RotateVectorByWXYZ(const float v[3], const float q[4], float r[3]);
1269 static void RotateVectorByWXYZ(const double v[3], const double q[4], double r[3]);
1271
1273
1278 static void Orthogonalize3x3(const float A[3][3], float B[3][3]);
1279 static void Orthogonalize3x3(const double A[3][3], double B[3][3]);
1281
1283
1289 static void Diagonalize3x3(const float A[3][3], float w[3], float V[3][3]);
1290 static void Diagonalize3x3(const double A[3][3], double w[3], double V[3][3]);
1292
1294
1304 const float A[3][3], float U[3][3], float w[3], float VT[3][3]);
1306 const double A[3][3], double U[3][3], double w[3], double VT[3][3]);
1308
1317 double a00, double a01, double a10, double a11, double b0, double b1, double& x0, double& x1);
1318
1327 static vtkTypeBool SolveLinearSystem(double** A, double* x, int size);
1328
1335 static vtkTypeBool InvertMatrix(double** A, double** AI, int size);
1336
1343 double** A, double** AI, int size, int* tmp1Size, double* tmp2Size);
1344
1367 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size);
1368
1374 static vtkTypeBool LUFactorLinearSystem(double** A, int* index, int size, double* tmpSize);
1375
1384 static void LUSolveLinearSystem(double** A, int* index, double* x, int size);
1385
1394 static double EstimateMatrixCondition(const double* const* A, int size);
1395
1397
1405 static vtkTypeBool Jacobi(float** a, float* w, float** v);
1406 static vtkTypeBool Jacobi(double** a, double* w, double** v);
1408
1410
1419 static vtkTypeBool JacobiN(float** a, int n, float* w, float** v);
1420 static vtkTypeBool JacobiN(double** a, int n, double* w, double** v);
1422
1437 int numberOfSamples, double** xt, int xOrder, double** mt);
1438
1453 static vtkTypeBool SolveLeastSquares(int numberOfSamples, double** xt, int xOrder, double** yt,
1454 int yOrder, double** mt, int checkHomogeneous = 1);
1455
1457
1464 static void RGBToHSV(const float rgb[3], float hsv[3])
1465 {
1466 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1467 }
1468 static void RGBToHSV(float r, float g, float b, float* h, float* s, float* v);
1469 static void RGBToHSV(const double rgb[3], double hsv[3])
1470 {
1471 RGBToHSV(rgb[0], rgb[1], rgb[2], hsv, hsv + 1, hsv + 2);
1472 }
1473 static void RGBToHSV(double r, double g, double b, double* h, double* s, double* v);
1475
1477
1484 static void HSVToRGB(const float hsv[3], float rgb[3])
1485 {
1486 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1487 }
1488 static void HSVToRGB(float h, float s, float v, float* r, float* g, float* b);
1489 static void HSVToRGB(const double hsv[3], double rgb[3])
1490 {
1491 HSVToRGB(hsv[0], hsv[1], hsv[2], rgb, rgb + 1, rgb + 2);
1492 }
1493 static void HSVToRGB(double h, double s, double v, double* r, double* g, double* b);
1495
1497
1501 static void ProLabToXYZ(const double prolab[3], double xyz[3])
1502 {
1503 ProLabToXYZ(prolab[0], prolab[1], prolab[2], xyz + 0, xyz + 1, xyz + 2);
1504 }
1505 static void ProLabToXYZ(double L, double a, double b, double* x, double* y, double* z);
1507
1509
1513 static void XYZToProLab(const double xyz[3], double prolab[3])
1514 {
1515 XYZToProLab(xyz[0], xyz[1], xyz[2], prolab + 0, prolab + 1, prolab + 2);
1516 }
1517 static void XYZToProLab(double x, double y, double z, double* L, double* a, double* b);
1519
1521
1524 static void LabToXYZ(const double lab[3], double xyz[3])
1525 {
1526 LabToXYZ(lab[0], lab[1], lab[2], xyz + 0, xyz + 1, xyz + 2);
1527 }
1528 static void LabToXYZ(double L, double a, double b, double* x, double* y, double* z);
1530
1532
1535 static void XYZToLab(const double xyz[3], double lab[3])
1536 {
1537 XYZToLab(xyz[0], xyz[1], xyz[2], lab + 0, lab + 1, lab + 2);
1538 }
1539 static void XYZToLab(double x, double y, double z, double* L, double* a, double* b);
1541
1543
1546 static void XYZToRGB(const double xyz[3], double rgb[3])
1547 {
1548 XYZToRGB(xyz[0], xyz[1], xyz[2], rgb + 0, rgb + 1, rgb + 2);
1549 }
1550 static void XYZToRGB(double x, double y, double z, double* r, double* g, double* b);
1552
1554
1557 static void RGBToXYZ(const double rgb[3], double xyz[3])
1558 {
1559 RGBToXYZ(rgb[0], rgb[1], rgb[2], xyz + 0, xyz + 1, xyz + 2);
1560 }
1561 static void RGBToXYZ(double r, double g, double b, double* x, double* y, double* z);
1563
1565
1571
1572 static void RGBToLab(const double rgb[3], double lab[3])
1573 {
1574 RGBToLab(rgb[0], rgb[1], rgb[2], lab + 0, lab + 1, lab + 2);
1575 }
1576 static void RGBToLab(double red, double green, double blue, double* L, double* a, double* b);
1578
1580
1583 static void ProLabToRGB(const double prolab[3], double rgb[3])
1584 {
1585 ProLabToRGB(prolab[0], prolab[1], prolab[2], rgb + 0, rgb + 1, rgb + 2);
1586 }
1587 static void ProLabToRGB(double L, double a, double b, double* red, double* green, double* blue);
1589
1591
1598 static void RGBToProLab(const double rgb[3], double prolab[3])
1599 {
1600 RGBToProLab(rgb[0], rgb[1], rgb[2], prolab + 0, prolab + 1, prolab + 2);
1601 }
1602 static void RGBToProLab(double red, double green, double blue, double* L, double* a, double* b);
1604
1606
1609 static void LabToRGB(const double lab[3], double rgb[3])
1610 {
1611 LabToRGB(lab[0], lab[1], lab[2], rgb + 0, rgb + 1, rgb + 2);
1612 }
1613 static void LabToRGB(double L, double a, double b, double* red, double* green, double* blue);
1615
1617
1620 static void UninitializeBounds(double bounds[6])
1621 {
1622 bounds[0] = 1.0;
1623 bounds[1] = -1.0;
1624 bounds[2] = 1.0;
1625 bounds[3] = -1.0;
1626 bounds[4] = 1.0;
1627 bounds[5] = -1.0;
1628 }
1629
1630
1632
1635 static vtkTypeBool AreBoundsInitialized(const double bounds[6])
1636 {
1637 if (bounds[1] - bounds[0] < 0.0)
1638 {
1639 return 0;
1640 }
1641 return 1;
1642 }
1643
1644
1649 template <class T>
1650 static T ClampValue(const T& value, const T& min, const T& max);
1651
1653
1657 static void ClampValue(double* value, const double range[2]);
1658 static void ClampValue(double value, const double range[2], double* clamped_value);
1659 static void ClampValues(double* values, int nb_values, const double range[2]);
1660 static void ClampValues(
1661 const double* values, int nb_values, const double range[2], double* clamped_values);
1663
1670 static double ClampAndNormalizeValue(double value, const double range[2]);
1671
1676 template <class T1, class T2>
1677 static void TensorFromSymmetricTensor(const T1 symmTensor[6], T2 tensor[9]);
1678
1684 template <class T>
1685 static void TensorFromSymmetricTensor(T tensor[9]);
1686
1696 double range_min, double range_max, double scale = 1.0, double shift = 0.0);
1697
1706 static vtkTypeBool GetAdjustedScalarRange(vtkDataArray* array, int comp, double range[2]);
1707
1712 static vtkTypeBool ExtentIsWithinOtherExtent(const int extent1[6], const int extent2[6]);
1713
1720 const double bounds1[6], const double bounds2[6], const double delta[3]);
1721
1728 const double point[3], const double bounds[6], const double delta[3]);
1729
1740 const double bounds[6], const double normal[3], const double point[3]);
1741
1751 static double Solve3PointCircle(
1752 const double p1[3], const double p2[3], const double p3[3], double center[3]);
1753
1757 static double Inf();
1758
1762 static double NegInf();
1763
1767 static double Nan();
1768
1772 static vtkTypeBool IsInf(double x);
1773
1777 static inline vtkTypeBool IsNan(double x);
1778
1783 static bool IsFinite(double x);
1784
1789 static int QuadraticRoot(double a, double b, double c, double min, double max, double* u);
1790
1796 static vtkIdType ComputeGCD(vtkIdType m, vtkIdType n) { return (n ? ComputeGCD(n, m % n) : m); }
1797
1802 {
1806 };
1807
1830 template <class Iter1, class Iter2, class Iter3>
1831 static void Convolve1D(Iter1 beginSample, Iter1 endSample, Iter2 beginKernel, Iter2 endKernel,
1832 Iter3 beginOut, Iter3 endOut, ConvolutionMode mode = ConvolutionMode::FULL)
1833 {
1834 int sampleSize = std::distance(beginSample, endSample);
1835 int kernelSize = std::distance(beginKernel, endKernel);
1836 int outSize = std::distance(beginOut, endOut);
1837
1838 if (sampleSize <= 0 || kernelSize <= 0 || outSize <= 0)
1839 {
1840 return;
1841 }
1842
1843 int begin = 0;
1844 int end = outSize;
1845
1846 switch (mode)
1847 {
1849 begin = static_cast<int>(std::ceil((std::min)(sampleSize, kernelSize) / 2.0)) - 1;
1850 end = begin + (std::max)(sampleSize, kernelSize);
1851 break;
1853 begin = (std::min)(sampleSize, kernelSize) - 1;
1854 end = begin + std::abs(sampleSize - kernelSize) + 1;
1855 break;
1857 default:
1858 break;
1859 }
1860
1861 for (int i = begin; i < end; i++)
1862 {
1863 Iter3 out = beginOut + i - begin;
1864 *out = 0;
1865 for (int j = (std::max)(i - sampleSize + 1, 0); j <= (std::min)(i, kernelSize - 1); j++)
1866 {
1867 *out += *(beginSample + (i - j)) * *(beginKernel + j);
1868 }
1869 }
1870 }
1871
1876 static void GetPointAlongLine(double result[3], double p1[3], double p2[3], const double offset)
1877 {
1878 double directionVector[3] = { p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2] };
1879 vtkMath::Normalize(directionVector);
1880 result[0] = p2[0] + (offset * directionVector[0]);
1881 result[1] = p2[1] + (offset * directionVector[1]);
1882 result[2] = p2[2] + (offset * directionVector[2]);
1883 }
1884
1885protected:
1886 vtkMath() = default;
1887 ~vtkMath() override = default;
1888
1890
1891private:
1892 vtkMath(const vtkMath&) = delete;
1893 void operator=(const vtkMath&) = delete;
1894};
1895
1896//----------------------------------------------------------------------------
1897inline float vtkMath::RadiansFromDegrees(float x)
1898{
1899 return x * 0.017453292f;
1900}
1901
1902//----------------------------------------------------------------------------
1903inline double vtkMath::RadiansFromDegrees(double x)
1904{
1905 return x * 0.017453292519943295;
1906}
1907
1908//----------------------------------------------------------------------------
1909inline float vtkMath::DegreesFromRadians(float x)
1910{
1911 return x * 57.2957795131f;
1912}
1913
1914//----------------------------------------------------------------------------
1915inline double vtkMath::DegreesFromRadians(double x)
1916{
1917 return x * 57.29577951308232;
1918}
1919
1920//----------------------------------------------------------------------------
1921inline bool vtkMath::IsPowerOfTwo(vtkTypeUInt64 x)
1922{
1923 return ((x != 0) & ((x & (x - 1)) == 0));
1924}
1925
1926//----------------------------------------------------------------------------
1927// Credit goes to Peter Hart and William Lewis on comp.lang.python 1997
1929{
1930 unsigned int z = static_cast<unsigned int>(((x > 0) ? x - 1 : 0));
1931 z |= z >> 1;
1932 z |= z >> 2;
1933 z |= z >> 4;
1934 z |= z >> 8;
1935 z |= z >> 16;
1936 return static_cast<int>(z + 1);
1937}
1938
1939//----------------------------------------------------------------------------
1940// Modify the trunc() operation provided by static_cast<int>() to get floor(),
1941// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1942inline int vtkMath::Floor(double x)
1943{
1944 int i = static_cast<int>(x);
1945 return i - (i > x);
1946}
1947
1948//----------------------------------------------------------------------------
1949// Modify the trunc() operation provided by static_cast<int>() to get ceil(),
1950// Note that in C++ conditions evaluate to values of 1 or 0 (true or false).
1951inline int vtkMath::Ceil(double x)
1952{
1953 int i = static_cast<int>(x);
1954 return i + (i < x);
1955}
1956
1957//----------------------------------------------------------------------------
1958template <class T>
1959inline T vtkMath::Min(const T& a, const T& b)
1960{
1961 return (b <= a ? b : a);
1962}
1963
1964//----------------------------------------------------------------------------
1965template <class T>
1966inline T vtkMath::Max(const T& a, const T& b)
1967{
1968 return (b > a ? b : a);
1969}
1970
1971//----------------------------------------------------------------------------
1972float vtkMath::Normalize(float v[3])
1973{
1974 float den = vtkMath::Norm(v);
1975 if (den != 0.0)
1976 {
1977 for (int i = 0; i < 3; ++i)
1978 {
1979 v[i] /= den;
1980 }
1981 }
1982 return den;
1983}
1984
1985//----------------------------------------------------------------------------
1986double vtkMath::Normalize(double v[3])
1987{
1988 double den = vtkMath::Norm(v);
1989 if (den != 0.0)
1990 {
1991 for (int i = 0; i < 3; ++i)
1992 {
1993 v[i] /= den;
1994 }
1995 }
1996 return den;
1997}
1998
1999//----------------------------------------------------------------------------
2000inline float vtkMath::Normalize2D(float v[2])
2001{
2002 float den = vtkMath::Norm2D(v);
2003 if (den != 0.0)
2004 {
2005 for (int i = 0; i < 2; ++i)
2006 {
2007 v[i] /= den;
2008 }
2009 }
2010 return den;
2011}
2012
2013//----------------------------------------------------------------------------
2014inline double vtkMath::Normalize2D(double v[2])
2015{
2016 double den = vtkMath::Norm2D(v);
2017 if (den != 0.0)
2018 {
2019 for (int i = 0; i < 2; ++i)
2020 {
2021 v[i] /= den;
2022 }
2023 }
2024 return den;
2025}
2026
2027//----------------------------------------------------------------------------
2028inline float vtkMath::Determinant3x3(const float c1[3], const float c2[3], const float c3[3])
2029{
2030 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
2031 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
2032}
2033
2034//----------------------------------------------------------------------------
2035inline double vtkMath::Determinant3x3(const double c1[3], const double c2[3], const double c3[3])
2036{
2037 return c1[0] * c2[1] * c3[2] + c2[0] * c3[1] * c1[2] + c3[0] * c1[1] * c2[2] -
2038 c1[0] * c3[1] * c2[2] - c2[0] * c1[1] * c3[2] - c3[0] * c2[1] * c1[2];
2039}
2040
2041//----------------------------------------------------------------------------
2043 double a1, double a2, double a3, double b1, double b2, double b3, double c1, double c2, double c3)
2044{
2045 return (a1 * vtkMath::Determinant2x2(b2, b3, c2, c3) -
2046 b1 * vtkMath::Determinant2x2(a2, a3, c2, c3) + c1 * vtkMath::Determinant2x2(a2, a3, b2, b3));
2047}
2048
2049//----------------------------------------------------------------------------
2050inline float vtkMath::Distance2BetweenPoints(const float p1[3], const float p2[3])
2051{
2052 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
2053 (p1[2] - p2[2]) * (p1[2] - p2[2]));
2054}
2055
2056//----------------------------------------------------------------------------
2057inline double vtkMath::Distance2BetweenPoints(const double p1[3], const double p2[3])
2058{
2059 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
2060 (p1[2] - p2[2]) * (p1[2] - p2[2]));
2061}
2062
2063//----------------------------------------------------------------------------
2064template <typename ReturnTypeT, typename TupleRangeT1, typename TupleRangeT2, typename EnableT>
2065inline ReturnTypeT vtkMath::Distance2BetweenPoints(const TupleRangeT1& p1, const TupleRangeT2& p2)
2066{
2067 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]) +
2068 (p1[2] - p2[2]) * (p1[2] - p2[2]));
2069}
2070
2071//------------------------------------------------------------------------------
2072inline double vtkMath::Distance2BetweenPoints2D(const double p1[2], const double p2[2])
2073{
2074 return ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]));
2075}
2076
2077//----------------------------------------------------------------------------
2078template <class VectorT1, class VectorT2, class VectorT3>
2079void vtkMath::Cross(VectorT1&& a, VectorT2&& b, VectorT3& c)
2080{
2082 ValueType Cx = a[1] * b[2] - a[2] * b[1];
2083 ValueType Cy = a[2] * b[0] - a[0] * b[2];
2084 ValueType Cz = a[0] * b[1] - a[1] * b[0];
2085 c[0] = Cx;
2086 c[1] = Cy;
2087 c[2] = Cz;
2088}
2089
2090//----------------------------------------------------------------------------
2091// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
2092inline void vtkMath::Cross(const float a[3], const float b[3], float c[3])
2093{
2094 float Cx = a[1] * b[2] - a[2] * b[1];
2095 float Cy = a[2] * b[0] - a[0] * b[2];
2096 float Cz = a[0] * b[1] - a[1] * b[0];
2097 c[0] = Cx;
2098 c[1] = Cy;
2099 c[2] = Cz;
2100}
2101
2102//----------------------------------------------------------------------------
2103// Cross product of two 3-vectors. Result (a x b) is stored in c[3].
2104inline void vtkMath::Cross(const double a[3], const double b[3], double c[3])
2105{
2106 double Cx = a[1] * b[2] - a[2] * b[1];
2107 double Cy = a[2] * b[0] - a[0] * b[2];
2108 double Cz = a[0] * b[1] - a[1] * b[0];
2109 c[0] = Cx;
2110 c[1] = Cy;
2111 c[2] = Cz;
2112}
2113
2114//----------------------------------------------------------------------------
2115template <class T>
2116inline double vtkDeterminant3x3(const T A[3][3])
2117{
2118 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] -
2119 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];
2120}
2121
2122//----------------------------------------------------------------------------
2123inline double vtkMath::Determinant3x3(const float A[3][3])
2124{
2125 return vtkDeterminant3x3(A);
2126}
2127
2128//----------------------------------------------------------------------------
2129inline double vtkMath::Determinant3x3(const double A[3][3])
2130{
2131 return vtkDeterminant3x3(A);
2132}
2133
2134//----------------------------------------------------------------------------
2135template <class T>
2136inline T vtkMath::ClampValue(const T& value, const T& min, const T& max)
2137{
2138 assert("pre: valid_range" && min <= max);
2139 return std::clamp(value, min, max);
2140}
2141
2142//----------------------------------------------------------------------------
2143inline void vtkMath::ClampValue(double* value, const double range[2])
2144{
2145 if (value && range)
2146 {
2147 assert("pre: valid_range" && range[0] <= range[1]);
2148
2149 *value = vtkMath::ClampValue(*value, range[0], range[1]);
2150 }
2151}
2152
2153//----------------------------------------------------------------------------
2154inline void vtkMath::ClampValue(double value, const double range[2], double* clamped_value)
2155{
2156 if (range && clamped_value)
2157 {
2158 assert("pre: valid_range" && range[0] <= range[1]);
2159
2160 *clamped_value = vtkMath::ClampValue(value, range[0], range[1]);
2161 }
2162}
2163
2164// ---------------------------------------------------------------------------
2165inline double vtkMath::ClampAndNormalizeValue(double value, const double range[2])
2166{
2167 assert("pre: valid_range" && range[0] <= range[1]);
2168
2169 double result;
2170 if (range[0] == range[1])
2171 {
2172 result = 0.0;
2173 }
2174 else
2175 {
2176 // clamp
2177 result = vtkMath::ClampValue(value, range[0], range[1]);
2178
2179 // normalize
2180 result = (result - range[0]) / (range[1] - range[0]);
2181 }
2182
2183 assert("post: valid_result" && result >= 0.0 && result <= 1.0);
2184
2185 return result;
2186}
2187
2188//-----------------------------------------------------------------------------
2189template <class T1, class T2>
2190inline void vtkMath::TensorFromSymmetricTensor(const T1 symmTensor[9], T2 tensor[9])
2191{
2192 for (int i = 0; i < 3; ++i)
2193 {
2194 tensor[4 * i] = symmTensor[i];
2195 }
2196 tensor[1] = tensor[3] = symmTensor[3];
2197 tensor[2] = tensor[6] = symmTensor[5];
2198 tensor[5] = tensor[7] = symmTensor[4];
2199}
2200
2201//-----------------------------------------------------------------------------
2202template <class T>
2204{
2205 tensor[6] = tensor[5]; // XZ
2206 tensor[7] = tensor[4]; // YZ
2207 tensor[8] = tensor[2]; // ZZ
2208 tensor[4] = tensor[1]; // YY
2209 tensor[5] = tensor[7]; // YZ
2210 tensor[2] = tensor[6]; // XZ
2211 tensor[1] = tensor[3]; // XY
2212}
2213VTK_ABI_NAMESPACE_END
2214
2215namespace
2216{
2217template <class QuaternionT, class MatrixT>
2218inline void vtkQuaternionToMatrix3x3(QuaternionT&& quat, MatrixT&& A)
2219{
2221
2222 Scalar ww = quat[0] * quat[0];
2223 Scalar wx = quat[0] * quat[1];
2224 Scalar wy = quat[0] * quat[2];
2225 Scalar wz = quat[0] * quat[3];
2226
2227 Scalar xx = quat[1] * quat[1];
2228 Scalar yy = quat[2] * quat[2];
2229 Scalar zz = quat[3] * quat[3];
2230
2231 Scalar xy = quat[1] * quat[2];
2232 Scalar xz = quat[1] * quat[3];
2233 Scalar yz = quat[2] * quat[3];
2234
2235 Scalar rr = xx + yy + zz;
2236 // normalization factor, just in case quaternion was not normalized
2237 Scalar f = 1 / (ww + rr);
2238 Scalar s = (ww - rr) * f;
2239 f *= 2;
2240
2242
2243 MatrixT& Ar = A;
2244 Wrapper::template Get<0, 0>(Ar) = xx * f + s;
2245 Wrapper::template Get<1, 0>(Ar) = (xy + wz) * f;
2246 Wrapper::template Get<2, 0>(Ar) = (xz - wy) * f;
2247
2248 Wrapper::template Get<0, 1>(Ar) = (xy - wz) * f;
2249 Wrapper::template Get<1, 1>(Ar) = yy * f + s;
2250 Wrapper::template Get<2, 1>(Ar) = (yz + wx) * f;
2251
2252 Wrapper::template Get<0, 2>(Ar) = (xz + wy) * f;
2253 Wrapper::template Get<1, 2>(Ar) = (yz - wx) * f;
2254 Wrapper::template Get<2, 2>(Ar) = zz * f + s;
2255}
2256} // anonymous namespace
2257
2258VTK_ABI_NAMESPACE_BEGIN
2259//------------------------------------------------------------------------------
2260inline void vtkMath::QuaternionToMatrix3x3(const float quat[4], float A[3][3])
2261{
2262 vtkQuaternionToMatrix3x3(quat, A);
2263}
2264
2265//------------------------------------------------------------------------------
2266inline void vtkMath::QuaternionToMatrix3x3(const double quat[4], double A[3][3])
2267{
2268 vtkQuaternionToMatrix3x3(quat, A);
2269}
2270
2271//-----------------------------------------------------------------------------
2272template <class QuaternionT, class MatrixT, class EnableT>
2273inline void vtkMath::QuaternionToMatrix3x3(QuaternionT&& q, MatrixT&& A)
2274{
2275 vtkQuaternionToMatrix3x3(std::forward<QuaternionT>(q), std::forward<MatrixT>(A));
2276}
2277VTK_ABI_NAMESPACE_END
2278
2279namespace
2280{
2281//------------------------------------------------------------------------------
2282// The solution is based on
2283// Berthold K. P. Horn (1987),
2284// "Closed-form solution of absolute orientation using unit quaternions,"
2285// Journal of the Optical Society of America A, 4:629-642
2286template <class MatrixT, class QuaternionT>
2287inline void vtkMatrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& quat)
2288{
2290
2291 Scalar N[4][4];
2292
2294
2295 MatrixT& Ar = A;
2296
2297 // on-diagonal elements
2298 N[0][0] = Wrapper::template Get<0, 0>(Ar) + Wrapper::template Get<1, 1>(Ar) +
2299 Wrapper::template Get<2, 2>(Ar);
2300 N[1][1] = Wrapper::template Get<0, 0>(Ar) - Wrapper::template Get<1, 1>(Ar) -
2301 Wrapper::template Get<2, 2>(Ar);
2302 N[2][2] = -Wrapper::template Get<0, 0>(Ar) + Wrapper::template Get<1, 1>(Ar) -
2303 Wrapper::template Get<2, 2>(Ar);
2304 N[3][3] = -Wrapper::template Get<0, 0>(Ar) - Wrapper::template Get<1, 1>(Ar) +
2305 Wrapper::template Get<2, 2>(Ar);
2306
2307 // off-diagonal elements
2308 N[0][1] = N[1][0] = Wrapper::template Get<2, 1>(Ar) - Wrapper::template Get<1, 2>(Ar);
2309 N[0][2] = N[2][0] = Wrapper::template Get<0, 2>(Ar) - Wrapper::template Get<2, 0>(Ar);
2310 N[0][3] = N[3][0] = Wrapper::template Get<1, 0>(Ar) - Wrapper::template Get<0, 1>(Ar);
2311
2312 N[1][2] = N[2][1] = Wrapper::template Get<1, 0>(Ar) + Wrapper::template Get<0, 1>(Ar);
2313 N[1][3] = N[3][1] = Wrapper::template Get<0, 2>(Ar) + Wrapper::template Get<2, 0>(Ar);
2314 N[2][3] = N[3][2] = Wrapper::template Get<2, 1>(Ar) + Wrapper::template Get<1, 2>(Ar);
2315
2316 Scalar eigenvectors[4][4], eigenvalues[4];
2317
2318 // convert into format that JacobiN can use,
2319 // then use Jacobi to find eigenvalues and eigenvectors
2320 Scalar *NTemp[4], *eigenvectorsTemp[4];
2321 for (int i = 0; i < 4; ++i)
2322 {
2323 NTemp[i] = N[i];
2324 eigenvectorsTemp[i] = eigenvectors[i];
2325 }
2326 vtkMath::JacobiN(NTemp, 4, eigenvalues, eigenvectorsTemp);
2327
2328 // the first eigenvector is the one we want
2329 quat[0] = eigenvectors[0][0];
2330 quat[1] = eigenvectors[1][0];
2331 quat[2] = eigenvectors[2][0];
2332 quat[3] = eigenvectors[3][0];
2333}
2334} // anonymous namespace
2335
2336VTK_ABI_NAMESPACE_BEGIN
2337//------------------------------------------------------------------------------
2338inline void vtkMath::Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
2339{
2340 vtkMatrix3x3ToQuaternion(A, quat);
2341}
2342
2343//------------------------------------------------------------------------------
2344inline void vtkMath::Matrix3x3ToQuaternion(const double A[3][3], double quat[4])
2345{
2346 vtkMatrix3x3ToQuaternion(A, quat);
2347}
2348
2349//-----------------------------------------------------------------------------
2350template <class MatrixT, class QuaternionT, class EnableT>
2351inline void vtkMath::Matrix3x3ToQuaternion(MatrixT&& A, QuaternionT&& q)
2352{
2353 vtkMatrix3x3ToQuaternion(std::forward<MatrixT>(A), std::forward<QuaternionT>(q));
2354}
2355VTK_ABI_NAMESPACE_END
2356
2357namespace vtk_detail
2358{
2359VTK_ABI_NAMESPACE_BEGIN
2360// Can't specialize templates inside a template class, so we move the impl here.
2361template <typename OutT>
2362void RoundDoubleToIntegralIfNecessary(double val, OutT* ret)
2363{ // OutT is integral -- clamp and round
2364 if (!vtkMath::IsNan(val))
2365 {
2366 double min = static_cast<double>(vtkTypeTraits<OutT>::Min());
2367 double max = static_cast<double>(vtkTypeTraits<OutT>::Max());
2368 val = vtkMath::ClampValue(val, min, max);
2369 *ret = static_cast<OutT>((val >= 0.0) ? (val + 0.5) : (val - 0.5));
2370 }
2371 else
2372 *ret = 0;
2373}
2374template <>
2375inline void RoundDoubleToIntegralIfNecessary(double val, double* retVal)
2376{ // OutT is double: passthrough
2377 *retVal = val;
2378}
2379template <>
2380inline void RoundDoubleToIntegralIfNecessary(double val, float* retVal)
2381{ // OutT is float -- just clamp (as doubles, then the cast to float is well-defined.)
2382 if (!vtkMath::IsNan(val))
2383 {
2384 double min = static_cast<double>(vtkTypeTraits<float>::Min());
2385 double max = static_cast<double>(vtkTypeTraits<float>::Max());
2386 val = vtkMath::ClampValue(val, min, max);
2387 }
2388
2389 *retVal = static_cast<float>(val);
2390}
2391VTK_ABI_NAMESPACE_END
2392} // end namespace vtk_detail
2393
2394VTK_ABI_NAMESPACE_BEGIN
2395//-----------------------------------------------------------------------------
2396#if defined(VTK_HAS_ISINF) || defined(VTK_HAS_STD_ISINF)
2397#define VTK_MATH_ISINF_IS_INLINE
2398inline vtkTypeBool vtkMath::IsInf(double x)
2399{
2400#if defined(VTK_HAS_STD_ISINF)
2401 return std::isinf(x);
2402#else
2403 return (isinf(x) != 0); // Force conversion to bool
2404#endif
2405}
2406#endif
2407
2408//-----------------------------------------------------------------------------
2409#if defined(VTK_HAS_ISNAN) || defined(VTK_HAS_STD_ISNAN)
2410#define VTK_MATH_ISNAN_IS_INLINE
2412{
2413#if defined(VTK_HAS_STD_ISNAN)
2414 return std::isnan(x);
2415#else
2416 return (isnan(x) != 0); // Force conversion to bool
2417#endif
2418}
2419#endif
2420
2421//-----------------------------------------------------------------------------
2422#if defined(VTK_HAS_ISFINITE) || defined(VTK_HAS_STD_ISFINITE) || defined(VTK_HAS_FINITE)
2423#define VTK_MATH_ISFINITE_IS_INLINE
2424inline bool vtkMath::IsFinite(double x)
2425{
2426#if defined(VTK_HAS_STD_ISFINITE)
2427 return std::isfinite(x);
2428#elif defined(VTK_HAS_ISFINITE)
2429 return (isfinite(x) != 0); // Force conversion to bool
2430#else
2431 return (finite(x) != 0); // Force conversion to bool
2432#endif
2433}
2434#endif
2435
2436VTK_ABI_NAMESPACE_END
2437#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:108
static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1 &p1, const TupleRangeT2 &p2)
Compute distance squared between two points p1 and p2.
Definition vtkMath.h:2065
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:1796
static void XYZToProLab(const double xyz[3], double prolab[3])
Convert Color from the CIE XYZ system to ProLAB.
Definition vtkMath.h:1513
static double Norm2D(const double x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:871
static void XYZToProLab(double x, double y, double z, double *L, double *a, double *b)
Convert Color from the CIE XYZ system to ProLAB.
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:1876
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:836
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:1889
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:1030
static void XYZToLab(const double xyz[3], double lab[3])
Convert Color from the CIE XYZ system to CIE-L*ab.
Definition vtkMath.h:1535
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:1117
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:1524
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:1484
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:898
static T Max(const T &a, const T &b)
Returns the maximum of the two arguments provided.
Definition vtkMath.h:1966
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:850
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:1951
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:1489
~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:1013
static void RGBToXYZ(const double rgb[3], double xyz[3])
Convert color from the RGB system to CIE XYZ.
Definition vtkMath.h:1557
static void QuaternionToMatrix3x3(const float quat[4], float A[3][3])
Convert a quaternion to a 3x3 rotation matrix.
Definition vtkMath.h:2260
static int NearestPowerOfTwo(int x)
Compute the nearest power of two that is not less than x.
Definition vtkMath.h:1928
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:1140
static float Normalize2D(float v[2])
Normalize (in place) a 2-vector.
Definition vtkMath.h:2000
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:1897
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:1831
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 void RGBToProLab(double red, double green, double blue, double *L, double *a, double *b)
Convert color from the RGB system to Prolab The input RGB must be values in the range [0,...
static void ProLabToXYZ(const double prolab[3], double xyz[3])
Convert color from the ProLAB system to CIE XYZ.
Definition vtkMath.h:1501
static vtkTypeBool AreBoundsInitialized(const double bounds[6])
Are the bounds initialized?
Definition vtkMath.h:1635
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 ProLabToXYZ(double L, double a, double b, double *x, double *y, double *z)
Convert color from the ProLAB system to CIE XYZ.
static void RGBToProLab(const double rgb[3], double prolab[3])
Convert color from the RGB system to Prolab The input RGB must be values in the range [0,...
Definition vtkMath.h:1598
static void ProLabToRGB(double L, double a, double b, double *red, double *green, double *blue)
Convert color from the ProLab system to RGB.
static void Matrix3x3ToQuaternion(const float A[3][3], float quat[4])
Convert a 3x3 matrix into a quaternion.
Definition vtkMath.h:2338
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:1546
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:2165
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:831
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:1572
static void ProLabToRGB(const double prolab[3], double rgb[3])
Convert color from the ProLab system to RGB.
Definition vtkMath.h:1583
static int Floor(double x)
Rounds a double to the nearest integer not greater than itself.
Definition vtkMath.h:1942
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:897
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:1469
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:1909
static float Determinant2x2(const float c1[2], const float c2[2])
Compute determinant of 2x2 matrix.
Definition vtkMath.h:888
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:1620
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:1921
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:1972
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:2123
static float Dot2D(const float x[2], const float y[2])
Dot product of two 2-vectors.
Definition vtkMath.h:826
ConvolutionMode
Support the convolution operations.
Definition vtkMath.h:1802
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:1464
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:1073
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:1609
static float Norm2D(const float x[2])
Compute the norm of a 2-vector.
Definition vtkMath.h:865
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:1959
static void InvertMatrix(MatrixT1 &&M1, MatrixT2 &&M2)
Computes the inverse of input matrix M1 into M2.
Definition vtkMath.h:1096
static void Cross(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Cross product of two 3-vectors.
Definition vtkMath.h:2079
static void MultiplyMatrix(MatrixT1 &&M1, MatrixT2 &&M2, MatrixT3 &&M3)
Multiply matrices such that M3 = M1 x M2.
Definition vtkMath.h:972
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:2136
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 double Distance2BetweenPoints2D(const double p1[2], const double p2[2])
Compute distance squared between two 2D points p1 and p2.
Definition vtkMath.h:2072
static vtkMatrixUtilities::ScalarTypeExtractor< VectorT >::value_type SquaredNorm(VectorT &&x)
Computes the dot product between 2 vectors x and y.
Definition vtkMath.h:1049
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:1000
Park and Miller Sequence of pseudo random numbers.
represent and manipulate 3D points
Definition vtkPoints.h:140
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:2362
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:2116
#define Distance2BetweenPoints2D(p1, p2)
int vtkIdType
Definition vtkType.h:363
#define max(a, b)