21#include <fmt/format.h>
38 typename = std::enable_if_t<std::is_floating_point_v<T>>>
39constexpr T
round(
const T& value,
size_t digits)
41 auto factor = std::pow(10, digits);
42 return std::round(value * factor) / factor;
50 typename = std::enable_if_t<std::is_floating_point_v<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<
typename Out,
size_t Bits,
typename In,
73 typename = std::enable_if_t<std::is_integral_v<Out>>,
74 typename = std::enable_if_t<std::is_integral_v<In>>>
77 static_assert(Bits <
sizeof(In) * 8);
78 static_assert(Bits <
sizeof(Out) * 8);
80 constexpr size_t N =
sizeof(Out) * 8 - Bits;
81 return static_cast<Out
>(
static_cast<Out
>((in &
static_cast<In
>(gcem::pow(2, Bits) - 1)) << N) >> N);
90template<
typename Derived>
96 Eigen::Matrix<typename Derived::Scalar, 3, 3> skewMat;
97 skewMat << 0, -a(2), a(1),
108template<
typename Derived>
114 Eigen::Matrix<typename Derived::Scalar, 3, 3> skewMat2;
115 skewMat2 << std::pow(a(2), 2) + std::pow(a(1), 2), a(0) * a(1), a(0) * a(2),
116 a(0) * a(1), std::pow(a(2), 2) + std::pow(a(0), 2), a(1) * a(2),
117 a(0) * a(2), a(1) * a(2), std::pow(a(0), 2) + std::pow(a(1), 2);
124 typename = std::enable_if_t<std::is_floating_point_v<T>>>
127 return 1.0 / std::cos(x);
132 typename = std::enable_if_t<std::is_floating_point_v<T>>>
135 return 1.0 / std::sin(x);
144 return (T(0) < val) - (val < T(0));
153template<
typename Derived>
154std::pair<Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>,
155 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime>>
160 auto n = Qmatrix.rows();
161 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> Q = Qmatrix.template triangularView<Eigen::Lower>();
162 Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime> L;
163 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime> D;
165 if constexpr (Derived::RowsAtCompileTime == Eigen::Dynamic)
167 L = Eigen::Matrix<typename Derived::Scalar, Eigen::Dynamic, Eigen::Dynamic>::Zero(n, n);
172 L = Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>::Zero();
176 for (Eigen::Index i = n - 1; i >= 0; i--)
179 L(i, seq(0, i)) = Q(i, seq(0, i)) / std::sqrt(Q(i, i));
181 for (Eigen::Index j = 0; j <= i - 1; j++)
183 Q(j, seq(0, j)) -= L(i, seq(0, j)) * L(i, j);
185 L(i, seq(0, i)) /= L(i, i);
196template<
typename Derived>
197std::pair<Eigen::Matrix<typename Derived::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime>,
198 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime>>
204 typename Derived::PlainObject L = Q.template triangularView<Eigen::Lower>();
205 Eigen::Vector<typename Derived::Scalar, Derived::RowsAtCompileTime> D;
208 if constexpr (Derived::RowsAtCompileTime == Eigen::Dynamic) { D.setZero(n); }
209 else { D.setZero(); }
211 for (Eigen::Index j = n - 1; j >= 0; j--)
213 for (Eigen::Index i = n - 1; i >= j + 1; i--)
215 L(i, j) = (L(i, j) - L(seq(i + 1, n - 1), j).dot(L(seq(i + 1, n - 1), i))) / L(i, i);
217 double t = L(j, j) - L(seq(j + 1, n - 1), j).dot(L(seq(j + 1, n - 1), j));
218 double c = t / L(j, j);
223 L(j, j) = std::sqrt(t);
225 for (Eigen::Index i = 0; i < n; i++)
227 L.row(i).leftCols(i) /= L(i, i);
228 D(i) = std::pow(L(i, i), 2.0);
243template<
typename DerivedA,
typename DerivedQ>
246 static_assert(DerivedA::ColsAtCompileTime == Eigen::Dynamic || DerivedA::ColsAtCompileTime == 1);
251 return a.transpose() * Q.inverse() * a;
294 return -1.0 * fabs(x);
302template<
typename Derived>
303typename Derived::PlainObject
lerp(
const Eigen::MatrixBase<Derived>& a,
const Eigen::MatrixBase<Derived>& b,
const typename Derived::Scalar& t)
305 return a + t * (b - a);
#define INS_ASSERT_USER_ERROR(_EXP, _MSG)
Assert function with message.
Definition Assert.h:21
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:199
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:244
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:75
T csc(const T &x)
Calculates the cosecant of a value (csc(x) = sec(pi/2 - x) = 1 / sin(x))
Definition Math.hpp:133
int sgn(const T &val)
Returns the sign of the given value.
Definition Math.hpp:142
T sign(const T &x, const T &y)
Change the sign of x according to the value of y.
Definition Math.hpp:287
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:109
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:91
constexpr T roundSignificantDigits(T value, size_t digits)
Round the number to the specified amount of significant digits.
Definition Math.hpp:51
T sec(const T &x)
Calculates the secant of a value (sec(x) = csc(pi/2 - x) = 1 / cos(x))
Definition Math.hpp:125
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.
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:303
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:156
constexpr T round(const T &value, size_t digits)
Round the number to the specified amount of digits.
Definition Math.hpp:39