0.3.0
Loading...
Searching...
No Matches
Eigen.hpp
Go to the documentation of this file.
1// This file is part of INSTINCT, the INS Toolkit for Integrated
2// Navigation Concepts and Training by the Institute of Navigation of
3// the University of Stuttgart, Germany.
4//
5// This Source Code Form is subject to the terms of the Mozilla Public
6// License, v. 2.0. If a copy of the MPL was not distributed with this
7// file, You can obtain one at https://mozilla.org/MPL/2.0/.
8
13
14#pragma once
15
16#include <util/Logger.hpp>
17#include <Eigen/Core>
18#include <Eigen/Dense>
19
20#include <nlohmann/json.hpp>
21using json = nlohmann::json;
22
23#include <fmt/ostream.h>
24
25#include "Assert.h"
26
27namespace Eigen
28{
29using Array5d = Array<double, 5, 1>;
30using Array6d = Array<double, 6, 1>;
31using Vector5d = Matrix<double, 5, 1>;
32using Vector6d = Matrix<double, 6, 1>;
33using Matrix5d = Matrix<double, 5, 5>;
34using Matrix6d = Matrix<double, 6, 6>;
35
36using Array3ld = Array<long double, 3, 1>;
37
38using Vector3ld = Matrix<long double, 3, 1>;
39using Vector4ld = Matrix<long double, 4, 1>;
40
41using Matrix3ld = Matrix<long double, 3, 3>;
42using Matrix4ld = Matrix<long double, 4, 4>;
43
44using Quaternionld = Quaternion<long double>;
45
46using AngleAxisld = AngleAxis<long double>;
47
54template<typename _Scalar, int _Rows, int _Cols>
55void to_json(json& j, const Matrix<_Scalar, _Rows, _Cols>& matrix)
56{
57 for (int r = 0; r < matrix.rows(); r++)
58 {
59 for (int c = 0; c < matrix.cols(); c++)
60 {
61 if (std::isnan(matrix(r, c))) { j[std::to_string(r)][std::to_string(c)] = "NaN"; }
62 else { j[std::to_string(r)][std::to_string(c)] = matrix(r, c); }
63 }
64 }
65}
66
73template<typename _Scalar, int _Rows, int _Cols>
74void from_json(const json& j, Matrix<_Scalar, _Rows, _Cols>& matrix)
75{
76 if constexpr (_Rows == -1 || _Cols == -1)
77 {
78 int rows = -1;
79 int cols = -1;
80 for (int r = 0; j.contains(std::to_string(r)); r++)
81 {
82 rows = std::max(r, rows);
83 for (int c = 0; j.at(std::to_string(r)).contains(std::to_string(c)); c++)
84 {
85 cols = std::max(c, cols);
86 }
87 }
88 matrix = Eigen::Matrix<_Scalar, _Rows, _Cols>::Zero(rows + 1, cols + 1);
89 }
90
91 for (int r = 0; j.contains(std::to_string(r)); r++) // NOLINT(readability-misleading-indentation)
92 {
93 for (int c = 0; j.at(std::to_string(r)).contains(std::to_string(c)); c++)
94 {
95 if (j.at(std::to_string(r)).at(std::to_string(c)).is_string())
96 {
97 auto str = j.at(std::to_string(r)).at(std::to_string(c)).get<std::string>();
98 if (str == "NaN") { matrix(r, c) = std::nan(""); }
99 else { LOG_WARN("Reading matrix value failed at position ({}, {}). Value is an unknown string '{}'.", r, c, str); }
100 }
101 else if (j.at(std::to_string(r)).at(std::to_string(c)).is_number())
102 {
103 j.at(std::to_string(r)).at(std::to_string(c)).get_to(matrix(r, c));
104 }
105 else
106 {
107 LOG_WARN("Reading matrix value failed at position ({}, {}). Value has the type '{}' which cannot be converted into a floating point number.", r, c, j.at(std::to_string(r)).at(std::to_string(c)).type_name());
108 }
109 }
110 }
111}
112
113} // namespace Eigen
114
115namespace NAV
116{
117
122template<typename Derived>
123void removeRows(Eigen::DenseBase<Derived>& matrix, size_t index, size_t length)
124{
125 INS_ASSERT_USER_ERROR(static_cast<size_t>(matrix.rows()) >= index + length, "Tried to remove rows which do not exist");
126
127 std::vector<int> indicesToKeep;
128 indicesToKeep.reserve(static_cast<size_t>(matrix.rows()) - length);
129 for (int i = 0; i < static_cast<int>(index); i++) { indicesToKeep.push_back(i); }
130 for (int i = static_cast<int>(index + length); i < matrix.rows(); i++) { indicesToKeep.push_back(i); }
131
132 matrix = matrix(indicesToKeep, Eigen::all).eval();
133}
134
138template<typename Derived>
139void removeRows(Eigen::DenseBase<Derived>& matrix, const std::vector<int>& rowIndices)
140{
141 std::vector<int> rowIndicesToKeep;
142 rowIndicesToKeep.reserve(static_cast<size_t>(matrix.rows()) - rowIndices.size());
143 for (int i = 0; i < matrix.rows(); i++)
144 {
145 if (std::ranges::find(rowIndices, i) == rowIndices.end())
146 {
147 rowIndicesToKeep.push_back(i);
148 }
149 }
150
151 matrix = matrix(rowIndicesToKeep, Eigen::all).eval();
152}
153
158template<typename Derived>
159void removeCols(Eigen::DenseBase<Derived>& matrix, size_t index, size_t length)
160{
161 INS_ASSERT_USER_ERROR(static_cast<size_t>(matrix.cols()) >= index + length, "Tried to remove cols which do not exist");
162
163 std::vector<int> indicesToKeep;
164 indicesToKeep.reserve(static_cast<size_t>(matrix.cols()) - length);
165 for (int i = 0; i < static_cast<int>(index); i++) { indicesToKeep.push_back(i); }
166 for (int i = static_cast<int>(index + length); i < matrix.cols(); i++) { indicesToKeep.push_back(i); }
167
168 matrix = matrix(Eigen::all, indicesToKeep).eval();
169}
170
174template<typename Derived>
175void removeCols(Eigen::DenseBase<Derived>& matrix, const std::vector<int>& colIndices)
176{
177 std::vector<int> colIndicesToKeep;
178 colIndicesToKeep.reserve(static_cast<size_t>(matrix.cols()) - colIndices.size());
179 for (int i = 0; i < matrix.cols(); i++)
180 {
181 if (std::ranges::find(colIndices, i) == colIndices.end())
182 {
183 colIndicesToKeep.push_back(i);
184 }
185 }
186
187 matrix = matrix(Eigen::all, colIndicesToKeep).eval();
188}
189
196template<typename Derived>
197void removeRowsAndCols(Eigen::DenseBase<Derived>& matrix, size_t row, size_t rows, size_t col, size_t cols)
198{
199 INS_ASSERT_USER_ERROR(static_cast<size_t>(matrix.rows()) >= row + rows, "Tried to remove rows which do not exist");
200 INS_ASSERT_USER_ERROR(static_cast<size_t>(matrix.cols()) >= col + cols, "Tried to remove cols which do not exist");
201
202 std::vector<int> rowsToKeep;
203 rowsToKeep.reserve(static_cast<size_t>(matrix.rows()) - rows);
204 for (int i = 0; i < static_cast<int>(row); i++) { rowsToKeep.push_back(i); }
205 for (int i = static_cast<int>(row + rows); i < matrix.rows(); i++) { rowsToKeep.push_back(i); }
206
207 std::vector<int> colsToKeep;
208 colsToKeep.reserve(static_cast<size_t>(matrix.cols()) - cols);
209 for (int i = 0; i < static_cast<int>(col); i++) { colsToKeep.push_back(i); }
210 for (int i = static_cast<int>(col + cols); i < matrix.cols(); i++) { colsToKeep.push_back(i); }
211
212 matrix = matrix(rowsToKeep, colsToKeep).eval();
213}
214
219template<typename Derived>
220void removeRowsAndCols(Eigen::DenseBase<Derived>& matrix, const std::vector<int>& rowIndices, const std::vector<int>& colIndices)
221{
222 std::vector<int> rowIndicesToKeep;
223 rowIndicesToKeep.reserve(static_cast<size_t>(matrix.rows()) - rowIndices.size());
224 for (int i = 0; i < matrix.rows(); i++)
225 {
226 if (std::ranges::find(rowIndices, i) == rowIndices.end())
227 {
228 rowIndicesToKeep.push_back(i);
229 }
230 }
231
232 std::vector<int> colIndicesToKeep;
233 colIndicesToKeep.reserve(static_cast<size_t>(matrix.cols()) - colIndices.size());
234 for (int i = 0; i < matrix.cols(); i++)
235 {
236 if (std::ranges::find(colIndices, i) == colIndices.end())
237 {
238 colIndicesToKeep.push_back(i);
239 }
240 }
241
242 matrix = matrix(rowIndicesToKeep, colIndicesToKeep).eval();
243}
244
245} // namespace NAV
246
247#ifndef DOXYGEN_IGNORE
248
249// clang-format off
250
251template<typename T>
252 requires std::is_base_of_v<Eigen::DenseBase<T>, T>
253struct fmt::formatter<T> : ostream_formatter
254{};
255
256// FIXME: This is not compiling with gcc 11.3 but with >12.1.
257// template<typename T>
258// requires std::is_base_of_v<Eigen::QuaternionBase<T>, T>
259// struct fmt::formatter<T> : ostream_formatter
260// {};
261template<>
262struct fmt::formatter<Eigen::Quaternionf> : ostream_formatter
263{};
264template<>
265struct fmt::formatter<Eigen::Quaterniond> : ostream_formatter
266{};
267template<>
268struct fmt::formatter<Eigen::Quaternionld> : ostream_formatter
269{};
270
271// clang-format on
272
273#endif
Assertion helpers.
#define INS_ASSERT_USER_ERROR(_EXP, _MSG)
Assert function with message.
Definition Assert.h:21
Quaternion< long double > Quaternionld
Long double Eigen::Quaternion.
Definition Eigen.hpp:44
Array< double, 5, 1 > Array5d
Double 5x1 Eigen::Array.
Definition Eigen.hpp:29
Array< double, 6, 1 > Array6d
Double 6x1 Eigen::Array.
Definition Eigen.hpp:30
Matrix< double, 6, 6 > Matrix6d
Double 6x6 Eigen::Matrix.
Definition Eigen.hpp:34
Matrix< double, 5, 5 > Matrix5d
Double 5x5 Eigen::Matrix.
Definition Eigen.hpp:33
Matrix< long double, 3, 3 > Matrix3ld
Long double 3x3 Eigen::Matrix.
Definition Eigen.hpp:41
void removeRows(Eigen::DenseBase< Derived > &matrix, size_t index, size_t length)
Removes rows from a matrix or vector.
Definition Eigen.hpp:123
Matrix< long double, 3, 1 > Vector3ld
Long double 3x1 Eigen::Vector.
Definition Eigen.hpp:38
Matrix< double, 5, 1 > Vector5d
Double 5x1 Eigen::Vector.
Definition Eigen.hpp:31
Matrix< double, 6, 1 > Vector6d
Double 6x1 Eigen::Vector.
Definition Eigen.hpp:32
Matrix< long double, 4, 4 > Matrix4ld
Long double 4x4 Eigen::Matrix.
Definition Eigen.hpp:42
Array< long double, 3, 1 > Array3ld
Long double 3x1 Eigen::Array.
Definition Eigen.hpp:36
void removeCols(Eigen::DenseBase< Derived > &matrix, size_t index, size_t length)
Removes columns from a matrix or vector.
Definition Eigen.hpp:159
Matrix< long double, 4, 1 > Vector4ld
Long double 3x1 Eigen::Vector.
Definition Eigen.hpp:39
AngleAxis< long double > AngleAxisld
Long double Eigen::AngleAxis.
Definition Eigen.hpp:46
void removeRowsAndCols(Eigen::DenseBase< Derived > &matrix, size_t row, size_t rows, size_t col, size_t cols)
Removes rows and columns from a matrix or vector.
Definition Eigen.hpp:197
nlohmann::json json
json namespace
Definition FlowManager.hpp:21
Utility class for logging to console and file.
#define LOG_WARN
Error occurred, but a fallback option exists and program continues to work normally.
Definition Logger.hpp:71