DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
EigenUtil.hpp
Go to the documentation of this file.
1#pragma once
2/// @file EigenUtil.hpp
3/// @brief Eigen extensions: `to_string`, an fmt-safe wrapper, and fmt formatter
4/// specialisations for dense Eigen matrices.
5///
6/// The @ref DNDS::MatrixFMTSafe "MatrixFMTSafe" / @ref VectorFMTSafe helpers exist because modern fmtlib
7/// (with `fmt/ranges.h`) will detect `Eigen::Matrix` as a range and format
8/// it as `[a, b, c, ...]`, overriding the Eigen stream formatting that
9/// DNDSR wants. Wrapping Eigen types with these aliases hides the iterator
10/// interface from fmt and keeps the matrix pretty-print.
11
12#include "EigenPCH.hpp"
13
14#include "Defines.hpp"
15
16#include <fmt/core.h>
17#include <fmt/ostream.h>
18#include <fmt/format.h>
19#include "DeviceStorage.hpp"
20#include "Vector.hpp"
21
22namespace DNDS
23{
24 /// @brief Render an `Eigen::DenseBase` to a string via `operator<<`.
25 /// @param precision Setprecision applied if `> 0`.
26 /// @param scientific Whether to use `std::scientific` notation.
27 // TODO: lessen copying chance?
28 template <class dir>
29 std::string to_string(const Eigen::DenseBase<dir> &v,
30 int precision = 5,
31 bool scientific = false)
32 {
33 std::stringstream ss;
34 if (precision > 0)
35 ss << std::setprecision(precision);
36 if (scientific)
37 ss << std::scientific;
38 ss << v;
39 return ss.str();
40 }
41
42}
43
44namespace Eigen
45{
46 /**
47 * @brief `Eigen::Matrix` wrapper that hides `begin`/`end` from fmt.
48 *
49 * @details When both `fmt/ranges.h` and Eigen are present, the fmt range
50 * formatter picks up `Eigen::Matrix` via its iterator interface and renders
51 * it as a bracketed list (`[1, 2, 3]`). This disables Eigen's stream
52 * formatting path. This wrapper inherits from `Eigen::Matrix` but deletes
53 * `begin` / `end`, so fmt falls back to the ostream formatter.
54 *
55 * Use this type (or its @ref VectorFMTSafe / @ref RowVectorFMTSafe aliases) wherever
56 * Eigen objects need to pass through `fmt::format`.
57 */
58 template <class T, int M, int N, int options = AutoAlign | ((M == 1 && N != 1) ? Eigen ::RowMajor : !(M == 1 && N != 1) ? Eigen ::ColMajor
59 : Eigen ::ColMajor),
60 int max_m = M, int max_n = N>
61 struct MatrixFMTSafe : public Matrix<T, M, N, options, max_m, max_n>
62 {
63 using Base = Matrix<T, M, N, options, max_m, max_n>;
64 using Base::Base;
65
66 /// @brief Deleted to hide range interface from fmt.
67 void begin() = delete;
68 /// @brief Deleted to hide range interface from fmt.
69 void end() = delete;
70 };
71
72 /// @brief Column-vector alias of @ref DNDS::MatrixFMTSafe "MatrixFMTSafe".
73 template <class T, int M>
75
76 /// @brief Row-vector alias of @ref DNDS::MatrixFMTSafe "MatrixFMTSafe".
77 template <class T, int N>
79}
80
81namespace DNDS::Meta
82{
83 /// @brief Type trait: is `T` a @ref DNDS::MatrixFMTSafe "MatrixFMTSafe" with real scalar? Used by
84 /// the fmt formatter below to catch both wrapped and unwrapped Eigen types.
85 template <class T>
86 struct is_real_eigen_fmt_safe_matrix : public std::false_type
87 {
88 };
89
90 template <int M, int N, int options, int max_m, int max_n>
94
95 template <class T>
97
98 const bool v = is_real_eigen_fmt_safe_matrix_v<Eigen::MatrixFMTSafe<real, 3, 3>>;
99}
100
101// formatter support for dense eigen matrices
102// ! is not compatible with fmt/ranges.h
103// ! use Eigen::MatrixFMTSafe if fmt/ranges.h is present
104// ! Eigen::Vector s would be fine if fmt/ranges.h is present, but using fmt/ranges.h syntax
105template <typename T, typename Char>
106struct fmt::formatter<T, Char,
107 std::enable_if_t<DNDS::Meta::is_eigen_dense_v<std::remove_cv_t<T>> ||
108 DNDS::Meta::is_real_eigen_fmt_safe_matrix_v<std::remove_cv_t<T>>>>
109// template <int M, int N, int options, int max_m, int max_n, class Char>
110// struct fmt::formatter<Eigen::Matrix<DNDS::real, M, N, options, max_m, max_n>, Char>
111{
112 // using TMat = Eigen::Matrix<DNDS::real, M, N, options, max_m, max_n>;
113 using TMat = std::remove_cv_t<T>;
114 char align = '>';
115 char sign = ' ';
116 int width = -1;
117 int precision = 16;
118 char type = 'g';
119 std::string formatSpecC = "{}";
120
121 auto parse(fmt::format_parse_context &ctx)
122 {
123 auto it = ctx.begin(), end = ctx.end();
124 bool afterDot = false;
125 while (it != end && *it != '}')
126 { // a home-cooked version of float point format parser
127 switch (*it)
128 {
129 case '<':
130 case '>':
131 case '^':
132 align = *it++;
133 break;
134 case '+':
135 case '-':
136 case ' ':
137 sign = *it++;
138 break;
139 case 'e':
140 case 'E':
141 case 'f':
142 case 'F':
143 case 'g':
144 case 'G':
145 type = *it++;
146 break;
147 case '.':
148 afterDot = true;
149 it++;
150 break;
151 default:
152 {
153 if (*it >= '0' && *it <= '9')
154 {
155 std::string v;
156 v.reserve(20);
157 while (it != end && *it >= '0' && *it <= '9')
158 v.push_back(*it++);
159 if (afterDot)
160 precision = std::stoi(v);
161 else
162 width = std::stoi(v);
163 }
164 else
165 DNDS_assert_info(false, fmt::format("invalid char {}", *it));
166 }
167 break;
168 }
169 }
170 if (width == -1)
171 formatSpecC = fmt::format(FMT_STRING("{{:{0}{1}.{3}{4}}}"), align, sign, width, precision, type);
172 else
173 formatSpecC = fmt::format(FMT_STRING("{{:{0}{1}{2}.{3}{4}}}"), align, sign, width, precision, type);
174 return it;
175 }
176
177 auto format(const TMat &mat, fmt::format_context &ctx) const
178 {
179 std::string buf;
180 buf.reserve(mat.size() * 10);
181 buf.push_back('[');
182 for (Eigen::Index i = 0; i < mat.rows(); ++i)
183 {
184 for (Eigen::Index j = 0; j < mat.cols(); ++j)
185 {
186 if (j > 0)
187 buf.append(",");
188 fmt::format_to(std::back_inserter(buf), formatSpecC, mat(i, j));
189 }
190 if (i < mat.rows() - 1)
191 buf.append(";\n");
192 }
193 buf.push_back(']');
194 return fmt::format_to(ctx.out(), "{}", buf);
195 }
196};
197
198namespace DNDS
199{
200
201 template <DeviceBackend B, typename t_scalar, int M, int N>
203 {
204 public:
205 static_assert(M >= 0 || M == Eigen::Dynamic, "M needs to be a valid eigen size");
206 static_assert(N >= 0 || N == Eigen::Dynamic, "N needs to be a valid eigen size");
207 using t_matrix = Eigen::Matrix<std::remove_cv_t<t_scalar>, M, N>;
208 using t_map_const = Eigen::Map<const t_matrix>;
209 using t_map = std::conditional_t<std::is_const_v<t_scalar>,
210 t_map_const, Eigen::Map<t_matrix>>;
212
213 protected:
214 t_scalar *ptr = nullptr;
215 std::conditional_t<M >= 0, EmptyNoDefault, rowsize> M_dynamic = 0;
216 std::conditional_t<N >= 0, EmptyNoDefault, rowsize> N_dynamic = 0;
217
218 public:
220
222 : ptr(n_ptr), M_dynamic(m), N_dynamic(n)
223 {
224 }
225
226 DNDS_DEVICE_CALLABLE [[nodiscard]] t_scalar *data() const
227 {
228 return ptr;
229 }
230
231 DNDS_DEVICE_CALLABLE [[nodiscard]] rowsize rows() const
232 {
233 if constexpr (M >= 0)
234 return M;
235 else
236 return M_dynamic;
237 }
238
239 DNDS_DEVICE_CALLABLE [[nodiscard]] rowsize cols() const
240 {
241 if constexpr (N >= 0)
242 return N;
243 else
244 return N_dynamic;
245 }
246
248 {
249 return {ptr, rows(), cols()};
250 }
251
253 {
254 return {ptr, rows(), cols()};
255 }
256 };
257
258 static_assert(std::is_trivially_copyable_v<EigenMatrixView<DeviceBackend::Host, real, Eigen::Dynamic, Eigen::Dynamic>>);
259
260 template <typename t_scalar, int M, int N>
262 {
263 public:
264 using t_matrix = Eigen::Matrix<t_scalar, M, N>;
265 using t_map = Eigen::Map<t_matrix>;
266
267 protected:
271
272 public:
274 {
275 h_data.resize(this->size());
276 }
277
278 rowsize rows() const
279 {
280 if constexpr (M >= 0)
281 return M;
282 else
283 return M_dynamic;
284 }
285
286 rowsize cols() const
287 {
288 if constexpr (N >= 0)
289 return N;
290 else
291 return N_dynamic;
292 }
293
294 rowsize size() const
295 {
296 return rows() * cols();
297 }
298
299 void resize(int m, int n)
300 {
301 if constexpr (M >= 0)
302 DNDS_assert(M == m);
303 else
304 M_dynamic = m;
305 if constexpr (N >= 0)
306 DNDS_assert(N == n);
307 else
308 N_dynamic = n;
309 h_data.resize(this->size());
310 }
311
313 {
314 return {h_data.data(), this->rows(), this->cols()};
315 }
316
317 template <DeviceBackend B>
319
320 void to_host()
321 {
322 h_data.to_host();
323 }
324
326 {
327 h_data.to_device(B);
328 }
329
330 template <DeviceBackend B>
332 {
333 DNDS_assert_info(h_data.device() == B || (h_data.device() == DeviceBackend::Unknown && B == DeviceBackend::Host),
334 "data not on this backend");
335 switch (B)
336 {
339 return t_deviceView<B>{h_data.data(), rows(), cols()};
340#ifdef DNDS_USE_CUDA
341 case DeviceBackend::CUDA:
342 return t_deviceView<B>{h_data.dataDevice(), rows(), cols()};
343#endif
344 default:
345 DNDS_assert_info(false, std::string("this device not implemented -- ") + device_backend_name(B));
346 return t_deviceView<B>{h_data.dataDevice(), rows(), cols()};
347 }
348 }
349 };
350}
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE(T, T_Self)
Definition Defines.hpp:83
#define DNDS_DEVICE_CALLABLE
Definition Defines.hpp:76
Device memory abstraction layer with backend-specific storage and factory creation.
Pre-compiled-header style shim that includes the heavy Eigen headers under DNDSR's warning suppressio...
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:113
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
Host-device vector types with optional GPU storage and device-side views.
std::conditional_t< std::is_const_v< t_scalar >, t_map_const, Eigen::Map< t_matrix > > t_map
DNDS_DEVICE_CALLABLE rowsize cols() const
DNDS_DEVICE_CALLABLE rowsize rows() const
DNDS_DEVICE_CALLABLE t_map_const map() const
Eigen::Map< const t_matrix > t_map_const
DNDS_DEVICE_CALLABLE t_scalar * data() const
Eigen::Matrix< std::remove_cv_t< t_scalar >, M, N > t_matrix
DNDS_DEVICE_CALLABLE t_map map()
void resize(int m, int n)
t_deviceView< B > deviceView()
Eigen::Map< t_matrix > t_map
void to_device(DeviceBackend B)
host_device_vector< t_scalar > h_data
Eigen::Matrix< t_scalar, M, N > t_matrix
EigenMatrixView< B, t_scalar, M, N > t_deviceView
constexpr bool is_real_eigen_fmt_safe_matrix_v
Definition EigenUtil.hpp:96
the host side operators are provided as implemented
DeviceBackend
Enumerates the backends a DeviceStorage / Array can live on.
@ Unknown
Unset / sentinel.
@ Host
Plain CPU memory.
const char * device_backend_name(DeviceBackend B)
Canonical string name for a DeviceBackend (used in log messages).
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:109
std::string to_string(const Eigen::DenseBase< dir > &v, int precision=5, bool scientific=false)
Render an Eigen::DenseBase to a string via operator<<.
Definition EigenUtil.hpp:29
Type trait: is T a MatrixFMTSafe with real scalar? Used by the fmt formatter below to catch both wrap...
Definition EigenUtil.hpp:87
Eigen::Matrix wrapper that hides begin/end from fmt.
Definition EigenUtil.hpp:62
void begin()=delete
Deleted to hide range interface from fmt.
Matrix< T, M, N, options, max_m, max_n > Base
Definition EigenUtil.hpp:63
void end()=delete
Deleted to hide range interface from fmt.
Eigen::Matrix< real, 5, 1 > v
constexpr DNDS::index N
Eigen::Vector3d n(1.0, 0.0, 0.0)