OpenKalman
eigen-library-interface.hpp
Go to the documentation of this file.
1 /* This file is part of OpenKalman, a header-only C++ library for
2  * Kalman filters and other recursive filters.
3  *
4  * Copyright (c) 2019-2024 Christopher Lee Ogden <ogden@gatech.edu>
5  *
6  * This Source Code Form is subject to the terms of the Mozilla Public
7  * License, v. 2.0. If a copy of the MPL was not distributed with this
8  * file, You can obtain one at https://mozilla.org/MPL/2.0/.
9  */
10 
16 #ifndef OPENKALMAN_EIGEN_LIBRARY_INTERFACE_HPP
17 #define OPENKALMAN_EIGEN_LIBRARY_INTERFACE_HPP
18 
19 #include <type_traits>
20 #include <tuple>
21 #include <random>
23 
24 namespace OpenKalman::interface
25 {
26 #ifdef __cpp_concepts
27  template<Eigen3::eigen_general<true> T>
28  struct library_interface<T>
29 #else
30  template<typename T>
31  struct library_interface<T, std::enable_if_t<Eigen3::eigen_general<T, true>>>
32 #endif
33  {
34  template<typename Derived>
35  using library_base = Eigen3::EigenAdapterBase<Derived,
36  std::conditional_t<Eigen3::eigen_array_general<T>, Eigen::ArrayBase<Derived>,
37  std::conditional_t<Eigen3::eigen_matrix_general<T>, Eigen::MatrixBase<Derived>, Eigen::EigenBase<Derived>>>>;
38 
39  private:
40 
41  template<typename Arg, std::size_t...I, typename...Ind>
42  static constexpr void
43  check_index_bounds(const Arg& arg, std::index_sequence<I...>, Ind...ind)
44  {
45  ([](auto max, auto ind){ if (ind < 0 or ind >= max) throw std::out_of_range {
46  ("Index " + std::to_string(I) + " is out of bounds: it is " + std::to_string(ind) +
47  " but should be in range [0..." + std::to_string(max - 1) + "].")};
48  }(get_index_dimension_of<I>(arg), ind),...);
49  }
50 
51 
52  template<typename Arg>
53  static constexpr decltype(auto)
54  get_coeff(Arg&& arg, Eigen::Index i, Eigen::Index j)
55  {
57  //check_index_bounds(arg, std::make_index_sequence<2> {}, i, j);
58  using Scalar = scalar_type_of_t<Arg>;
59 
60  if constexpr (Eigen3::eigen_DiagonalMatrix<Arg> or Eigen3::eigen_DiagonalWrapper<Arg>)
61  {
62  if (i == j)
63  return static_cast<Scalar>(get_coeff(nested_object(std::forward<Arg>(arg)), i, 0));
64  else
65  return static_cast<Scalar>(0);
66  }
67  else if constexpr (Eigen3::eigen_TriangularView<Arg>)
68  {
69  constexpr int Mode {std::decay_t<Arg>::Mode};
70  if ((i > j and (Mode & int{Eigen::Upper}) != 0x0) or (i < j and (Mode & int{Eigen::Lower}) != 0x0))
71  return static_cast<Scalar>(0);
72  else
73  return static_cast<Scalar>(get_coeff(nested_object(std::forward<Arg>(arg)), i, j));
74  }
75  else if constexpr (Eigen3::eigen_SelfAdjointView<Arg>)
76  {
77  return get_coeff(nested_object(std::forward<Arg>(arg)), i, j);
78  }
79  else
80  {
81  auto evaluator = Eigen::internal::evaluator<std::decay_t<Arg>>(arg);
82  constexpr bool use_coeffRef = (Eigen::internal::traits<std::decay_t<Arg>>::Flags & Eigen::LvalueBit) != 0x0 and
83  std::is_lvalue_reference_v<Arg> and not std::is_const_v<std::remove_reference_t<Arg>>;
84  if constexpr (use_coeffRef) return evaluator.coeffRef(i, j);
85  else return evaluator.coeff(i, j);
86  }
87  }
88 
89 
90  template<typename Indices>
91  static constexpr std::tuple<Eigen::Index, Eigen::Index>
92  extract_indices(const Indices& indices)
93  {
94  auto it = stdex::ranges::begin(indices);
95  auto e = stdex::ranges::end(indices);
96  if (it == e) return {0, 0};
97  auto i = static_cast<std::size_t>(*it);
98  if (++it == e) return {i, 0};
99  auto j = static_cast<std::size_t>(*it);
100  if (++it == e) return {i, j};
101  throw std::logic_error("Wrong number of indices on component access");
102  }
103 
104  public:
105 
106 #ifdef __cpp_lib_ranges
107  template<typename Arg, std::ranges::input_range Indices> requires
108  std::convertible_to<std::ranges::range_value_t<Indices>, const typename std::decay_t<Arg>::Index> and
109  (collections::size_of_v<Indices> == stdex::dynamic_extent or collections::size_of_v<Indices> <= 2)
110  static constexpr values::scalar decltype(auto)
111 #else
112  template<typename Arg, typename Indices, std::enable_if_t<
113  (collections::size_of_v<Indices> == stdex::dynamic_extent or collections::size_of_v<Indices> <= 2), int> = 0>
114  static constexpr decltype(auto)
115 #endif
116  get_component(Arg&& arg, const Indices& indices)
117  {
118  using Scalar = scalar_type_of_t<Arg>;
119  auto [i, j] = extract_indices(indices);
120  if constexpr (Eigen3::eigen_SelfAdjointView<Arg>)
121  {
122  constexpr int Mode {std::decay_t<Arg>::Mode};
123  bool transp = (i > j and (Mode & int{Eigen::Upper}) != 0) or (i < j and (Mode & int{Eigen::Lower}) != 0);
124  if constexpr (values::complex<Scalar>)
125  {
126  if (transp) return static_cast<Scalar>(values::conj(get_coeff(std::as_const(arg), j, i)));
127  else return static_cast<Scalar>(get_coeff(std::as_const(arg), i, j));
128  }
129  else
130  {
131  if (transp) return get_coeff(std::forward<Arg>(arg), j, i);
132  else return get_coeff(std::forward<Arg>(arg), i, j);
133  }
134  }
135  else
136  {
137  return get_coeff(std::forward<Arg>(arg), i, j);
138  }
139  }
140 
141 
142 #ifdef __cpp_lib_ranges
143  template<typename Arg, std::ranges::input_range Indices> requires (not std::is_const_v<Arg>) and
144  std::convertible_to<std::ranges::range_value_t<Indices>, const typename Arg::Index> and
145  (collections::size_of_v<Indices> == stdex::dynamic_extent or collections::size_of_v<Indices> <= 2) and
146  std::assignable_from<decltype(get_coeff(std::declval<Arg&>(), 0, 0)), const scalar_type_of_t<Arg>&> and
147  ((Eigen::internal::traits<std::decay_t<Arg>>::Flags & Eigen::LvalueBit) != 0) and
148  (not Eigen3::eigen_DiagonalWrapper<Arg>) and (not Eigen3::eigen_TriangularView<Arg>)
149 #else
150  template<typename Arg, typename Indices, std::enable_if_t<(not std::is_const_v<Arg>) and
151  (collections::size_of_v<Indices> == stdex::dynamic_extent or collections::size_of_v<Indices> <= 2) and
152  std::is_assignable<decltype(get_coeff(std::declval<Arg&>(), 0, 0)), const scalar_type_of_t<Arg>&>::value and
153  (Eigen::internal::traits<std::decay_t<Arg>>::Flags & Eigen::LvalueBit) != 0 and
154  (not Eigen3::eigen_DiagonalWrapper<Arg>) and (not Eigen3::eigen_TriangularView<Arg>), int> = 0>
155 #endif
156  static void
157  set_component(Arg& arg, const scalar_type_of_t<Arg>& s, const Indices& indices)
158  {
159  using Scalar = scalar_type_of_t<Arg>;
160  auto [i, j] = extract_indices(indices);
161  if constexpr (Eigen3::eigen_SelfAdjointView<Arg>)
162  {
163  constexpr int Mode {std::decay_t<Arg>::Mode};
164  bool transp = (i > j and (Mode & int{Eigen::Upper}) != 0) or (i < j and (Mode & int{Eigen::Lower}) != 0);
165  if constexpr (values::complex<Scalar>)
166  {
167  if (transp) get_coeff(arg, j, i) = values::conj(s);
168  else get_coeff(arg, i, j) = s;
169  }
170  else
171  {
172  if (transp) get_coeff(arg, j, i) = s;
173  else get_coeff(arg, i, j) = s;
174  }
175  }
176  else
177  {
178  get_coeff(arg, i, j) = s;
179  }
180  }
181 
182  private:
183 
184  template<typename Arg>
185  static decltype(auto)
186  wrap_if_nests_by_reference(Arg&& arg)
187  {
188  if constexpr (Eigen3::eigen_general<Arg>)
189  {
190  constexpr auto Flags = Eigen::internal::traits<std::remove_reference_t<Arg>>::Flags;
191  if constexpr (std::is_lvalue_reference_v<Arg> and static_cast<bool>(Flags & Eigen::NestByRefBit))
192  return std::forward<Arg>(arg);
193  else
194  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg));
195  }
196  else
197  {
198  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg));
199  }
200  }
201 
202  public:
203 
204  template<typename Arg>
205  static decltype(auto)
206  to_native_matrix(Arg&& arg)
207  {
208  if constexpr (Eigen3::eigen_wrapper<Arg>)
209  {
210  return std::forward<Arg>(arg);
211  }
212  else if constexpr (internal::library_wrapper<Arg>)
213  {
214  return to_native_matrix(OpenKalman::nested_object(std::forward<Arg>(arg)));
215  }
216  else if constexpr (Eigen3::eigen_ArrayWrapper<Arg>)
217  {
218  return to_native_matrix(nested_object(std::forward<Arg>(arg)));
219  }
220  else if constexpr (Eigen3::eigen_array_general<Arg>)
221  {
222  return wrap_if_nests_by_reference(std::forward<Arg>(arg)).matrix();
223  }
224  else if constexpr (Eigen3::eigen_matrix_general<Arg>)
225  {
226  return wrap_if_nests_by_reference(std::forward<Arg>(arg));
227  }
228  else if constexpr (not Eigen3::eigen_general<Arg> and directly_accessible<Arg> and std::is_lvalue_reference_v<Arg>)
229  {
230  using Scalar = scalar_type_of_t<Arg>;
231  constexpr int rows = dynamic_dimension<Arg, 0> ? Eigen::Dynamic : index_dimension_of_v<Arg, 0>;
232  constexpr int cols = dynamic_dimension<Arg, 1> ? Eigen::Dynamic : index_dimension_of_v<Arg, 1>;
233 
234  if constexpr (layout_of_v<Arg> == data_layout::stride)
235  {
236  auto [s0, s1] = internal::strides(arg);
237  using S0 = decltype(s0);
238  using S1 = decltype(s1);
239  constexpr int es0 = []() -> int {
240  if constexpr (values::fixed<S0>) return static_cast<std::ptrdiff_t>(S0{});
241  else return Eigen::Dynamic;
242  }();
243  constexpr int es1 = []() -> int {
244  if constexpr (values::fixed<S1>) return static_cast<std::ptrdiff_t>(S1{});
245  else return Eigen::Dynamic;
246  }();
247  using IndexType = typename std::decay_t<Arg>::Index;
248  auto is0 = static_cast<IndexType>(static_cast<std::ptrdiff_t>(s0));
249  auto is1 = static_cast<IndexType>(static_cast<std::ptrdiff_t>(s1));
250 
251  if constexpr (values::fixed<S0> and
252  (es0 == 1 or (values::fixed<S1> and es0 < es1)))
253  {
254  using M = Eigen::Matrix<Scalar, rows, cols, Eigen::ColMajor>;
255  Eigen::Stride<es1, es0> strides {is1, is0};
256  return Eigen::Map<const M, Eigen::Unaligned, decltype(strides)> {internal::raw_data(arg), rows, cols, strides};
257  }
258  else if constexpr (values::fixed<S1> and
259  (es1 == 1 or (values::fixed<S0> and es1 < es0)))
260  {
261  using M = Eigen::Matrix<Scalar, rows, cols, Eigen::RowMajor>;
262  Eigen::Stride<es0, es1> strides {is0, is1};
263  return Eigen::Map<const M, Eigen::Unaligned, decltype(strides)> {internal::raw_data(arg), rows, cols, strides};
264  }
265  else
266  {
267  if (is1 < is0)
268  {
269  using M = Eigen::Matrix<Scalar, rows, cols, Eigen::RowMajor>;
270  Eigen::Stride<es0, es1> strides {is0, is1};
271  return Eigen::Map<const M, Eigen::Unaligned, decltype(strides)> {internal::raw_data(arg), rows, cols, strides};
272  }
273  else
274  {
275  using M = Eigen::Matrix<Scalar, rows, cols, Eigen::ColMajor>;
276  Eigen::Stride<es1, es0> strides {is1, is0};
277  return Eigen::Map<const M, Eigen::Unaligned, decltype(strides)> {internal::raw_data(arg), rows, cols, strides};
278  }
279  }
280  }
281  else
282  {
283  constexpr auto l = layout_of_v<Arg> == data_layout::right ? Eigen::RowMajor : Eigen::ColMajor;
284  using M = Eigen::Matrix<Scalar, rows, cols, l>;
285  return Eigen::Map<const M> {internal::raw_data(arg), rows, cols};
286  }
287  }
288  else
289  {
290  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg));
291  }
292  }
293 
294  private:
295 
296  template<typename Scalar, int rows, int cols, auto options>
297  using dense_type = std::conditional_t<Eigen3::eigen_array_general<T, true>,
298  Eigen::Array<Scalar, rows, cols, options>, Eigen::Matrix<Scalar, rows, cols, options>>;
299 
300  template<typename Scalar, std::size_t rows, std::size_t cols, auto options>
301  using writable_type = dense_type<Scalar,
302  (rows == stdex::dynamic_extent ? Eigen::Dynamic : static_cast<int>(rows)),
303  (cols == stdex::dynamic_extent ? Eigen::Dynamic : static_cast<int>(cols)), options>;
304 
305  public:
306 
307 #ifdef __cpp_concepts
308  template<typename To, Eigen3::eigen_general From> requires (std::assignable_from<To&, From&&>)
309 #else
310  template<typename To, typename From, std::enable_if_t<Eigen3::eigen_general<From> and std::is_assignable_v<To&, From&&>, int> = 0>
311 #endif
312  static void copy(To& a, From&& b)
313  {
314  if constexpr (Eigen3::eigen_DiagonalWrapper<From>)
315  {
316  if constexpr (vector<nested_object_of_t<From>>)
317  a = std::forward<From>(b);
318  else
319  a = diagonal_of(std::forward<From>(b)).asDiagonal();
320  }
321  else
322  {
323  a = std::forward<From>(b);
324  }
325  }
326 
327  private:
328 
329  template<typename Descriptors>
330  static constexpr decltype(auto)
331  extract_descriptors(Descriptors&& descriptors)
332  {
333  if constexpr (pattern_collection<Descriptors>)
334  {
335  constexpr auto dim = collections::size_of_v<Descriptors>;
336  static_assert(dim <= 2);
337  if constexpr (dim == 0) return std::tuple {patterns::Axis{}, patterns::Axis{}};
338  else if constexpr (dim == 1) return std::tuple {std::get<0>(std::forward<Descriptors>(descriptors)), patterns::Axis{}};
339  else return std::forward<Descriptors>(descriptors);
340  }
341  else
342  {
343  auto it = stdex::ranges::begin(descriptors);
344  auto e = stdex::ranges::end(descriptors);
345  if (it == e) return std::tuple {patterns::Axis{}, patterns::Axis{}};
346  auto i = *it;
347  if (++it == e) return std::tuple {i, patterns::Axis{}};
348  auto j = *it;
349  if (++it == e) return std::tuple {i, j};
350  throw std::logic_error("Wrong number of vector space descriptors");
351  }
352  }
353 
354  public:
355 
356 #ifdef __cpp_concepts
357  template<data_layout layout, values::number Scalar, patterns::euclidean_pattern_collection Ds>
358 #else
359  template<data_layout layout, typename Scalar, typename Ds, std::enable_if_t<values::number<Scalar> and
360  patterns::euclidean_pattern_collection<Ds>, int> = 0>
361 #endif
362  static auto make_default(Ds&& ds)
363  {
364  using IndexType = typename Eigen::Index;
365 
366  if constexpr (pattern_collection<Ds> and collections::size_of_v<Ds> > 2)
367  {
368  constexpr auto options = layout == data_layout::right ? Eigen::RowMajor : Eigen::ColMajor;
369 
370  if constexpr (fixed_pattern_collection<Ds>)
371  {
372  auto sizes = std::apply([](auto&&...d){
373  return Eigen::Sizes<static_cast<std::ptrdiff_t>(patterns::dimension_of_v<decltype(d)>)...> {};
374  }, std::forward<Ds>(ds));
375 
376  return Eigen::TensorFixedSize<Scalar, decltype(sizes), options, IndexType> {};
377  }
378  else
379  {
380  return std::apply([](auto&&...d){
381  using Ten = Eigen::Tensor<Scalar, collections::size_of_v<Ds>, options, IndexType>;
382  return Ten {static_cast<IndexType>(get_dimension(d))...};
383  }, std::forward<Ds>(ds));
384  }
385  }
386  else
387  {
388  auto [d0, d1] = extract_descriptors(std::forward<Ds>(ds));
389  constexpr auto dim0 = patterns::dimension_of_v<decltype(d0)>;
390  constexpr auto dim1 = patterns::dimension_of_v<decltype(d1)>;
391 
392  static_assert(layout != data_layout::right or dim0 == 1 or dim1 != 1,
393  "Eigen does not allow creation of a row-major column vector.");
394  static_assert(layout != data_layout::left or dim0 != 1 or dim1 == 1,
395  "Eigen does not allow creation of a column-major row vector.");
396 
397  constexpr auto options =
398  layout == data_layout::right or (layout == data_layout::none and dim0 == 1 and dim1 != 1) ?
399  Eigen::RowMajor : Eigen::ColMajor;
400 
401  using M = writable_type<Scalar, dim0, dim1, options>;
402 
403  if constexpr (dim0 == stdex::dynamic_extent or dim1 == stdex::dynamic_extent)
404  return M {static_cast<IndexType>(get_dimension(d0)), static_cast<IndexType>(get_dimension(d1))};
405  else
406  return M {};
407  }
408 
409 
410  }
411 
412 
413 #ifdef __cpp_concepts
414  template<data_layout layout, writable Arg, std::convertible_to<scalar_type_of_t<Arg>> S, std::convertible_to<scalar_type_of_t<Arg>>...Ss>
415  requires (layout == data_layout::right) or (layout == data_layout::left)
416 #else
417  template<data_layout layout, typename Arg, typename S, typename...Ss, std::enable_if_t<writable<Arg> and
418  (layout == data_layout::right or layout == data_layout::left) and
419  std::conjunction<std::is_convertible<S, typename scalar_type_of<Arg>::type>,
420  std::is_convertible<Ss, typename scalar_type_of<Arg>::type>...>::value, int> = 0>
421 #endif
422  static void fill_components(Arg& arg, const S s, const Ss ... ss)
423  {
424  if constexpr (layout == data_layout::left) ((arg.transpose() << s), ... , ss);
425  else ((arg << s), ... , ss);
426  }
427 
428 
429 #ifdef __cpp_lib_ranges
430  template<values::dynamic C, patterns::euclidean_pattern_collection Ds> requires
431  (collections::size_of_v<Ds> != stdex::dynamic_extent) and (collections::size_of_v<Ds> <= 2)
432  static constexpr constant_matrix auto
433 #else
434  template<typename C, typename Ds, std::enable_if_t<values::dynamic<C> and
435  patterns::euclidean_pattern_collection<Ds> and
436  (collections::size_of_v<Ds> != stdex::dynamic_extent) and (collections::size_of_v<Ds> <= 2)), int> = 0>
437  static constexpr auto
438 #endif
439  make_constant(C&& c, Ds&& ds)
440  {
441  auto [d0, d1] = extract_descriptors(std::forward<Ds>(ds));
442  constexpr auto dim0 = patterns::dimension_of_v<decltype(d0)>;
443  constexpr auto dim1 = patterns::dimension_of_v<decltype(d1)>;
444 
445  auto value = values::to_value_type(std::forward<C>(c));
446  using Scalar = std::decay_t<decltype(value)>;
447  constexpr auto options = (dim0 == 1 and dim1 != 1) ? Eigen::RowMajor : Eigen::ColMajor;
448  using M = writable_type<Scalar, dim0, dim1, options>;
449 
450  using IndexType = typename Eigen::Index;
451 
452  if constexpr (fixed_pattern_collection<Ds>)
453  return M::Constant(value);
454  else
455  return M::Constant(static_cast<IndexType>(dim0), static_cast<IndexType>(get_dimension(d1)), value);
456  // Note: We don't address orders greater than 2 because Eigen::Tensor's constant has substantial limitations.
457  }
458 
459 
460 #ifdef __cpp_concepts
461  template<typename Scalar, patterns::euclidean_pattern_collection Ds> requires
462  (collections::size_of_v<Ds> != stdex::dynamic_extent) and (collections::size_of_v<Ds> <= 2)
463  static constexpr identity_matrix auto
464 #else
465  template<typename Scalar, typename...Ds, std::enable_if_t<patterns::euclidean_pattern_collection<Ds> and
466  (collections::size_of_v<Ds> != stdex::dynamic_extent) and (collections::size_of_v<Ds> <= 2), int> = 0>
467  static constexpr auto
468 #endif
469  make_identity_matrix(Ds&& ds)
470  {
471  auto [d0, d1] = extract_descriptors(std::forward<Ds>(ds));
472  constexpr auto dim0 = patterns::dimension_of_v<decltype(d0)>;
473  constexpr auto dim1 = patterns::dimension_of_v<decltype(d1)>;
474 
475  constexpr auto options = (dim0 == 1 and dim1 != 1) ? Eigen::RowMajor : Eigen::ColMajor;
476  using M = writable_type<Scalar, dim0, dim1, options>;
477 
478  using IndexType = typename Eigen::Index;
479 
480  if constexpr (fixed_pattern_collection<Ds>)
481  return M::Identity();
482  else
483  return M::Identity(static_cast<IndexType>(get_dimension(d0)), static_cast<IndexType>(get_dimension(d1)));
484  }
485 
486 
487 #ifdef __cpp_concepts
488  template<triangle_type t, Eigen3::eigen_dense_general Arg> requires std::is_lvalue_reference_v<Arg> or
489  (not std::is_lvalue_reference_v<typename Eigen::internal::ref_selector<std::decay_t<Arg>>::non_const_type>)
490 #else
491  template<triangle_type t, typename Arg, std::enable_if_t<Eigen3::eigen_dense_general<Arg> and (std::is_lvalue_reference_v<Arg> or
492  not std::is_lvalue_reference_v<typename Eigen::internal::ref_selector<std::remove_reference_t<Arg>>::non_const_type>), int> = 0>
493 #endif
494  static constexpr auto
495  make_triangular_matrix(Arg&& arg)
496  {
497  constexpr auto Mode = t == triangle_type::upper ? Eigen::Upper : Eigen::Lower;
498  return arg.template triangularView<Mode>();
499  }
500 
501 
502 #ifdef __cpp_concepts
503  template<HermitianAdapterType t, Eigen3::eigen_dense_general Arg> requires std::is_lvalue_reference_v<Arg>
504 #else
505  template<HermitianAdapterType t, typename Arg, std::enable_if_t<Eigen3::eigen_dense_general<Arg> and std::is_lvalue_reference_v<Arg>, int> = 0>
506 #endif
507  static constexpr auto
508  make_hermitian_adapter(Arg&& arg)
509  {
510  constexpr auto Mode = t == HermitianAdapterType::upper ? Eigen::Upper : Eigen::Lower;
511  return arg.template selfadjointView<Mode>();
512  }
513 
514 
515  // to_euclidean not defined--rely on default
516 
517  // from_euclidean not defined--rely on default
518 
519  // wrap_angles not defined--rely on default
520 
521 
522 #ifdef __cpp_concepts
523  template<typename Arg, typename...Begin, typename...Size> requires (sizeof...(Begin) <= 2)
524 #else
525  template<typename Arg, typename...Begin, typename...Size, std::enable_if_t<(sizeof...(Begin) <= 2), int> = 0>
526 #endif
527  static auto
528  get_slice(Arg&& arg, const std::tuple<Begin...>& begin, const std::tuple<Size...>& size)
529  {
530  auto b0 = [](const auto& begin){
531  using Begin0 = std::decay_t<decltype(std::get<0>(begin))>;
532  if constexpr (values::fixed<Begin0>) return std::integral_constant<Eigen::Index, Begin0::value>{};
533  else return static_cast<Eigen::Index>(std::get<0>(begin));
534  }(begin);
535 
536  auto b1 = [](const auto& begin){
537  if constexpr (sizeof...(Begin) < 2) return std::integral_constant<Eigen::Index, 0>{};
538  else
539  {
540  using Begin1 = std::decay_t<decltype(std::get<1>(begin))>;
541  if constexpr (values::fixed<Begin1>) return std::integral_constant<Eigen::Index, Begin1::value>{};
542  else return static_cast<Eigen::Index>(std::get<1>(begin));
543  }
544  }(begin);
545 
546  auto s0 = [](const auto& size){
547  using Size0 = std::decay_t<decltype(std::get<0>(size))>;
548  if constexpr (values::fixed<Size0>) return std::integral_constant<Eigen::Index, Size0::value>{};
549  else return static_cast<Eigen::Index>(std::get<0>(size));
550  }(size);
551 
552  auto s1 = [](const auto& size){
553  if constexpr (sizeof...(Size) < 2) return std::integral_constant<Eigen::Index, 1>{};
554  else
555  {
556  using Size1 = std::decay_t<decltype(std::get<1>(size))>;
557  if constexpr (values::fixed<Size1>) return std::integral_constant<Eigen::Index, Size1::value>{};
558  else return static_cast<Eigen::Index>(std::get<1>(size));
559  }
560  }(size);
561 
562  constexpr int S0 = static_cast<int>(values::fixed<decltype(s0)> ? static_cast<Eigen::Index>(s0) : Eigen::Dynamic);
563  constexpr int S1 = static_cast<int>(values::fixed<decltype(s1)> ? static_cast<Eigen::Index>(s1) : Eigen::Dynamic);
564 
565  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
566  using M = decltype(m);
567 
568  constexpr auto Flags = Eigen::internal::traits<std::remove_reference_t<M>>::Flags;
569 
570  if constexpr (directly_accessible<M> and not (std::is_lvalue_reference_v<M> and static_cast<bool>(Flags & Eigen::NestByRefBit)))
571  {
572  // workaround for Eigen::Block's special treatment of directly accessible nested types:
573  auto rep {std::forward<M>(m).template replicate<1,1>()};
574  using B = Eigen::Block<const decltype(rep), S0, S1>;
575  if constexpr ((values::fixed<Size> and ...))
576  return B {std::move(rep), std::move(b0), std::move(b1)};
577  else
578  return B {std::move(rep), std::move(b0), std::move(b1), std::move(s0), std::move(s1)};
579  }
580  else
581  {
582  using M_noref = std::remove_reference_t<M>;
583  using XprType = std::conditional_t<std::is_lvalue_reference_v<M>, M_noref, const M_noref>;
584  using B = Eigen::Block<XprType, S0, S1>;
585  if constexpr ((values::fixed<Size> and ...))
586  return B {std::forward<M>(m), std::move(b0), std::move(b1)};
587  else
588  return B {std::forward<M>(m), std::move(b0), std::move(b1), std::move(s0), std::move(s1)};
589  }
590  }
591 
592 
593 #ifdef __cpp_concepts
594  template<Eigen3::eigen_dense_general Arg, Eigen3::eigen_general Block, typename...Begin> requires (sizeof...(Begin) <= 2)
595 #else
596  template<typename Arg, typename Block, typename...Begin, std::enable_if_t<
597  Eigen3::eigen_dense_general<Arg> and Eigen3::eigen_general<Block> and (sizeof...(Begin) <= 2), int> = 0>
598 #endif
599  static void
600  set_slice(Arg& arg, Block&& block, const Begin&...begin)
601  {
602  auto [b0, b1] = [](Eigen::Index bs0, Eigen::Index bs1, const auto&...){ return std::tuple {bs0, bs1}; }
603  (static_cast<std::size_t>(begin)..., 0_uz, 0_uz);
604 
605  if constexpr (Eigen3::eigen_Block<Block>)
606  {
607  if (std::addressof(arg) == std::addressof(block.nestedExpression()) and b0 == block.startRow() and b1 == block.startCol())
608  return;
609  }
610 
611  constexpr auto Bx0 = static_cast<int>(index_dimension_of_v<Block, 0>);
612  constexpr auto Bx1 = static_cast<int>(index_dimension_of_v<Block, 1>);
613  using Bk = Eigen::Block<std::remove_reference_t<Arg>, Bx0, Bx1>;
614 
615  if constexpr (not has_dynamic_dimensions<Block>)
616  {
617  Bk {arg, b0, b1} = std::forward<Block>(block);
618  }
619  else
620  {
621  auto s0 = static_cast<Eigen::Index>(get_index_dimension_of<0>(block));
622  auto s1 = static_cast<Eigen::Index>(get_index_dimension_of<1>(block));
623  Bk {arg, b0, b1, s0, s1} = std::forward<Block>(block);
624  }
625  }
626 
627 
628 #ifdef __cpp_concepts
629  template<triangle_type t, Eigen3::eigen_SelfAdjointView A, Eigen3::eigen_general B> requires (not hermitian_matrix<A>) and
630  set_triangle_defined_for<T, t, decltype(OpenKalman::nested_object(std::declval<A>())), B&&>
631 #else
632  template<triangle_type t, typename A, typename B, std::enable_if_t<
633  Eigen3::eigen_general<B> and Eigen3::eigen_SelfAdjointView<A> and (not hermitian_matrix<A>) and
634  set_triangle_defined_for<T, t, decltype(OpenKalman::nested_object(std::declval<A>())), B&&>, int> = 0>
635 #endif
636  static void
637  set_triangle(A&& a, B&& b)
638  {
639  // A SelfAdjointView won't always be a hermitian_adapter
642  internal::set_triangle<t>(OpenKalman::nested_object(std::forward<A>(a)), OpenKalman::conjugate_transpose(std::forward<B>(b)));
643  else
644  internal::set_triangle<t>(OpenKalman::nested_object(std::forward<A>(a)), std::forward<B>(b));
645  }
646 
647 
648 #ifdef __cpp_concepts
649  template<triangle_type t, typename A, Eigen3::eigen_general B> requires
650  ((diagonal_matrix<A> and internal::has_nested_vector<A>) or t == triangle_type::diagonal) and
651  std::assignable_from<decltype(OpenKalman::diagonal_of(std::declval<A&&>())), decltype(OpenKalman::diagonal_of(std::declval<B&&>()))>
652 #else
653  template<triangle_type t, typename A, typename B, std::enable_if_t<
654  Eigen3::eigen_general<B> and ((diagonal_matrix<A> and internal::has_nested_vector<A>) or t == triangle_type::diagonal) and
655  std::is_assignable_v<decltype(OpenKalman::diagonal_of(std::declval<A&&>())), decltype(OpenKalman::diagonal_of(std::declval<B&&>()))>, int> = 0>
656 #endif
657  static void
658  set_triangle(A&& a, B&& b)
659  {
660  copy(OpenKalman::diagonal_of(std::forward<A>(a)), OpenKalman::diagonal_of(std::forward<B>(b)));
661  }
662 
663 
664 #ifdef __cpp_concepts
665  template<triangle_type t, typename A, Eigen3::eigen_general B> requires
666  (Eigen3::eigen_MatrixWrapper<A> or Eigen3::eigen_ArrayWrapper<A>) and
667  set_triangle_defined_for<T, t, decltype(OpenKalman::nested_object(std::declval<A>())), B&&>
668 #else
669  template<triangle_type t, typename A, typename B, std::enable_if_t<Eigen3::eigen_general<B> and
670  (Eigen3::eigen_MatrixWrapper<A> or Eigen3::eigen_ArrayWrapper<A>) and
671  set_triangle_defined_for<T, t, decltype(OpenKalman::nested_object(std::declval<A>())), B&&>, int> = 0>
672 #endif
673  static void
674  set_triangle(A&& a, B&& b)
675  {
676  internal::set_triangle<t>(OpenKalman::nested_object(std::forward<A>(a)), std::forward<B>(b));
677  }
678 
679 
680 #ifdef __cpp_concepts
681  template<triangle_type t, Eigen3::eigen_dense_general A, Eigen3::eigen_general B> requires
682  (not Eigen3::eigen_MatrixWrapper<A>) and (not Eigen3::eigen_ArrayWrapper<A>) and
683  writable<A&&> and (t != triangle_type::diagonal)
684 #else
685  template<triangle_type t, typename A, typename B, std::enable_if_t<
686  Eigen3::eigen_dense_general<A> and Eigen3::eigen_general<B> and
687  (not Eigen3::eigen_MatrixWrapper<A>) and (not Eigen3::eigen_ArrayWrapper<A>) and
688  writable<A&&> and (t != triangle_type::diagonal), int> = 0>
689 #endif
690  static void
691  set_triangle(A&& a, B&& b)
692  {
693  a.template triangularView<t == triangle_type::upper ? Eigen::Upper : Eigen::Lower>() = std::forward<B>(b);
694  }
695 
696 
697 #ifdef __cpp_concepts
698  template<Eigen3::eigen_dense_general Arg> requires square_shaped<Arg> or dimension_size_of_index_is<Arg, 0, 1> or
699  std::is_lvalue_reference_v<Arg> or (not std::is_lvalue_reference_v<typename std::decay_t<Arg>::Nested>)
700 #else
701  template<typename Arg, std::enable_if_t<Eigen3::eigen_dense_general<Arg> and
702  (square_shaped<Arg> or dimension_size_of_index_is<Arg, 0, 1> or std::is_lvalue_reference_v<Arg> or
703  not std::is_lvalue_reference_v<typename std::decay_t<Arg>::Nested>), int> = 0>
704 #endif
705  static constexpr auto
706  to_diagonal(Arg&& arg)
707  {
708  if constexpr (not vector<Arg>) if (not is_vector(arg)) throw std::invalid_argument {
709  "Argument of to_diagonal must have 1 column; instead it has " + std::to_string(get_index_dimension_of<1>(arg))};
710 
711  if constexpr (square_shaped<Arg> or dimension_size_of_index_is<Arg, 0, 1>)
712  {
713  return internal::make_fixed_size_adapter(std::forward<Arg>(arg)); // Make one-dimensional matrix
714  }
715  else if constexpr (Eigen3::eigen_array_general<Arg>)
716  {
717  return arg.matrix().asDiagonal();
718  }
719  else
720  {
721  return arg.asDiagonal();
722  }
723  }
724 
725 
726 #ifdef __cpp_concepts
727  template<Eigen3::eigen_SelfAdjointView Arg> requires (hermitian_adapter_type_of_v<Arg> == HermitianAdapterType::lower)
728 #else
729  template<typename Arg, std::enable_if_t<Eigen3::eigen_SelfAdjointView<Arg>, int> = 0>
730 #endif
731  static constexpr auto
732  to_diagonal(Arg&& arg)
733  {
734  // If it is a column vector, a "lower SelfAdjointView wrapper is irrelevant
735  return OpenKalman::to_diagonal(nested_object(std::forward<Arg>(arg)));
736  }
737 
738 
739 #ifdef __cpp_concepts
740  template<Eigen3::eigen_TriangularView Arg> requires (triangle_type_of_v<Arg> == triangle_type::lower)
741 #else
742  template<typename Arg, std::enable_if_t<Eigen3::eigen_TriangularView<Arg> and (triangle_type_of_v<Arg> == triangle_type::lower), int> = 0>
743 #endif
744  static constexpr auto
745  to_diagonal(Arg&& arg)
746  {
747  // If it is a column vector, a "lower" TriangularView wrapper is irrelevant
748  return OpenKalman::to_diagonal(nested_object(std::forward<Arg>(arg)));
749  }
750 
751 
752  template<typename Arg>
753 #ifdef __cpp_concepts
754  static constexpr vector decltype(auto)
755 #else
756  static constexpr decltype(auto)
757 #endif
758  diagonal_of(Arg&& arg)
759  {
760  using Scalar = scalar_type_of_t<Arg>;
761 
762  if constexpr (Eigen3::eigen_DiagonalWrapper<Arg>)
763  {
764  using Diag = decltype(nested_object(std::forward<Arg>(arg))); //< must be nested_object(...) rather than .diagonal() because of const_cast
765  using EigenTraits = Eigen::internal::traits<std::decay_t<Diag>>;
766  constexpr auto rows = EigenTraits::RowsAtCompileTime;
767  constexpr auto cols = EigenTraits::ColsAtCompileTime;
768 
769  static_assert(cols != 1, "For Eigen::DiagonalWrapper<T> interface, T should never be a column vector "
770  "because diagonal_of function handles this case.");
771  if constexpr (cols == 0)
772  {
773  auto ret {nested_object(std::forward<Arg>(arg))};
774  return ret;
775  }
776  else if constexpr (rows == 1 or rows == 0)
777  {
778  auto ret {OpenKalman::transpose(nested_object(std::forward<Arg>(arg)))};
779  return ret;
780  }
781  else if constexpr (rows == Eigen::Dynamic or cols == Eigen::Dynamic)
782  {
783  decltype(auto) diag = nested_object(std::forward<Arg>(arg));
784  using M = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
785  return M {M::Map(to_dense_object(std::forward<Diag>(diag)).data(),
786  get_index_dimension_of<0>(diag) * get_index_dimension_of<1>(diag))};
787  }
788  else // rows > 1 and cols > 1
789  {
790  using M = Eigen::Matrix<Scalar, rows * cols, 1>;
791  return M {M::Map(to_dense_object(nested_object(std::forward<Arg>(arg))).data())};
792  }
793  }
794  else if constexpr (Eigen3::eigen_SelfAdjointView<Arg> or Eigen3::eigen_TriangularView<Arg>)
795  {
796  return OpenKalman::diagonal_of(nested_object(std::forward<Arg>(arg)));
797  }
798  else if constexpr (Eigen3::eigen_Identity<Arg>)
799  {
800  struct Op
801  {
802  constexpr auto operator()(std::size_t a, std::size_t b) const { return std::min(a, b); }
803  };
804  auto dim = values::operation(Op{}, get_index_dimension_of<0>(arg), get_index_dimension_of<1>(arg));
805  return make_constant<Arg, Scalar, 1>(dim);
806  }
807  else if constexpr (Eigen3::eigen_dense_general<Arg>)
808  {
809  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
810  using M = decltype(m);
811  using M_noref = std::remove_reference_t<M>;
812  using MatrixType = std::conditional_t<std::is_lvalue_reference_v<M>, M_noref, const M_noref>;
813  return Eigen::Diagonal<MatrixType, 0> {std::forward<M>(m)};
814  }
815  else
816  {
817  return OpenKalman::diagonal_of(to_native_matrix(std::forward<Arg>(arg)));
818  }
819  }
820 
821 
822  template<typename Arg, typename Factor0 = std::integral_constant<std::size_t, 1>, typename Factor1 = std::integral_constant<std::size_t, 1>>
823  static auto
824  broadcast(Arg&& arg, const Factor0& factor0 = Factor0{}, const Factor1& factor1 = Factor1{})
825  {
826  constexpr int F0 = []{
827  if constexpr (values::fixed<Factor0>) return static_cast<std::size_t>(Factor0{});
828  else return Eigen::Dynamic;
829  }();
830  constexpr int F1 = []{
831  if constexpr (values::fixed<Factor1>) return static_cast<std::size_t>(Factor1{});
832  else return Eigen::Dynamic;
833  }();
834 
835  using IndexType = typename std::decay_t<Arg>::Index;
836 
837  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
838  using M = decltype(m);
839  using R = Eigen::Replicate<std::remove_reference_t<M>, F0, F1>;
840  if constexpr (values::fixed<Factor0> and values::fixed<Factor1>)
841  return R {std::forward<M>(m)};
842  else
843  return R {std::forward<M>(m), static_cast<IndexType>(factor0), static_cast<IndexType>(factor1)};
844  }
845 
846  private:
847 
848  // Only to be used in a non-evaluated context
849  template<typename Op, typename...S>
850  static constexpr auto dummy_op(Op op, S...s)
851  {
852  if constexpr (std::is_invocable_v<Op, S...>) return op(s...);
853  else if constexpr (std::is_invocable_v<Op, std::size_t, std::size_t>) return op(std::size_t{0}, std::size_t{0});
854  else return op(std::size_t{0});
855  }
856 
857 
858  template<typename Op, typename...Args>
859  struct EigenNaryOp;
860 
861  template<typename Op, typename Arg>
862  struct EigenNaryOp<Op, Arg> { using type = Eigen::CwiseUnaryOp<Op, Arg>; };
863 
864  template<typename Op, typename Arg1, typename Arg2>
865  struct EigenNaryOp<Op, Arg1, Arg2> { using type = Eigen::CwiseBinaryOp<Op, Arg1, Arg2>; };
866 
867  template<typename Op, typename Arg1, typename Arg2, typename Arg3>
868  struct EigenNaryOp<Op, Arg1, Arg2, Arg3> { using type = Eigen::CwiseTernaryOp<Op, Arg1, Arg2, Arg3>; };
869 
870  public:
871 
872 #ifdef __cpp_concepts
873  template<patterns::pattern...Ds, typename Operation, indexible...Args> requires
874  (sizeof...(Ds) <= 2) and (sizeof...(Args) <= 3) and
875  (values::number<std::invoke_result_t<Operation, scalar_type_of_t<Args>...>> or
876  (sizeof...(Args) == 0 and
877  (values::number<std::invoke_result_t<Operation, std::conditional_t<true, std::size_t, Ds>...>> or
878  values::number<std::invoke_result_t<Operation, std::size_t>>)))
879 #else
880  template<typename...Ds, typename Operation, typename...Args, std::enable_if_t<sizeof...(Ds) <= 2 and sizeof...(Args) <= 3 and
881  (patterns::pattern<Ds> and ...) and (indexible<Args> and ...) and
882  (values::number<typename std::invoke_result<Operation, typename scalar_type_of<Args>::type...>::type> or
883  (sizeof...(Args) == 0 and
884  (values::number<typename std::invoke_result<Operation, std::conditional_t<true, std::size_t, Ds>...>::type> or
885  values::number<typename std::invoke_result<Operation, std::size_t>::type>))), int> = 0>
886 #endif
887  static auto
888  n_ary_operation(const std::tuple<Ds...>& tup, Operation&& operation, Args&&...args)
889  {
890  decltype(auto) op = Eigen3::native_operation(std::forward<Operation>(operation));
891  using Op = decltype(op);
892  using Scalar = decltype(dummy_op(operation, std::declval<scalar_type_of_t<Args>>()...));
893 
894  if constexpr (sizeof...(Args) == 0)
895  {
896  using P = dense_writable_matrix_t<T, data_layout::none, Scalar, std::tuple<Ds...>>;
897  return Eigen::CwiseNullaryOp<std::remove_reference_t<Op>, P> {
898  static_cast<typename P::Index>(get_dimension(std::get<0>(tup))),
899  static_cast<typename P::Index>(get_dimension(std::get<1>(tup))),
900  std::forward<Op>(op)};
901  }
902  else
903  {
904  auto seq = std::index_sequence_for<Ds...>{};
905  using CW = typename EigenNaryOp<std::decay_t<Op>, std::remove_reference_t<Args>...>::type;
906  return CW {std::forward<Args>(args)..., std::forward<Op>(op)};
907  }
908  }
909 
910 
911 #ifdef __cpp_concepts
912  template<std::size_t...indices, typename BinaryFunction, typename Arg> requires ((indices <= 1) and ...)
913 #else
914  template<std::size_t...indices, typename BinaryFunction, typename Arg, std::enable_if_t<((indices <= 1) and ...), int> = 0>
915 #endif
916  static auto
917  reduce(BinaryFunction&& b, Arg&& arg)
918  {
919  if constexpr (Eigen3::eigen_dense_general<Arg>)
920  {
921  auto&& op = Eigen3::native_operation(std::forward<BinaryFunction>(b));
922  using Op = decltype(op);
923 
924  if constexpr (((indices == 0) or ...) and ((indices == 1) or ...))
925  {
926  return std::forward<Arg>(arg).redux(std::forward<Op>(op)); // Complete reduction, which will be a scalar.
927  }
928  else
929  {
930  constexpr auto dir = ((indices == 0) and ...) ? Eigen::Vertical : Eigen::Horizontal;
931  using ROp = Eigen::internal::member_redux<std::decay_t<Op>, scalar_type_of_t<Arg>>;
932  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
933  using M = decltype(m);
934  using M_noref = std::remove_reference_t<M>;
935  using MatrixType = std::conditional_t<std::is_lvalue_reference_v<M>, M_noref, const M_noref>;
936  return Eigen::PartialReduxExpr<MatrixType, ROp, dir> {std::forward<M>(m), ROp{std::forward<Op>(op)}};
937  }
938  }
939  else
940  {
941  return reduce<indices...>(std::forward<BinaryFunction>(b), to_native_matrix(std::forward<Arg>(arg)));
942  }
943  }
944 
945 
946  template<typename Arg>
947  static constexpr decltype(auto)
948  conjugate(Arg&& arg)
949  {
950  if constexpr (Eigen3::eigen_dense_general<Arg>)
951  {
952  using UnaryOp = Eigen::internal::scalar_conjugate_op<scalar_type_of_t<Arg>>;
953  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
954  using M = decltype(m);
955  return Eigen::CwiseUnaryOp<UnaryOp, std::remove_reference_t<M>> {std::forward<M>(m)};
956  }
957  else if constexpr (Eigen3::eigen_TriangularView<Arg> or Eigen3::eigen_SelfAdjointView<Arg>)
958  {
959  return std::forward<Arg>(arg).conjugate();
960  }
961  else if constexpr (triangular_matrix<Arg>)
962  {
963  return OpenKalman::conjugate(TriangularAdapter {std::forward<Arg>(arg)});
964  }
965  else
966  {
967  return conjugate(to_native_matrix(std::forward<Arg>(arg)));
968  }
969  // Note: the global conjugate function already handles Eigen::DiagonalMatrix and Eigen::DiagonalWrapper
970  }
971 
972 
973  template<typename Arg>
974  static constexpr decltype(auto)
975  transpose(Arg&& arg)
976  {
977  if constexpr (Eigen3::eigen_matrix_general<Arg>)
978  {
979  decltype(auto) m = to_native_matrix(std::forward<Arg>(arg));
980  using M = decltype(m);
981  using M_noref = std::remove_reference_t<M>;
982  using MatrixType = std::conditional_t<std::is_lvalue_reference_v<M>, M_noref, const M_noref>;
983  return Eigen::Transpose<MatrixType> {std::forward<M>(m)};
984  }
985  else if constexpr (Eigen3::eigen_TriangularView<Arg> or Eigen3::eigen_SelfAdjointView<Arg>)
986  {
987  return std::forward<Arg>(arg).transpose();
988  }
989  else if constexpr (triangular_matrix<Arg>)
990  {
991  static_assert(triangle_type_of_v<Arg> == triangle_type::upper or triangle_type_of_v<Arg> == triangle_type::lower);
992  constexpr auto t = triangular_matrix<Arg, triangle_type::upper> ? triangle_type::lower : triangle_type::upper;
993  return OpenKalman::make_triangular_matrix<t>(transpose(to_native_matrix(std::forward<Arg>(arg))));
994  }
995  else
996  {
997  return transpose(to_native_matrix(std::forward<Arg>(arg)));
998  }
999  // Note: the global transpose function already handles zero, constant, diagonal, constant-diagonal, and symmetric cases.
1000  }
1001 
1002 
1003 #ifdef __cpp_concepts
1004  template<typename Arg> requires (not Eigen3::eigen_dense_general<Arg>)
1005 #else
1006  template<typename Arg, std::enable_if_t<not Eigen3::eigen_dense_general<Arg>, int> = 0>
1007 #endif
1008  static constexpr decltype(auto)
1009  conjugate_transpose(Arg&& arg)
1010  {
1011  if constexpr (Eigen3::eigen_TriangularView<Arg> or Eigen3::eigen_SelfAdjointView<Arg>)
1012  {
1013  return std::forward<Arg>(arg).conjugate_transpose();
1014  }
1015  else if constexpr (triangular_matrix<Arg>)
1016  {
1017  static_assert(triangle_type_of_v<Arg> == triangle_type::upper or triangle_type_of_v<Arg> == triangle_type::lower);
1018  constexpr auto t = triangular_matrix<Arg, triangle_type::upper> ? triangle_type::lower : triangle_type::upper;
1019  return OpenKalman::make_triangular_matrix<t>(OpenKalman::conjugate_transpose(to_native_matrix(std::forward<Arg>(arg))));
1020  }
1021  else
1022  {
1023  return OpenKalman::conjugate_transpose(to_native_matrix(std::forward<Arg>(arg)));
1024  }
1025  // Note: the global adjoint function already handles zero, constant, diagonal, non-complex, and hermitian cases.
1026  }
1027 
1028 
1029  template<typename Arg>
1030  static constexpr auto
1031  determinant(Arg&& arg)
1032  {
1033  if constexpr (Eigen3::eigen_matrix_general<Arg, true>)
1034  return std::forward<Arg>(arg).determinant();
1035  else if constexpr (Eigen3::eigen_array_general<Arg, true>)
1036  return std::forward<Arg>(arg).matrix().determinant();
1037  else
1038  return to_native_matrix(std::forward<Arg>(arg)).determinant();
1039  // Note: the global determinant function already handles Eigen::TriangularView, Eigen::DiagonalMatrix, and Eigen::DiagonalWrapper
1040  }
1041 
1042 
1043 #ifdef __cpp_concepts
1044  template<typename A, Eigen3::eigen_general B>
1045 #else
1046  template<typename A, typename B, std::enable_if_t<Eigen3::eigen_general<B>, int> = 0>
1047 #endif
1048  static constexpr auto
1049  sum(A&& a, B&& b)
1050  {
1051  constexpr HermitianAdapterType h = hermitian_adapter_type_of_v<A, B> == HermitianAdapterType::any ?
1052  HermitianAdapterType::lower : hermitian_adapter_type_of_v<A, B>;
1053 
1054  auto f = [](auto&& x) -> decltype(auto) {
1055  constexpr bool herm = hermitian_matrix<A> and hermitian_matrix<B>;
1056  if constexpr ((triangle_type_of_v<A, B> != triangle_type::any and triangular_adapter<decltype(x)>) or (herm and hermitian_adapter<decltype(x), h>))
1057  return nested_object(std::forward<decltype(x)>(x));
1058  else if constexpr (herm and hermitian_adapter<decltype(x)>)
1059  return OpenKalman::transpose(nested_object(std::forward<decltype(x)>(x)));
1060  else
1061  return std::forward<decltype(x)>(x);
1062  };
1063 
1064  auto op = Eigen::internal::scalar_sum_op<scalar_type_of_t<A>, scalar_type_of_t<B>>{};
1065 
1066  decltype(auto) a_wrap = to_native_matrix(f(std::forward<A>(a)));
1067  using AWrap = decltype(a_wrap);
1068  decltype(auto) b_wrap = to_native_matrix(f(std::forward<B>(b)));
1069  using BWrap = decltype(b_wrap);
1070  using CW = Eigen::CwiseBinaryOp<decltype(op), std::remove_reference_t<AWrap>, std::remove_reference_t<BWrap>>;
1071  CW s {std::forward<AWrap>(a_wrap), std::forward<BWrap>(b_wrap), std::move(op)};
1072 
1073  if constexpr (triangle_type_of_v<A, B> != triangle_type::any) return OpenKalman::make_triangular_matrix<triangle_type_of_v<A, B>>(std::move(s));
1074  else if constexpr (hermitian_matrix<A> and hermitian_matrix<B>) return make_hermitian_matrix<h>(std::move(s));
1075  else return s;
1076  }
1077 
1078  private:
1079 
1080  template<typename Arg, typename S, typename Op>
1081  static constexpr auto
1082  scalar_op_impl(Arg&& arg, S&& s, Op&& op)
1083  {
1084  using Scalar = scalar_type_of_t<Arg>;
1085  using ConstOp = Eigen::internal::scalar_constant_op<Scalar>;
1086  decltype(auto) m {to_native_matrix(std::forward<Arg>(arg))};
1087  using M = decltype(m);
1088  using PlainObjectType = dense_writable_matrix_t<M>;
1089  using Index = typename PlainObjectType::Index;
1090  auto c {Eigen::CwiseNullaryOp<ConstOp, PlainObjectType> {
1091  static_cast<Index>(get_index_dimension_of<0>(m)),
1092  static_cast<Index>(get_index_dimension_of<1>(m)),
1093  ConstOp {static_cast<Scalar>(values::to_value_type(s))}}};
1094  using CW = Eigen::CwiseBinaryOp<Op, std::remove_reference_t<M>, decltype(c)>;
1095  return CW {std::forward<M>(m), c, std::forward<Op>(op)};
1096  }
1097 
1098  public:
1099 
1100 #ifdef __cpp_concepts
1101  template<indexible Arg, values::scalar S>
1102  static constexpr vector_space_descriptors_may_match_with<Arg> auto
1103 #else
1104  template<typename Arg, typename S>
1105  static constexpr auto
1106 #endif
1107  scalar_product(Arg&& arg, S&& s)
1108  {
1109  using Scalar = scalar_type_of_t<Arg>;
1110  using Op = Eigen::internal::scalar_product_op<Scalar, Scalar>;
1111  return scalar_op_impl(std::forward<Arg>(arg), std::forward<S>(s), Op{});
1112  }
1113 
1114 
1115 #ifdef __cpp_concepts
1116  template<indexible Arg, values::scalar S>
1117  static constexpr vector_space_descriptors_may_match_with<Arg> auto
1118 #else
1119  template<typename Arg, typename S>
1120  static constexpr auto
1121 #endif
1122  scalar_quotient(Arg&& arg, S&& s)
1123  {
1124  using Scalar = scalar_type_of_t<Arg>;
1125  using Op = Eigen::internal::scalar_quotient_op<Scalar, Scalar>;
1126  return scalar_op_impl(std::forward<Arg>(arg), std::forward<S>(s), Op{});
1127  }
1128 
1129 
1130 #ifdef __cpp_concepts
1131  template<Eigen3::eigen_general A, Eigen3::eigen_general B>
1132 #else
1133  template<typename A, typename B, std::enable_if_t<(Eigen3::eigen_general<A> and Eigen3::eigen_general<B>), int> = 0>
1134 #endif
1135  static constexpr auto
1136  contract(A&& a, B&& b)
1137  {
1138  if constexpr (diagonal_matrix<A> and not Eigen3::eigen_DiagonalWrapper<A> and not Eigen3::eigen_DiagonalMatrix<A>)
1139  {
1140  return contract(to_native_matrix(diagonal_of(std::forward<A>(a))).asDiagonal(), std::forward<B>(b));
1141  }
1142  else if constexpr (diagonal_matrix<B> and not Eigen3::eigen_DiagonalWrapper<B> and not Eigen3::eigen_DiagonalMatrix<B>)
1143  {
1144  return contract(std::forward<A>(a), to_native_matrix(diagonal_of(std::forward<B>(b))).asDiagonal());
1145  }
1146  else if constexpr (diagonal_matrix<A> or diagonal_matrix<B>)
1147  {
1148  decltype(auto) a_wrap = to_native_matrix(std::forward<A>(a));
1149  using AWrap = decltype(a_wrap);
1150  decltype(auto) b_wrap = to_native_matrix(std::forward<B>(b));
1151  using BWrap = decltype(b_wrap);
1152  using Prod = Eigen::Product<std::remove_reference_t<AWrap>, std::remove_reference_t<BWrap>, Eigen::LazyProduct>;
1153  return Prod {std::forward<AWrap>(a_wrap), std::forward<BWrap>(b_wrap)};
1154  }
1155  else if constexpr (not Eigen3::eigen_matrix_general<B>)
1156  {
1157  return contract(std::forward<A>(a), to_native_matrix(std::forward<B>(b)));
1158  }
1159  else if constexpr (Eigen3::eigen_matrix_general<A> or Eigen3::eigen_TriangularView<A> or Eigen3::eigen_SelfAdjointView<A>)
1160  {
1161  using Prod = Eigen::Product<std::remove_reference_t<A>, std::remove_reference_t<B>, Eigen::DefaultProduct>;
1162  return to_dense_object(Prod {std::forward<A>(a), std::forward<B>(b)});
1163  }
1164  else
1165  {
1166  return contract(to_native_matrix(std::forward<A>(a)), std::forward<B>(b));
1167  }
1168  }
1169 
1170 
1171 #ifdef __cpp_concepts
1172  template<bool on_the_right, writable A, indexible B> requires Eigen3::eigen_dense_general<A> or
1173  Eigen3::eigen_DiagonalMatrix<A> or Eigen3::eigen_DiagonalWrapper<A>
1174 #else
1175  template<bool on_the_right, typename A, typename B, std::enable_if_t<writable<A> and
1176  (Eigen3::eigen_dense_general<A> or Eigen3::eigen_DiagonalMatrix<A> or (Eigen3::eigen_DiagonalWrapper<A> and diagonal_matrix<A> and internal::has_nested_vector<A, 0>)), int> = 0>
1177 #endif
1178  static A&
1179  contract_in_place(A& a, B&& b)
1180  {
1181  if constexpr (Eigen3::eigen_DiagonalMatrix<A> or Eigen3::eigen_DiagonalWrapper<A>)
1182  {
1183  static_assert(diagonal_matrix<B>);
1184  nested_object(a) = diagonal_of(a).array() * diagonal_of(std::forward<B>(b)).array();
1185  }
1186  else if constexpr (Eigen3::eigen_TriangularView<A>)
1187  {
1188  static_assert(triangular_matrix<B, triangle_type_of_v<A>>);
1189  nested_object(a) = contract(a, std::forward<B>(b));
1190  }
1191  else
1192  {
1193  auto&& ma = [](A& a) -> decltype(auto) {
1194  if constexpr (Eigen3::eigen_array_general<A>) return a.matrix();
1195  else return (a);
1196  }(a);
1197 
1198  if constexpr (on_the_right)
1199  return ma.applyOnTheRight(OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
1200  else
1201  return ma.applyOnTheLeft(OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
1202  return a;
1203  }
1204  }
1205 
1206 
1207 #ifdef __cpp_concepts
1208  template<triangle_type tri, Eigen3::eigen_SelfAdjointView A>
1209 #else
1210  template<triangle_type tri, typename A, std::enable_if_t<Eigen3::eigen_SelfAdjointView<A>, int> = 0>
1211 #endif
1212  static constexpr auto
1213  cholesky_factor(A&& a)
1214  {
1215  using NestedMatrix = std::decay_t<nested_object_of_t<A>>;
1216  using Scalar = scalar_type_of_t<A>;
1217  auto dim = *is_square_shaped(a);
1218 
1219  if constexpr (constant_matrix<NestedMatrix>)
1220  {
1221  // If nested matrix is a positive constant matrix, construct the Cholesky factor using a shortcut.
1222 
1223  auto s = constant_value {a};
1224 
1225  if (values::to_value_type(s) < Scalar(0))
1226  {
1227  // Cholesky factor elements are complex, so throw an exception.
1228  throw (std::runtime_error("cholesky_factor of constant HermitianAdapter: result is indefinite"));
1229  }
1230 
1231  if constexpr(tri == triangle_type::diagonal)
1232  {
1233  static_assert(diagonal_matrix<A>);
1234  return OpenKalman::to_diagonal(make_constant<A>(values::sqrt(s), dim, Dimensions<1>{}));
1235  }
1236  else if constexpr(tri == triangle_type::lower)
1237  {
1238  auto euc_dim = get_dimension(dim);
1239  auto col0 = make_constant<A>(values::sqrt(s), euc_dim, Dimensions<1>{});
1240  auto othercols = make_zero<A>(euc_dim, euc_dim - Dimensions<1>{});
1241  return attach_patterns(OpenKalman::make_triangular_matrix<tri>(concatenate_horizontal(col0, othercols)), dim, dim);
1242  }
1243  else
1244  {
1245  static_assert(tri == triangle_type::upper);
1246  auto euc_dim = get_dimension(dim);
1247  auto row0 = make_constant<A>(values::sqrt(s), Dimensions<1>{}, dim);
1248  auto otherrows = make_zero<A>(euc_dim - Dimensions<1>{}, euc_dim);
1249  return attach_patterns(OpenKalman::make_triangular_matrix<tri>(concatenate_vertical(row0, otherrows)), dim, dim);
1250  }
1251  }
1252  else
1253  {
1254  // For the general case, perform an LLT Cholesky decomposition.
1255  using M = dense_writable_matrix_t<A>;
1256  M b;
1257  auto LL_x = a.llt();
1258  if (LL_x.info() == Eigen::Success)
1259  {
1260  if constexpr((tri == triangle_type::upper and hermitian_adapter_type_of_v<A> == HermitianAdapterType::upper) or
1261  tri == triangle_type::lower and hermitian_adapter_type_of_v<A> == HermitianAdapterType::lower)
1262  {
1263  b = std::move(LL_x.matrixLLT());
1264  }
1265  else
1266  {
1267  constexpr unsigned int uplo = tri == triangle_type::upper ? Eigen::Upper : Eigen::Lower;
1268  b.template triangularView<uplo>() = LL_x.matrixLLT().conjugate_transpose();
1269  }
1270  }
1271  else [[unlikely]]
1272  {
1273  // If covariance is not positive definite, use the more robust LDLT decomposition.
1274  auto LDL_x = nested_object(std::forward<A>(a)).ldlt();
1275  if ((not LDL_x.isPositive() and not LDL_x.isNegative()) or LDL_x.info() != Eigen::Success) [[unlikely]]
1276  {
1277  if (LDL_x.isPositive() and LDL_x.isNegative()) // Covariance is zero, even though decomposition failed.
1278  {
1279  if constexpr(tri == triangle_type::lower)
1280  b.template triangularView<Eigen::Lower>() = make_zero(nested_object(a));
1281  else
1282  b.template triangularView<Eigen::Upper>() = make_zero(nested_object(a));
1283  }
1284  else // Covariance is indefinite, so throw an exception.
1285  {
1286  throw (std::runtime_error("cholesky_factor of HermitianAdapter: covariance is indefinite"));
1287  }
1288  }
1289  else if constexpr(tri == triangle_type::lower)
1290  {
1291  b.template triangularView<Eigen::Lower>() =
1292  LDL_x.matrixL().toDenseMatrix() * LDL_x.vectorD().cwiseSqrt().asDiagonal();
1293  }
1294  else
1295  {
1296  b.template triangularView<Eigen::Upper>() =
1297  LDL_x.vectorD().cwiseSqrt().asDiagonal() * LDL_x.matrixU().toDenseMatrix();
1298  }
1299  }
1300  return TriangularAdapter<M, tri> {std::move(b)};
1301  }
1302  }
1303 
1304 
1305  template<HermitianAdapterType significant_triangle, typename A, typename U, typename Alpha>
1306  static decltype(auto)
1307  rank_update_hermitian(A&& a, U&& u, const Alpha alpha)
1308  {
1309  if constexpr (Eigen3::eigen_matrix_general<A>)
1310  {
1311  static_assert(writable<A>);
1312  constexpr auto s = significant_triangle == HermitianAdapterType::lower ? Eigen::Lower : Eigen::Upper;
1313  a.template selfadjointView<s>().template rankUpdate(std::forward<U>(u), alpha);
1314  return std::forward<A>(a);
1315  }
1316  else
1317  {
1318  return rank_update_hermitian<significant_triangle>(to_native_matrix(std::forward<A>(a)), std::forward<U>(u), alpha);
1319  }
1320  }
1321 
1322 
1323  template<triangle_type triangle, typename A, typename U, typename Alpha>
1324  static decltype(auto)
1325  rank_update_triangular(A&& a, U&& u, const Alpha alpha)
1326  {
1327  if constexpr (Eigen3::eigen_matrix_general<A>)
1328  {
1329  static_assert(writable<A>);
1330  constexpr auto t = triangle == triangle_type::lower ? Eigen::Lower : Eigen::Upper;
1331  using Scalar = scalar_type_of_t<A>;
1332  for (std::size_t i = 0; i < get_index_dimension_of<1>(u); i++)
1333  {
1334  if (Eigen::internal::llt_inplace<Scalar, t>::rankUpdate(a, get_chip<1>(u, i), alpha) >= 0)
1335  throw (std::runtime_error("rank_update_triangular: product is not positive definite"));
1336  }
1337  return std::forward<A>(a);
1338  }
1339  else
1340  {
1341  return rank_update_triangular<triangle>(to_native_matrix(std::forward<A>(a)), std::forward<U>(u), alpha);
1342  }
1343  }
1344 
1345 
1346 #ifdef __cpp_concepts
1347  template<bool must_be_unique, bool must_be_exact, typename A, typename B> requires Eigen3::eigen_matrix_general<B>
1348 #else
1349  template<bool must_be_unique, bool must_be_exact, typename A, typename B, std::enable_if_t<Eigen3::eigen_matrix_general<B>, int> = 0>
1350 #endif
1351  static constexpr auto
1352  solve(A&& a, B&& b)
1353  {
1354  using Scalar = scalar_type_of_t<A>;
1355 
1356  constexpr std::size_t a_rows = dynamic_dimension<A, 0> ? index_dimension_of_v<B, 0> : index_dimension_of_v<A, 0>;
1357  constexpr std::size_t a_cols = index_dimension_of_v<A, 1>;
1358  constexpr std::size_t b_cols = index_dimension_of_v<B, 1>;
1359 
1360  if constexpr (Eigen3::eigen_TriangularView<A>)
1361  {
1362  auto ret {Eigen::Solve {std::forward<A>(a), std::forward<B>(b)}};
1363  if constexpr (std::is_lvalue_reference_v<A> and std::is_lvalue_reference_v<B>) return ret;
1364  else return to_dense_object(std::move(ret));
1365  }
1366  else if constexpr (Eigen3::eigen_SelfAdjointView<A>)
1367  {
1368  constexpr auto uplo = hermitian_adapter_type_of_v<A> == triangle_type::upper ? Eigen::Upper : Eigen::Lower;
1369  auto v {std::forward<A>(a).template selfadjointView<uplo>()};
1370  auto llt {v.llt()};
1371 
1373  if (llt.info() == Eigen::Success)
1374  {
1375  ret = Eigen::Solve {llt, std::forward<B>(b)};
1376  }
1377  else [[unlikely]]
1378  {
1379  // A is semidefinite. Use LDLT decomposition instead.
1380  auto ldlt {v.ldlt()};
1381  if ((not ldlt.isPositive() and not ldlt.isNegative()) or ldlt.info() != Eigen::Success)
1382  {
1383  throw (std::runtime_error("Eigen solve (hermitian case): A is indefinite"));
1384  }
1385  ret = Eigen::Solve {ldlt, std::forward<B>(b)};
1386  }
1387  return ret;
1388  }
1389  else if constexpr (Eigen3::eigen_matrix_general<A>)
1390  {
1392  if constexpr (must_be_exact or must_be_unique)
1393  {
1394  auto a_cols_rt = get_index_dimension_of<1>(a);
1395  auto qr {Eigen::ColPivHouseholderQR<Mat> {std::forward<A>(a)}};
1396 
1397  if constexpr (must_be_unique)
1398  {
1399  if (qr.rank() < a_cols_rt) throw std::runtime_error {"solve function requests a "
1400  "unique solution, but A is rank-deficient, so result X is not unique"};
1401  }
1402 
1403  auto res {to_dense_object(Eigen::Solve {std::move(qr), std::forward<B>(b)})};
1404 
1405  if constexpr (must_be_exact)
1406  {
1407  bool a_solution_exists = (a * res).isApprox(b, a_cols_rt * std::numeric_limits<scalar_type_of_t<A>>::epsilon());
1408 
1409  if (a_solution_exists) return res;
1410  else throw std::runtime_error {"solve function requests an exact solution, "
1411  "but the solution is only an approximation"};
1412  }
1413  else
1414  {
1415  return res;
1416  }
1417  }
1418  else
1419  {
1420  return to_dense_object(Eigen::Solve {Eigen::HouseholderQR<Mat> {std::forward<A>(a)}, std::forward<B>(b)});
1421  }
1422  }
1423  else
1424  {
1425  return solve<must_be_unique, must_be_exact>(to_native_matrix(std::forward<A>(a)), std::forward<B>(b));
1426  }
1427  }
1428 
1429  private:
1430 
1431  template<typename A>
1432  static constexpr auto
1433  QR_decomp_impl(A&& a)
1434  {
1435  using Scalar = scalar_type_of_t<A>;
1436  constexpr auto rows = index_dimension_of_v<A, 0>;
1437  constexpr auto cols = index_dimension_of_v<A, 1>;
1438  using MatrixType = Eigen3::eigen_matrix_t<Scalar, rows, cols>;
1439  using ResultType = Eigen3::eigen_matrix_t<Scalar, cols, cols>;
1440 
1441  Eigen::HouseholderQR<MatrixType> QR {std::forward<A>(a)};
1442 
1443  if constexpr (dynamic_dimension<A, 1>)
1444  {
1445  auto rt_cols = get_index_dimension_of<1>(a);
1446 
1447  ResultType ret {rt_cols, rt_cols};
1448 
1449  if constexpr (dynamic_dimension<A, 0>)
1450  {
1451  auto rt_rows = get_index_dimension_of<0>(a);
1452 
1453  if (rt_rows < rt_cols)
1454  ret << QR.matrixQR().topRows(rt_rows),
1456  else
1457  ret = QR.matrixQR().topRows(rt_cols);
1458  }
1459  else
1460  {
1461  if (rows < rt_cols)
1462  ret << QR.matrixQR().template topRows<rows>(),
1464  else
1465  ret = QR.matrixQR().topRows(rt_cols);
1466  }
1467 
1468  return ret;
1469  }
1470  else
1471  {
1472  ResultType ret;
1473 
1474  if constexpr (dynamic_dimension<A, 0>)
1475  {
1476  auto rt_rows = get_index_dimension_of<0>(a);
1477 
1478  if (rt_rows < cols)
1479  ret << QR.matrixQR().topRows(rt_rows),
1481  else
1482  ret = QR.matrixQR().template topRows<cols>();
1483  }
1484  else
1485  {
1486  if constexpr (rows < cols)
1487  ret << QR.matrixQR().template topRows<rows>(), Eigen3::eigen_matrix_t<Scalar, cols - rows, cols>::Zero();
1488  else
1489  ret = QR.matrixQR().template topRows<cols>();
1490  }
1491 
1492  return ret;
1493  }
1494  }
1495 
1496  public:
1497 
1498  template<typename A>
1499  static constexpr auto
1500  LQ_decomposition(A&& a)
1501  {
1502  return OpenKalman::make_triangular_matrix<triangle_type::lower>(OpenKalman::conjugate_transpose(QR_decomp_impl(OpenKalman::conjugate_transpose(std::forward<A>(a)))));
1503  }
1504 
1505 
1506  template<typename A>
1507  static constexpr auto
1508  QR_decomposition(A&& a)
1509  {
1510  return OpenKalman::make_triangular_matrix<triangle_type::upper>(QR_decomp_impl(std::forward<A>(a)));
1511  }
1512 
1513  };
1514 
1515 
1516 }
1517 
1518 #endif
decltype(auto) constexpr concatenate_vertical(V &&v, Vs &&... vs)
Concatenate one or more typed matrices objects vertically.
Definition: typed-matrix-overloads.hpp:270
constexpr auto is_square_shaped(const T &t)
At least 2 and at most N indices have the same extent.
Definition: is_square_shaped.hpp:43
constexpr auto n_ary_operation(const std::tuple< Ds... > &d_tup, Operation &&operation, Args &&...args)
Perform a component-wise n-ary operation, using broadcasting to match the size of a pattern matrix...
Definition: n_ary_operation.hpp:325
decltype(auto) constexpr contract(A &&a, B &&b)
Matrix multiplication of A * B.
Definition: contract.hpp:54
Definition: basics.hpp:41
constexpr bool hermitian_adapter
Specifies that a type is a hermitian matrix adapter of a particular type.
Definition: hermitian_adapter.hpp:55
decltype(auto) rank_update_hermitian(A &&a, U &&u, scalar_type_of_t< A > alpha=1)
Do a rank update on a hermitian matrix.
Definition: rank_update_hermitian.hpp:45
Definition: eigen-forward-declarations.hpp:90
triangle_type
The type of a triangular matrix.
Definition: enumerations.hpp:26
decltype(auto) constexpr attach_patterns(Arg &&arg, P &&p)
Attach a pattern_collection to an indexible object.
Definition: attach_patterns.hpp:74
A lower-left triangular matrix.
decltype(auto) constexpr conjugate(Arg &&arg)
Take the complex conjugate of an indexible object.
Definition: conjugate.hpp:44
decltype(auto) constexpr get_slice(Arg &&arg, const std::tuple< Offset... > &offsets, const std::tuple< Extent... > &extents)
Extract a slice from a matrix or tensor.
Definition: get_slice.hpp:103
Forward declarations for OpenKalman&#39;s Eigen interface.
A triangular_adapter, where components above or below the diagonal (or both) are zero.
Definition: forward-class-declarations.hpp:147
constexpr bool pattern
An object describing the characteristics (e.g., dimensions, wrapping structure) of an index...
Definition: pattern.hpp:31
constexpr bool indexible
T is a multidimensional array type.
Definition: indexible.hpp:32
decltype(auto) constexpr QR_decomposition(A &&a)
Perform a QR decomposition of matrix A=Q[U,0], U is a upper-triangular matrix, and Q is orthogonal...
Definition: QR_decomposition.hpp:33
decltype(auto) constexpr to_value_type(Arg &&arg)
Convert, if necessary, a fixed or dynamic value to its underlying base type.
Definition: to_value_type.hpp:28
constexpr bool number
T is a numerical field type.
Definition: number.hpp:41
decltype(auto) constexpr conjugate_transpose(Arg &&arg)
Take the conjugate-transpose of an indexible_object.
Definition: conjugate_transpose.hpp:35
HermitianAdapterType
The type of a hermitian adapter, indicating which triangle of the nested matrix is used...
Definition: enumerations.hpp:79
constexpr bool triangular_adapter
Specifies that a type is a triangular adapter of triangle type tri.
Definition: triangular_adapter.hpp:47
constexpr bool value
T is a fixed or dynamic value that is reducible to a number.
Definition: value.hpp:45
decltype(auto) constexpr concatenate_horizontal(V &&v, Vs &&... vs)
Concatenate one or more matrix objects vertically.
Definition: typed-matrix-overloads.hpp:308
A diagonal matrix (both a lower-left and an upper-right triangular matrix).
decltype(auto) constexpr to_dense_object(Arg &&arg)
Convert the argument to a dense, writable matrix of a particular scalar type.
Definition: to_dense_object.hpp:37
decltype(auto) constexpr reduce(BinaryFunction &&b, Arg &&arg)
Perform a partial reduction based on an associative binary function, across one or more indices...
Definition: reduce.hpp:143
A structure representing the dimensions associated with of a particular index.
Definition: Dimensions.hpp:42
std::conditional_t< sizeof...(dims)==1, Eigen::Matrix< Scalar, detail::eigen_index_convert< dims >..., detail::eigen_index_convert< 1 > >, Eigen::Matrix< Scalar, detail::eigen_index_convert< dims >... > > eigen_matrix_t
An alias for a self-contained, writable, native Eigen matrix.
Definition: eigen-forward-declarations.hpp:491
decltype(auto) constexpr broadcast(Arg &&arg, const Factors &...factors)
Broadcast an object by replicating it by factors specified for each index.
Definition: broadcast.hpp:49
decltype(auto) constexpr apply(F &&f, T &&t)
A generalization of std::apply.
Definition: apply.hpp:49
constexpr bool triangular_matrix
Specifies that an argument is an indexible object having a given triangle_type (e.g., upper, lower, or diagonal).
Definition: triangular_matrix.hpp:36
Lower, upper, or diagonal matrix.
constexpr auto make_identity_matrix(P &&p)
Make an identity_matrix with a given shape pattern.
Definition: make_identity_matrix.hpp:37
An interface to various routines from the linear algebra library associated with indexible object T...
Definition: library_interface.hpp:42
Definition: object_traits.hpp:38
constexpr bool identity_matrix
Specifies that a type is known at compile time to be a rank-2 or lower identity matrix.
Definition: identity_matrix.hpp:50
decltype(auto) constexpr transpose(Arg &&arg)
Swap any two indices of an indexible_object.
Definition: transpose.hpp:49
constexpr A && contract_in_place(A &&a, B &&b)
In-place matrix multiplication of A * B, storing the result in A.
Definition: contract_in_place.hpp:38
decltype(auto) constexpr diagonal_of(Arg &&arg)
Extract a column vector (or column slice for rank>2 tensors) comprising the diagonal elements...
Definition: diagonal_of.hpp:44
constexpr auto solve(A &&a, B &&b)
Solve the equation AX = B for X, which may or may not be a unique solution.
Definition: solve.hpp:87
constexpr auto conj(const Arg &arg)
A constexpr function for the complex conjugate of a (complex) number.
Definition: conj.hpp:39
constexpr auto sqrt(const Arg &arg)
A constexpr alternative to std::sqrt.
Definition: sqrt.hpp:46
constexpr auto determinant(Arg &&arg)
Take the determinant of a matrix.
Definition: determinant.hpp:44
constexpr bool size
T is either an index representing a size, or unbounded_size_t, which indicates that the size is unbou...
Definition: size.hpp:71
Arg && fill_components(Arg &&arg, S...s)
This is an overloaded member function, provided for convenience. It differs from the above function o...
Definition: fill_components.hpp:67
constexpr auto is_vector(const T &t)
Return true if T is a vector at runtime.
Definition: is_vector.hpp:60
An upper-right triangular matrix.
decltype(auto) constexpr sum(Ts &&...ts)
Element-by-element sum of one or more objects.
Definition: sum.hpp:112
constexpr auto constant_value(T &&t)
The constant value associated with a constant_object or constant_diagonal_object. ...
Definition: constant_value.hpp:37
decltype(auto) constexpr LQ_decomposition(A &&a)
Perform an LQ decomposition of matrix A=[L,0]Q, L is a lower-triangular matrix, and Q is orthogonal...
Definition: LQ_decomposition.hpp:33
decltype(auto) rank_update_triangular(A &&a, U &&u, scalar_type_of_t< A > alpha=1)
Do a rank update on triangular matrix.
Definition: rank_update_triangular.hpp:48
std::decay_t< decltype(make_dense_object< T, layout, S >(std::declval< D >()))> dense_writable_matrix_t
An alias for a dense, writable matrix, patterned on parameter T.
Definition: dense_writable_matrix_t.hpp:38
Definition: trait_backports.hpp:64
constexpr auto make_zero(P &&p)
Make an indexible object in which every element is 0.
Definition: make_zero.hpp:37
decltype(auto) constexpr to_diagonal(Arg &&arg, P &&p)
Convert a column vector (or any other array with a 1D second index) into a diagonal_matrix.
Definition: to_diagonal.hpp:46
constexpr bool vector
T is a vector (e.g., column or row vector).
Definition: vector.hpp:60
constexpr auto get_dimension(const Arg &arg)
Get the vector dimension of patterns::pattern Arg.
Definition: get_dimension.hpp:53
A matrix with typed rows and columns.
Definition: forward-class-declarations.hpp:292
decltype(auto) constexpr make_triangular_matrix(Arg &&arg)
Create a triangular_matrix from a general matrix.
Definition: make_triangular_matrix.hpp:35
constexpr auto operation(Operation &&op, Args &&...args)
A potentially constant-evaluated operation involving some number of values.
Definition: operation.hpp:98
decltype(auto) constexpr cholesky_factor(A &&a)
Take the Cholesky factor of a matrix.
Definition: cholesky_factor.hpp:38
constexpr Arg && set_slice(Arg &&arg, Block &&block, const Begin &...begin)
Assign an object to a particular slice of a matrix or tensor.
Definition: set_slice.hpp:56