23#include <fmt/format.h>
39template<std::
floating_po
int T>
40constexpr T
round(
const T& value,
size_t digits)
42 auto factor = std::pow(10, digits);
43 return std::round(value * factor) / factor;
50template<std::
floating_po
int T>
53 if (value == 0.0) {
return 0.0; }
55 auto absVal = gcem::abs(value);
56 auto log10 =
static_cast<int32_t
>(gcem::log10(absVal));
57 auto exp = log10 + (log10 > 0 || (log10 == 0 && absVal >= 1.0));
58 auto fac =
static_cast<T
>(digits) -
static_cast<T
>(exp);
60 auto factor =
static_cast<T
>(gcem::pow(10.0, fac));
63 return static_cast<T
>(gcem::round(value * factor) / factor);
72template<std::
integral Out,
size_t Bits, std::
integral In>
75 static_assert(Bits <
sizeof(In) * 8);
76 static_assert(Bits <
sizeof(Out) * 8);
78 constexpr size_t N =
sizeof(Out) * 8 - Bits;
79 return static_cast<Out
>(
static_cast<Out
>((in &
static_cast<In
>(gcem::pow(2, Bits) - 1)) << N) >> N);
88template<
typename Derived>
94 Eigen::Matrix<typename Derived::Scalar, 3, 3> skewMat;
95 skewMat << 0, -a(2), a(1),
106template<
typename Derived>
112 Eigen::Matrix<typename Derived::Scalar, 3, 3> skewMat2;
113 skewMat2 << std::pow(a(2), 2) + std::pow(a(1), 2), a(0) * a(1), a(0) * a(2),
114 a(0) * a(1), std::pow(a(2), 2) + std::pow(a(0), 2), a(1) * a(2),
115 a(0) * a(2), a(1) * a(2), std::pow(a(0), 2) + std::pow(a(1), 2);
121template<std::
floating_po
int T>
124 return 1.0 / std::cos(x);
128template<std::
floating_po
int T>
131 return 1.0 / std::sin(x);
140 return (T(0) < val) - (val < T(0));
149template<
typename Derived>
150std::optional<std::pair<Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>,
151 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime>>>
156 auto n = Qmatrix.rows();
157 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> Q = Qmatrix.template triangularView<Eigen::Lower>();
158 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> L;
159 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime> D;
161 if constexpr (Derived::RowsAtCompileTime == Eigen::Dynamic)
163 L = Eigen::Matrix<typename Derived::Scalar, Eigen::Dynamic, Eigen::Dynamic>::Zero(n, n);
168 L = Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>::Zero();
172 for (Eigen::Index i = n - 1; i >= 0; i--)
175 if (Q(i, i) <= 0.0) {
return {}; }
176 L(i, seq(0, i)) = Q(i, seq(0, i)) / std::sqrt(Q(i, i));
178 for (Eigen::Index j = 0; j <= i - 1; j++)
180 Q(j, seq(0, j)) -= L(i, seq(0, j)) * L(i, j);
182 L(i, seq(0, i)) /= L(i, i);
185 return std::make_pair(L, D);
193template<
typename Derived>
194std::optional<std::pair<Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>,
195 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime>>>
201 typename Derived::PlainObject L = Q.template triangularView<Eigen::Lower>();
202 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime> D;
205 if constexpr (Derived::RowsAtCompileTime == Eigen::Dynamic) { D.setZero(n); }
206 else { D.setZero(); }
208 for (Eigen::Index j = n - 1; j >= 0; j--)
210 for (Eigen::Index i = n - 1; i >= j + 1; i--)
212 L(i, j) = (L(i, j) - L(seq(i + 1, n - 1), j).dot(L(seq(i + 1, n - 1), i))) / L(i, i);
214 double t = L(j, j) - L(seq(j + 1, n - 1), j).dot(L(seq(j + 1, n - 1), j));
215 if (t <= 0.0) {
return {}; }
216 double c = t / L(j, j);
217 cmin = std::min(c, cmin);
218 L(j, j) = std::sqrt(t);
220 for (Eigen::Index i = 0; i < n; i++)
222 L.row(i).leftCols(i) /= L(i, i);
223 D(i) = std::pow(L(i, i), 2.0);
227 return std::make_pair(L, D);
238template<
typename DerivedA,
typename DerivedQ>
241 static_assert(DerivedA::ColsAtCompileTime == Eigen::Dynamic || DerivedA::ColsAtCompileTime == 1);
246 return a.transpose() * Q.inverse() * a;
289 return -1.0 * fabs(x);
297template<
typename Derived>
298typename Derived::PlainObject
lerp(
const Eigen::MatrixBase<Derived>& a,
const Eigen::MatrixBase<Derived>& b,
const typename Derived::Scalar& t)
300 return a + t * (b - a);
316 auto i =
static_cast<size_t>(std::distance(data.begin(), std::upper_bound(data.begin(), data.end(), value)));
318 if (i == data.size() - 1) { i--; }
319 const auto& lb = data.at(i);
320 const auto& ub = data.at(i + 1);
321 double t = (value - lb) / (ub - lb);
323 return { .l = i, .u = i + 1, .t = t };
342 const auto& c00,
const auto& c10,
343 const auto& c01,
const auto& c11)
345 auto a = c00 * (1 - tx) + c10 * tx;
346 auto b = c01 * (1 - tx) + c11 * tx;
347 return a * (1 - ty) + b * ty;
#define INS_ASSERT_USER_ERROR(_EXP, _MSG)
Assert function with message.
Definition Assert.h:21
T sec(const T &x)
Calculates the secant of a value (sec(x) = csc(pi/2 - x) = 1 / cos(x))
Definition Math.hpp:122
LerpSearchResult lerpSearch(const auto &data, const auto &value)
Searches the value in the data container.
Definition Math.hpp:314
T csc(const T &x)
Calculates the cosecant of a value (csc(x) = sec(pi/2 - x) = 1 / sin(x))
Definition Math.hpp:129
DerivedA::Scalar squaredNormVectorMatrix(const Eigen::MatrixBase< DerivedA > &a, const Eigen::MatrixBase< DerivedQ > &Q)
Calculates the squared norm of the vector and matrix.
Definition Math.hpp:239
int sgn(const T &val)
Returns the sign of the given value.
Definition Math.hpp:138
std::optional< std::pair< Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime >, Eigen::Vector< typename Derived::Scalar, Derived::RowsAtCompileTime > > > LtDLdecomp_outerProduct(const Eigen::MatrixBase< Derived > &Qmatrix)
Find (L^T D L)-decomposition of Q-matrix via outer product method.
Definition Math.hpp:152
T sign(const T &x, const T &y)
Change the sign of x according to the value of y.
Definition Math.hpp:282
Eigen::Matrix< typename Derived::Scalar, 3, 3 > skewSymmetricMatrixSquared(const Eigen::MatrixBase< Derived > &a)
Calculates the square of a skew symmetric matrix of the given vector.
Definition Math.hpp:107
auto bilinearInterpolation(const double &tx, const double &ty, const auto &c00, const auto &c10, const auto &c01, const auto &c11)
Bilinear interpolation.
Definition Math.hpp:341
constexpr T roundSignificantDigits(T value, size_t digits)
Round the number to the specified amount of significant digits.
Definition Math.hpp:51
double calcEllipticalIntegral(double phi, double m)
Calculates the incomplete elliptical integral of the second kind.
Eigen::Matrix< typename Derived::Scalar, 3, 3 > skewSymmetricMatrix(const Eigen::MatrixBase< Derived > &a)
Calculates the skew symmetric matrix of the given vector. This is needed to perform the cross product...
Definition Math.hpp:89
uint64_t factorial(uint64_t n)
Calculates the factorial of an unsigned integer.
double normalCDF(double value)
Calculates the cumulative distribution function (CDF) of the standard normal distribution.
constexpr T round(const T &value, size_t digits)
Round the number to the specified amount of digits.
Definition Math.hpp:40
Derived::PlainObject lerp(const Eigen::MatrixBase< Derived > &a, const Eigen::MatrixBase< Derived > &b, const typename Derived::Scalar &t)
Linear interpolation between vectors.
Definition Math.hpp:298
constexpr Out interpretAs(In in)
Interprets the input integer with certain amount of Bits as Output type. Takes care of sign extension...
Definition Math.hpp:73
std::optional< std::pair< Eigen::Matrix< typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime >, Eigen::Vector< typename Derived::Scalar, Derived::RowsAtCompileTime > > > LtDLdecomp_choleskyFact(const Eigen::MatrixBase< Derived > &Q)
Find (L^T D L)-decomposition of Q-matrix via a backward Cholesky factorization in a bordering method ...
Definition Math.hpp:196
Lerp Search Result.
Definition Math.hpp:305
size_t u
Upper bound index (l + 1)
Definition Math.hpp:307
double t
[0, 1] for Interpolation, otherwise Extrapolation
Definition Math.hpp:308
size_t l
Lower bound index.
Definition Math.hpp:306