OpenKalman
eigen-tensor-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) 2023 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_TENSOR_LIBRARY_INTERFACE_HPP
17 #define OPENKALMAN_EIGEN_TENSOR_LIBRARY_INTERFACE_HPP
18 
19 #include <type_traits>
20 #include <tuple>
21 #include <random>
23 
24 
25 // ------------------- //
26 // library_interface //
27 // ------------------- //
28 
29 namespace OpenKalman::interface
30 {
31 #ifdef __cpp_concepts
32  template<Eigen3::eigen_tensor_general<true> T>
33  struct library_interface<T>
34 #else
35  template<typename T>
36  struct library_interface<T, std::enable_if_t<Eigen3::eigen_tensor_general<T, true>>>
37 #endif
38  {
39  private:
40 
41  using IndexType = typename Eigen::internal::traits<T>::Index;
42 
43  public:
44 
45  template<typename Derived>
47 
48 
49 #ifdef __cpp_lib_concepts
50  template<indexible Arg, std::size_t NumIndices> requires (NumIndices == index_count_v<Arg>)
51  static constexpr values::scalar decltype(auto)
52 #else
53  template <typename Arg, std::size_t NumIndices, std::enable_if_t<NumIndices == index_count_v<Arg>, int> = 0>
54  static constexpr decltype(auto)
55 #endif
56  get_component(Arg&& arg, const std::array<IndexType, NumIndices>& indices)
57  {
58  if constexpr ((Eigen::internal::traits<T>::Flags & Eigen::LvalueBit) != 0)
59  return Eigen::TensorEvaluator<std::remove_reference_t<Arg>, Eigen::DefaultDevice> {std::forward<Arg>(arg), Eigen::DefaultDevice{}}.coeffRef(indices);
60  else
61  return Eigen::TensorEvaluator<std::remove_reference_t<Arg>, Eigen::DefaultDevice> {std::forward<Arg>(arg), Eigen::DefaultDevice{}}.coeff(indices);
62  }
63 
64  private:
65 
66  template<std::size_t ix_count, typename I, typename End, typename...Ixs>
67  static constexpr decltype(auto)
68  make_ix_array(I i, const End& end, Ixs...indices)
69  {
70  if constexpr (ix_count == 0)
71  {
72  if (i != end) throw std::logic_error("Too many indices on component access");
73  return std::array<IndexType, sizeof...(Ixs)> {indices...};
74  }
75  else
76  {
77  if (i == end) throw std::logic_error("Not enough indices on component access");
78  auto this_ix {*i};
79  return make_ix_array<ix_count - 1>(++i, end, indices..., static_cast<IndexType>(this_ix));
80  }
81  }
82 
83  public:
84 
85 #ifdef __cpp_lib_ranges
86  template<indexible Arg, std::ranges::input_range Indices> requires
87  std::convertible_to<std::ranges::range_value_t<Indices>, const typename std::decay_t<Arg>::Index>
88  static constexpr values::scalar decltype(auto)
89 #else
90  template<typename Arg, typename Indices>
91  static constexpr decltype(auto)
92 #endif
93  get_component(Arg&& arg, const Indices& indices)
94  {
95  constexpr std::size_t ix_count = Eigen::internal::traits<T>::NumDimensions;
96  return get_component(std::forward<Arg>(arg), make_ix_array<ix_count>(stdex::ranges::begin(indices), stdex::ranges::end(indices)));
97  }
98 
99 
100 #ifdef __cpp_lib_concepts
101  template<indexible Arg, std::size_t NumIndices> requires (NumIndices == index_count_v<Arg>) and (not std::is_const_v<Arg>)
102 #else
104  not std::is_const<Arg>::value, int> = 0>
105 #endif
106  static void
107  set_component(Arg& arg, const scalar_type_of_t<T>& s, const std::array<IndexType, NumIndices>& indices)
108  {
109  Eigen::TensorEvaluator<Arg, Eigen::DefaultDevice> {arg, Eigen::DefaultDevice{}}.coeffRef(indices) = s;
110  }
111 
112 
113 #ifdef __cpp_lib_ranges
114  template<indexible Arg, std::ranges::input_range Indices> requires
115  std::convertible_to<std::ranges::range_value_t<Indices>, const typename std::decay_t<Arg>::Index> and
116  ((Eigen::internal::traits<std::decay_t<Arg>>::Flags & Eigen::LvalueBit) != 0x0) and (not std::is_const_v<Arg>)
117 #else
118  template<typename Arg, typename Indices, std::enable_if_t<
119  ((Eigen::internal::traits<std::decay_t<Arg>>::Flags & Eigen::LvalueBit) != 0x0) and (not std::is_const_v<Arg>), int> = 0>
120 #endif
121  static void
122  set_component(Arg& arg, const scalar_type_of_t<T>& s, const Indices& indices)
123  {
124  constexpr std::size_t ix_count = Eigen::internal::traits<T>::NumDimensions;
125  set_component(arg, s, make_ix_array<ix_count>(stdex::ranges::begin(indices), stdex::ranges::end(indices)));
126  }
127 
128  private:
129 
130  template<auto l, typename Arg, std::size_t...Is>
131  static constexpr auto
132  make_TensorMap_impl(const Arg& arg, std::index_sequence<Is...>)
133  {
134  using M = Eigen::TensorFixedSize<scalar_type_of_t<Arg>, Eigen::Sizes<static_cast<std::ptrdiff_t>(index_dimension_of_v<Arg, Is>)...>, l, IndexType>;
135  return Eigen::TensorMap<const M, l> {internal::raw_data(arg), index_dimension_of_v<Arg, Is>...};
136  }
137 
138 
139  template<auto l, typename Arg>
140  static constexpr auto
141  make_TensorMap(const Arg& arg)
142  {
143  if constexpr (has_dynamic_dimensions<Arg>)
144  {
145  using M = Eigen::Tensor<scalar_type_of_t<Arg>, index_count_v<Arg>, l, IndexType>;
146  return Eigen::TensorMap<const M, l>{internal::raw_data(arg)};
147  }
148  else
149  {
150  return make_TensorMap_impl<l>(arg, std::make_index_sequence<index_count_v<Arg>>{});
151  }
152  }
153 
154  public:
155 
156  template<typename Arg>
157  static decltype(auto)
158  to_native_matrix(Arg&& arg)
159  {
160  if constexpr (Eigen3::eigen_tensor_wrapper<Arg>)
161  {
162  return std::forward<Arg>(arg);
163  }
164  else if constexpr (internal::library_wrapper<Arg>)
165  {
166  return to_native_matrix(OpenKalman::nested_object(std::forward<Arg>(arg)));
167  }
168  else if constexpr (not Eigen3::eigen_tensor_general<Arg> and directly_accessible<Arg> and std::is_lvalue_reference_v<Arg>)
169  {
170  if constexpr (layout_of_v<Arg> == data_layout::stride and internal::has_static_strides<Arg>)
171  {
172  auto strides = internal::strides(arg);
173  constexpr std::ptrdiff_t stride0 = std::get<0>(strides);
174  constexpr std::ptrdiff_t strideN = std::get<index_count_v<Arg> - 1>(strides);
175  constexpr auto l = stride0 > strideN ? Eigen::RowMajor : Eigen::ColMajor;
176  return std::apply(
177  [](auto&& map, auto&&...s){
178  return Eigen::TensorStridingOp {map, std::array<std::size_t, index_count_v<Arg>>{s...}};
179  },
180  std::tuple_cat(make_TensorMap<l>(arg), std::move(strides)));
181  }
182  else
183  {
184  constexpr auto l = layout_of_v<Arg> == data_layout::right ? Eigen::RowMajor : Eigen::ColMajor;
185  return make_TensorMap<l>(arg);
186  }
187  }
188  else
189  {
190  return Eigen3::make_eigen_tensor_wrapper(std::forward<Arg>(arg));
191  }
192  }
193 
194 
195 #ifdef __cpp_concepts
196  template<typename To, Eigen3::eigen_tensor_general From> requires (std::assignable_from<To&, From&&>)
197 #else
198  template<typename To, typename From, std::enable_if_t<Eigen3::eigen_tensor_general<From> and std::is_assignable_v<To&, From&&>, int> = 0>
199 #endif
200  static void copy(To& a, From&& b)
201  {
202  a = std::forward<From>(b);
203  }
204 
205 
206  template<data_layout layout, typename Scalar, typename...D>
207  static auto make_default(D&&...d)
208  {
209  constexpr auto options = layout == data_layout::right ? Eigen::RowMajor : Eigen::ColMajor;
210  if constexpr (((dynamic_pattern<D>) or ...))
211  return Eigen::Tensor<Scalar, sizeof...(D), options, IndexType>(static_cast<IndexType>(get_dimension(d))...);
212  else
213  return Eigen::TensorFixedSize<Scalar, Eigen::Sizes<static_cast<std::ptrdiff_t>(coordinates::dimension_of_v<D>)...>, options, IndexType> {};
214  }
215 
216 
217 #ifdef __cpp_concepts
218  template<data_layout layout, writable Arg, std::convertible_to<scalar_type_of_t<Arg>> ... Scalars>
219  requires (layout == data_layout::right) or (layout == data_layout::left)
220 #else
221  template<data_layout layout, typename Arg, typename...Scalars, std::enable_if_t<writable<Arg> and
222  (layout == data_layout::right or or layout == data_layout::left) and
223  std::conjunction<std::is_convertible<Scalars, typename scalar_type_of<Arg>::type>::value...>::value, int> = 0>
224 #endif
225  static void fill_components(Arg& arg, const Scalars ... scalars)
226  {
227  if constexpr (layout == data_layout::left)
228  arg.swap_layout().setValues({scalars...});
229  else
230  arg.setValues({scalars...});
231  }
232 
233 
234 #ifdef __cpp_concepts
235  template<values::dynamic C, typename...Ds> requires (... and (not dynamic_pattern<Ds>))
236  static constexpr constant_matrix auto
237 #else
238  template<typename C, typename...Ds, std::enable_if_t<values::dynamic<C> and
239  (... and (not dynamic_pattern<Ds>)), int> = 0>
240  static constexpr auto
241 #endif
242  make_constant(C&& c, Ds&&...ds)
243  {
244  auto value = values::to_value_type(std::forward<C>(c));
245  using Scalar = std::decay_t<decltype(value)>;
246  auto m = make_default<data_layout::none, Scalar>(std::forward<Ds>(ds)...);
247  // m will be a dangling reference to TensorFixedSize, but Eigen only references its static dimensions, so there is no bug
248  return Eigen::TensorCwiseNullaryOp {m, Eigen::internal::scalar_constant_op<Scalar>(value)};
249  }
250 
251 
252  // make_identity_matrix not defined
253 
254 
255  template<typename Arg, typename...Begin, typename...Size>
256  static auto
257  get_slice(Arg&& arg, std::tuple<Begin...> begin, std::tuple<Size...> size)
258  {
259  auto offsets = std::apply([](auto&&...a) {
260  return std::array<std::size_t, sizeof...(Begin)> {std::forward<decltype(a)>(a)...};
261  }, begin);
262  auto extents = std::apply([](auto&&...a) {
263  return std::array<std::size_t, sizeof...(Size)> {std::forward<decltype(a)>(a)...};
264  }, size);
265  return to_native_matrix(std::forward<Arg>(arg)).slice(offsets, extents);
266  }
267 
268 
269 #ifdef __cpp_concepts
270  template<Eigen3::eigen_tensor_general<true> Arg, Eigen3::eigen_tensor_general Block, typename...Begin>
271 #else
272  template<typename Arg, typename Block, typename...Begin, std::enable_if_t<
273  Eigen3::eigen_tensor_general<Arg, true> and Eigen3::eigen_tensor_general<Block>, int> = 0>
274 #endif
275  static constexpr void
276  set_slice(Arg& arg, Block&& block, Begin...begin)
277  {
278  auto offsets = std::array {static_cast<std::size_t>(begin)...};
279  auto extents = std::apply([](auto&&...a) {
280  return std::array {get_dimension(std::forward<decltype(a)>(a))...}; }, all_vector_space_descriptors(block));
281  arg.slice(offsets, extents) = std::forward<Block>(block);
282  }
283 
284  /*
285  private:
286 
287 #ifdef __cpp_concepts
288  template<typename A>
289 #else
290  template<typename A, typename = void>
291 #endif
292  struct pass_through_eigenwrapper : std::false_type {};
293 
294 #ifdef __cpp_concepts
295  template<Eigen3::eigen_wrapper A>
296  struct pass_through_eigenwrapper<A>
297 #else
298  template<typename A>
299  struct pass_through_eigenwrapper<A, std::enable_if_t<Eigen3::eigen_wrapper<A>>>
300 #endif
301  : std::bool_constant<Eigen3::eigen_dense_general<nested_object_of_t<A>> or (diagonal_matrix<nested_object_of_t<A>> and internal::has_nested_vector<nested_object_of_t<A>>) or
302  triangular_adapter<nested_object_of_t<A>> or hermitian_adapter<nested_object_of_t<A>>> {};
303 
304  public:
305 
306  template<triangle_type t, typename A, typename B>
307  static decltype(auto) set_triangle(A&& a, B&& b)
308  {
309  if constexpr (Eigen3::eigen_MatrixWrapper<A> or Eigen3::eigen_ArrayWrapper<A> or pass_through_eigenwrapper<A>::value)
310  {
311  return internal::set_triangle<t>(nested_object(std::forward<A>(a)), std::forward<B>(b));
312  }
313  else if constexpr (not Eigen3::eigen_dense_general<A>)
314  {
315  return set_triangle<t>(Eigen3::make_eigen_wrapper(std::forward<A>(a)), std::forward<B>(b));
316  }
317  else
318  {
319  if constexpr (t == triangle_type::diagonal)
320  {
321  if constexpr (writable<A> and std::is_lvalue_reference_v<A>)
322  {
323  a.diagonal() = OpenKalman::diagonal_of(std::forward<B>(b));
324  return std::forward<A>(a);
325  }
326  else
327  {
328  auto aw = to_dense_object(std::forward<A>(a));
329  aw.diagonal() = OpenKalman::diagonal_of(std::forward<B>(b));
330  return aw;
331  }
332  }
333  else
334  {
335  decltype(auto) aw = [](A&& a) -> decltype(auto) {
336  if constexpr (writable<A> and std::is_lvalue_reference_v<A>) return std::forward<A>(a);
337  else return to_dense_object(std::forward<A>(a));
338  }(std::forward<A>(a));
339 
340  auto tview = aw.template triangularView<t == triangle_type::upper ? Eigen::Upper : Eigen::Lower>();
341  if constexpr (std::is_assignable_v<decltype(tview), B&&>)
342  tview = std::forward<B>(b);
343  else
344  tview = Eigen3::make_eigen_wrapper(std::forward<B>(b));
345  return std::forward<decltype(aw)>(aw);
346  }
347  }
348  }
349 
350 
351 #ifdef __cpp_concepts
352  template<Eigen3::eigen_dense_general Arg> requires
353  std::is_lvalue_reference_v<nested_object_of_t<Eigen::DiagonalWrapper<std::remove_reference_t<Arg>>>>
354 #else
355  template<typename Arg, std::enable_if_t<Eigen3::eigen_dense_general<Arg> and
356  std::is_lvalue_reference_v<typename nested_object_of<Eigen::DiagonalWrapper<std::remove_reference_t<Arg>>>::type>, int> = 0>
357 #endif
358  static constexpr auto
359  to_diagonal(Arg& arg)
360  {
361  if constexpr (not vector<Arg>) if (not is_vector(arg)) throw std::invalid_argument {
362  "Argument of to_diagonal must have 1 column; instead it has " + std::to_string(get_index_dimension_of<1>(arg))};
363  return Eigen::DiagonalWrapper<std::remove_reference_t<Arg>> {arg};
364  }
365 
366 
367 #ifdef __cpp_concepts
368  template<Eigen3::eigen_SelfAdjointView Arg>
369 #else
370  template<typename Arg, std::enable_if_t<Eigen3::eigen_SelfAdjointView<Arg>, int> = 0>
371 #endif
372  static constexpr auto
373  to_diagonal(Arg&& arg)
374  {
375  // If it is a column vector, the SelfAdjointView wrapper doesn't matter, and otherwise, the following will throw an exception:
376  return OpenKalman::to_diagonal(nested_object(std::forward<Arg>(arg)));
377  }
378 
379 
380 #ifdef __cpp_concepts
381  template<Eigen3::eigen_TriangularView Arg>
382 #else
383  template<typename Arg, std::enable_if_t<Eigen3::eigen_TriangularView<Arg>, int> = 0>
384 #endif
385  static constexpr auto
386  to_diagonal(Arg&& arg)
387  {
388  // If it is a column vector, the TriangularView wrapper doesn't matter, and otherwise, the following will thow an exception:
389  return OpenKalman::to_diagonal(nested_object(std::forward<Arg>(arg)));
390  }
391 
392  private:
393 
394  template<typename Arg, typename Id>
395  static constexpr decltype(auto)
396  diagonal_of_impl(Arg&& arg, Id&& id)
397  {
398  auto diag {[](Arg&& arg){
399  if constexpr (Eigen3::eigen_array_general<Arg, true>)
400  return make_self_contained<Arg>(std::forward<Arg>(arg).matrix().diagonal());
401  else // eigen_matrix_general<Arg, true>
402  return make_self_contained<Arg>(std::forward<Arg>(arg).diagonal());
403  }(std::forward<Arg>(arg))};
404  using Diag = const std::decay_t<decltype(diag)>;
405  static_assert(vector<Diag>);
406 
407  using D = std::conditional_t<dynamic_dimension<Diag, 0> and fixed_pattern<Id>, std::decay_t<Id>, vector_space_descriptor_of_t<Diag, 0>>;
408 
409  if constexpr (not dynamic_dimension<Diag, 0> or dynamic_pattern<D>) return diag;
410  else return internal::FixedSizeAdapter<Diag, D, Dimensions<1>> {std::move(diag)};
411  }
412 
413  public:
414 
415 #ifdef __cpp_concepts
416  template<Eigen3::eigen_general<true> Arg>
417 #else
418  template<typename Arg, std::enable_if_t<Eigen3::eigen_general<Arg, true>, int> = 0>
419 #endif
420  static constexpr decltype(auto)
421  diagonal_of(Arg&& arg)
422  {
423  auto d = is_square_shaped(arg);
424  if (not d) throw std::invalid_argument {"Argument of diagonal_of must be a square matrix; instead it has " +
425  std::to_string(get_index_dimension_of<0>(arg)) + " rows and " +
426  std::to_string(get_index_dimension_of<1>(arg)) + " columns"};
427 
428  using Scalar = scalar_type_of_t<Arg>;
429 
430  if constexpr (Eigen3::eigen_DiagonalWrapper<Arg>)
431  {
432  using Scalar = scalar_type_of_t<Arg>;
433  decltype(auto) diag {nested_object(std::forward<Arg>(arg))}; //< must be nested_object(...) rather than .diagonal() because of const_cast
434  using Diag = decltype(diag);
435  using EigenTraits = Eigen::internal::traits<std::decay_t<Diag>>;
436  constexpr auto rows = EigenTraits::RowsAtCompileTime;
437  constexpr auto cols = EigenTraits::ColsAtCompileTime;
438 
439  static_assert(cols != 1, "For Eigen::DiagonalWrapper<T> interface, T should never be a column vector "
440  "because diagonal_of function handles this case.");
441  if constexpr (cols == 0)
442  {
443  return std::forward<Diag>(diag);
444  }
445  else if constexpr (rows == 1 or rows == 0)
446  {
447  return OpenKalman::transpose(std::forward<Diag>(diag));
448  }
449  else if constexpr (rows == Eigen::Dynamic or cols == Eigen::Dynamic)
450  {
451  auto d {to_dense_object(std::forward<Diag>(diag))};
452  using M = Eigen::Matrix<Scalar, Eigen::Dynamic, 1>;
453  return M {M::Map(to_dense_object(std::forward<Diag>(diag)).data(),
454  get_index_dimension_of<0>(diag) * get_index_dimension_of<1>(diag))};
455  }
456  else // rows > 1 and cols > 1
457  {
458  using M = Eigen::Matrix<Scalar, rows * cols, 1>;
459  return M {M::Map(to_dense_object(std::forward<Diag>(diag)).data())};
460  }
461  }
462  else if constexpr (Eigen3::eigen_SelfAdjointView<Arg> or Eigen3::eigen_TriangularView<Arg>)
463  {
464  // Assume there are no dangling references
465  return OpenKalman::diagonal_of(nested_object(std::forward<Arg>(arg)));
466  }
467  else if constexpr (Eigen3::eigen_wrapper<Arg>)
468  {
469  if constexpr (Eigen3::eigen_general<nested_object_of_t<Arg>, true>)
470  return OpenKalman::diagonal_of(nested_object(std::forward<Arg>(arg)));
471  else
472  return diagonal_of_impl(std::forward<Arg>(arg), *d);
473  }
474  else if constexpr (Eigen3::eigen_Identity<Arg>)
475  {
476  constexpr std::size_t dim = dynamic_dimension<Arg, 0> ? index_dimension_of_v<Arg, 1> : index_dimension_of_v<Arg, 0>;
477  if constexpr (dim == stdex::dynamic_extent) return make_constant<Arg, Scalar, 1>(*d, Dimensions<1>{});
478  else return make_constant<Arg, Scalar, 1>(Dimensions<dim>{}, Dimensions<1>{});
479  }
480  else
481  {
482  return diagonal_of_impl(std::forward<Arg>(arg), *d);
483  };
484  }
485 
486  private:
487 
488  template<typename...Ds, typename...ArgDs, typename Arg, std::size_t...I>
489  static decltype(auto)
490  replicate_arg_impl(const std::tuple<Ds...>& p_tup, const std::tuple<ArgDs...>& arg_tup, Arg&& arg, std::index_sequence<I...>)
491  {
492  using R = Eigen::Replicate<std::decay_t<Arg>,
493  (coordinates::dimension_of_v<Ds> == stdex::dynamic_extent or coordinates::dimension_of_v<ArgDs> == stdex::dynamic_extent ?
494  Eigen::Dynamic : static_cast<IndexType>(coordinates::dimension_of_v<Ds> / coordinates::dimension_of_v<ArgDs>))...>;
495 
496  if constexpr (((coordinates::dimension_of_v<Ds> != stdex::dynamic_extent) and ...) and
497  ((coordinates::dimension_of_v<ArgDs> != stdex::dynamic_extent) and ...))
498  {
499  if constexpr (((coordinates::dimension_of_v<Ds> == coordinates::dimension_of_v<ArgDs>) and ...))
500  return std::forward<Arg>(arg);
501  else
502  return R {std::forward<Arg>(arg)};
503  }
504  else
505  {
506  auto ret = R {std::forward<Arg>(arg), static_cast<IndexType>(
507  get_dimension(std::get<I>(p_tup)) / get_dimension(std::get<I>(arg_tup)))...};
508  return ret;
509  }
510  }
511 
512 
513  template<typename...Ds, typename Arg>
514  static decltype(auto)
515  replicate_arg(const std::tuple<Ds...>& p_tup, Arg&& arg)
516  {
517  return replicate_arg_impl(p_tup, all_vector_space_descriptors(arg), std::forward<Arg>(arg), std::index_sequence_for<Ds...> {});
518  }
519 
520  public:
521 
522 #ifdef __cpp_concepts
523  template<typename...Ds, typename Operation, typename...Args> requires (sizeof...(Args) <= 3) and
524  std::invocable<Operation&&, scalar_type_of_t<Args>...> and
525  values::number<std::invoke_result_t<Operation&&, scalar_type_of_t<Args>...>>
526 #else
527  template<typename...Ds, typename Operation, typename...Args, std::enable_if_t<(sizeof...(Args) <= 3) and
528  std::is_invocable<Operation&&, typename scalar_type_of<Args>::type...>::value and
529  values::number<typename std::invoke_result<Operation&&, typename scalar_type_of<Args>::type...>::type>, int> = 0>
530 #endif
531  static auto
532  n_ary_operation(const std::tuple<Ds...>& tup, Operation&& operation, Args&&...args)
533  {
534  auto&& op = nat_op(std::forward<Operation>(operation));
535  using Op = decltype(op);
536 
537  if constexpr (sizeof...(Args) == 0)
538  {
539  using P = dense_writable_matrix_t<T, data_layout::none, scalar_type_of_t<T>, std::tuple<Ds...>>;
540  IndexType r = get_dimension(std::get<0>(tup));
541  IndexType c = get_dimension(std::get<1>(tup));
542  return Eigen::CwiseNullaryOp<std::decay_t<Op>, P> {r, c, std::forward<Op>(op)};
543  }
544  else if constexpr (sizeof...(Args) == 1)
545  {
546  return make_self_contained<Args...>(Eigen::CwiseUnaryOp {
547  replicate_arg(tup, std::forward<Args>(args))..., std::forward<Op>(op)});
548  }
549  else if constexpr (sizeof...(Args) == 2)
550  {
551  return make_self_contained<Args...>(Eigen::CwiseBinaryOp<std::decay_t<Op>,
552  std::decay_t<decltype(replicate_arg(tup, std::forward<Args>(args)))>...> {
553  replicate_arg(tup, std::forward<Args>(args))..., std::forward<Op>(op)});
554  }
555  else
556  {
557  return make_self_contained<Args...>(Eigen::CwiseTernaryOp<std::decay_t<Op>,
558  std::decay_t<decltype(replicate_arg(tup, std::forward<Args>(args)))>...> {
559  replicate_arg(tup, std::forward<Args>(args))..., std::forward<Op>(op)});
560  }
561  }
562  */
563 
564 
565  /*
566  template<std::size_t...indices, typename BinaryFunction, typename Arg>
567  static constexpr auto
568  reduce(BinaryFunction&& b, Arg&& arg)
569  {
570  if constexpr (Eigen3::eigen_dense_general<Arg>)
571  {
572  auto&& op = nat_op(std::forward<BinaryFunction>(b));
573  using Op = decltype(op);
574 
575  if constexpr (sizeof...(indices) == 2) // reduce in both directions
576  {
577  return std::forward<Arg>(arg).redux(std::forward<Op>(op));
578  }
579  else
580  {
581  using OpWrapper = Eigen::internal::member_redux<std::decay_t<Op>, scalar_type_of_t<Arg>>;
582  constexpr auto dir = ((indices == 0) and ...) ? Eigen::Vertical : Eigen::Horizontal;
583  using P = Eigen::PartialReduxExpr<std::decay_t<Arg>, OpWrapper, dir>;
584  return make_self_contained<Arg>(P {std::forward<Arg>(arg), OpWrapper {std::forward<Op>(op)}});
585  }
586  }
587  else
588  {
589  return reduce<indices...>(std::forward<BinaryFunction>(b), Eigen3::make_eigen_wrapper(std::forward<Arg>(arg)));
590  }
591  }
592 
593  // to_euclidean not defined--rely on default
594 
595  // from_euclidean not defined--rely on default
596 
597  // wrap_angles not defined--rely on default
598 
599 
600  template<typename Arg>
601  static constexpr decltype(auto)
602  conjugate(Arg&& arg)
603  {
604  // The global conjugate function already handles Eigen::DiagonalMatrix and Eigen::DiagonalWrapper
605  return std::forward<Arg>(arg).conjugate();
606  }
607 
608 
609  template<typename Arg>
610  static constexpr decltype(auto)
611  transpose(Arg&& arg)
612  {
613  if constexpr (Eigen3::eigen_wrapper<Arg>)
614  {
615  if constexpr (Eigen3::eigen_general<nested_object_of_t<Arg>, true>)
616  return transpose(nested_object(std::forward<Arg>(arg)));
617  else
618  return std::forward<Arg>(arg).transpose(); // Rely on inherited Eigen transpose method
619  }
620  else if constexpr (Eigen3::eigen_matrix_general<Arg, true> or Eigen3::eigen_TriangularView<Arg> or Eigen3::eigen_SelfAdjointView<Arg>)
621  return std::forward<Arg>(arg).transpose();
622  else if constexpr (Eigen3::eigen_array_general<Arg, true>)
623  return std::forward<Arg>(arg).matrix().transpose();
624  else if constexpr (triangular_matrix<Arg>)
625  return OpenKalman::transpose(TriangularAdapter {std::forward<Arg>(arg)});
626  else
627  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg)).transpose();
628  // Note: the global transpose function already handles zero, constant, constant-diagonal, and symmetric cases.
629  }
630 
631 
632  template<typename Arg>
633  static constexpr decltype(auto)
634  adjoint(Arg&& arg)
635  {
636  if constexpr (Eigen3::eigen_wrapper<Arg>)
637  {
638  if constexpr (Eigen3::eigen_general<nested_object_of_t<Arg>, true>)
639  return adjoint(nested_object(std::forward<Arg>(arg)));
640  else
641  return std::forward<Arg>(arg).adjoint(); // Rely on inherited Eigen adjoint method
642  }
643  else if constexpr (Eigen3::eigen_matrix_general<Arg, true> or Eigen3::eigen_TriangularView<Arg> or Eigen3::eigen_SelfAdjointView<Arg>)
644  return std::forward<Arg>(arg).adjoint();
645  else if constexpr (Eigen3::eigen_array_general<Arg, true>)
646  return std::forward<Arg>(arg).matrix().adjoint();
647  else if constexpr (triangular_matrix<Arg>)
648  return OpenKalman::adjoint(TriangularAdapter {std::forward<Arg>(arg)});
649  else
650  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg)).adjoint();
651  // Note: the global adjoint function already handles zero, constant, diagonal, non-complex, and hermitian cases.
652  }
653 
654 
655  template<typename Arg>
656  static constexpr auto
657  determinant(Arg&& arg)
658  {
659  if constexpr (Eigen3::eigen_matrix_general<Arg, true>)
660  return std::forward<Arg>(arg).determinant();
661  else if constexpr (Eigen3::eigen_array_general<Arg, true>)
662  return std::forward<Arg>(arg).matrix().determinant();
663  else
664  return Eigen3::make_eigen_wrapper(std::forward<Arg>(arg)).determinant();
665  // Note: the global determinant function already handles Eigen::TriangularView, Eigen::DiagonalMatrix, and Eigen::DiagonalWrapper
666  }
667 
668 
669  template<typename A, typename B>
670  static constexpr auto
671  sum(A&& a, B&& b)
672  {
673  auto s = make_self_contained<A, B>(std::forward<A>(a) + std::forward<B>(b));
674  if constexpr ((dynamic_dimension<decltype(s), 0> and (not dynamic_dimension<A, 0> or not dynamic_dimension<B, 0>)) or
675  (dynamic_dimension<decltype(s), 1> and (not dynamic_dimension<A, 1> or not dynamic_dimension<B, 1>)))
676  {
677  using S0 = vector_space_descriptor_of_t<decltype(s), 0>; using S1 = vector_space_descriptor_of_t<decltype(s), 1>;
678  using A0 = vector_space_descriptor_of_t<A, 0>; using A1 = vector_space_descriptor_of_t<A, 1>;
679  using B0 = vector_space_descriptor_of_t<B, 0>; using B1 = vector_space_descriptor_of_t<B, 1>;
680 
681  using R = std::conditional_t<dynamic_pattern<S0>,
682  std::conditional_t<dynamic_pattern<A0>, std::conditional_t<dynamic_pattern<B0>, S0, B0>, A0>, S0>;
683  using C = std::conditional_t<dynamic_pattern<S1>,
684  std::conditional_t<dynamic_pattern<A1>, std::conditional_t<dynamic_pattern<B1>, S1, B1>, A1>, S1>;
685  return internal::FixedSizeAdapter<std::decay_t<decltype(s)>, R, C> {std::move(s)};
686  }
687  else return s;
688  }
689 
690 
691  template<typename A, typename B>
692  static constexpr auto
693  contract(A&& a, B&& b)
694  {
695  if constexpr (diagonal_matrix<A> and internal::has_nested_vector<A>)
696  {
697  if constexpr (Eigen3::eigen_DiagonalWrapper<A> or Eigen3::eigen_DiagonalMatrix<A>)
698  return make_self_contained<A, B>(std::forward<A>(a) * OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
699  else
700  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(nested_object(std::forward<A>(a))).asDiagonal() *
701  OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
702  }
703  else if constexpr (diagonal_matrix<B> and internal::has_nested_vector<B>)
704  {
705  if constexpr (Eigen3::eigen_DiagonalWrapper<B> or Eigen3::eigen_DiagonalMatrix<B>)
706  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(std::forward<A>(a)) * std::forward<B>(b));
707  else
708  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(std::forward<A>(a)) *
709  OpenKalman::to_native_matrix<T>(nested_object(std::forward<B>(b))).asDiagonal());
710  }
711  else if constexpr (triangular_adapter<A>)
712  {
713  constexpr auto uplo = triangular_matrix<A, triangle_type::upper> ? Eigen::Upper : Eigen::Lower;
714  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(nested_object(std::forward<A>(a))).template triangularView<uplo>() *
715  OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
716  }
717  else if constexpr (triangular_adapter<B>)
718  {
719  constexpr auto uplo = triangular_matrix<A, triangle_type::upper> ? Eigen::Upper : Eigen::Lower;
720  auto prod = OpenKalman::to_native_matrix<T>(std::forward<A>(a));
721  prod.applyOnTheRight(OpenKalman::to_native_matrix<T>(nested_object(std::forward<B>(b))).template triangularView<uplo>());
722  return prod;
723  }
724  else if constexpr (hermitian_adapter<A>)
725  {
726  constexpr auto uplo = hermitian_adapter<A, HermitianAdapterType::upper> ? Eigen::Upper : Eigen::Lower;
727  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(nested_object(std::forward<A>(a))).template selfadjointView<uplo>() *
728  OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
729  }
730  else if constexpr (hermitian_adapter<B>)
731  {
732  constexpr auto uplo = triangular_matrix<A, triangle_type::upper> ? Eigen::Upper : Eigen::Lower;
733  auto prod = OpenKalman::to_native_matrix<T>(std::forward<A>(a));
734  prod.applyOnTheRight(OpenKalman::to_native_matrix<T>(nested_object(std::forward<B>(b))).template selfadjointView<uplo>());
735  return prod;
736  }
737  else
738  {
739  return make_self_contained<A, B>(OpenKalman::to_native_matrix<T>(std::forward<A>(a)) *
740  OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
741  }
742  }
743 
744 
745 #ifdef __cpp_concepts
746  template<bool on_the_right, writable A, indexible B> requires Eigen3::eigen_dense_general<A>
747 #else
748  template<bool on_the_right, typename A, typename B, std::enable_if_t<writable<A> and Eigen3::eigen_dense_general<A>,int> = 0>
749 #endif
750  static A&
751  contract_in_place(A& a, B&& b)
752  {
753  auto&& ma = [](A& a){
754  if constexpr (Eigen3::eigen_array_general<A, true>) return a.matrix();
755  else return (a);
756  }(a);
757 
758  if constexpr (on_the_right)
759  return ma.applyOnTheRight(OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
760  else
761  return ma.applyOnTheLeft(OpenKalman::to_native_matrix<T>(std::forward<B>(b)));
762  return a;
763  }
764 
765 
766  template<triangle_type tri, typename A>
767  static constexpr auto
768  cholesky_factor(A&& a)
769  {
770  using NestedMatrix = std::decay_t<nested_object_of_t<A>>;
771  using Scalar = scalar_type_of_t<A>;
772  constexpr auto dim = index_dimension_of_v<A, 0>;
773  using M = dense_writable_matrix_t<A>;
774 
775  if constexpr (std::is_same_v<
776  const NestedMatrix, const typename Eigen::MatrixBase<NestedMatrix>::ConstantReturnType>)
777  {
778  // If nested matrix is a positive constant matrix, construct the Cholesky factor using a shortcut.
779 
780  auto s = nested_object(std::forward<A>(a)).functor()();
781 
782  if (s < Scalar(0))
783  {
784  // Cholesky factor elements are complex, so throw an exception.
785  throw (std::runtime_error("cholesky_factor of constant HermitianAdapter: covariance is indefinite"));
786  }
787 
788  if constexpr(tri == triangle_type::diagonal)
789  {
790  static_assert(diagonal_matrix<A>);
791  auto vec = make_constant<A>(square_root(s), Dimensions<dim>{}, Dimensions<1>{});
792  return diagonal_adapter<decltype(vec)> {vec};
793  }
794  else if constexpr(tri == triangle_type::lower)
795  {
796  auto col0 = make_constant<A>(square_root(s), Dimensions<dim>{}, Dimensions<1>{});
797  auto othercols = make_zero<A>(get_pattern_collection<0>(a), get_pattern_collection<0>(a) - 1);
798  return TriangularAdapter<M, tri> {concatenate_horizontal(col0, othercols)};
799  }
800  else
801  {
802  static_assert(tri == triangle_type::upper);
803  auto row0 = make_constant<A>(square_root(s), Dimensions<1>{}, Dimensions<dim>{});
804  auto otherrows = make_zero<A>(get_pattern_collection<0>(a) - 1, get_pattern_collection<0>(a));
805  return TriangularAdapter<M, tri> {concatenate_vertical(row0, otherrows)};
806  }
807  }
808  else
809  {
810  // For the general case, perform an LLT Cholesky decomposition.
811  M b;
812  auto LL_x = a.view().llt();
813  if (LL_x.info() == Eigen::Success)
814  {
815  if constexpr(tri == hermitian_adapter_type_of_v<A>)
816  {
817  b = std::move(LL_x.matrixLLT());
818  }
819  else
820  {
821  constexpr unsigned int uplo = tri == triangle_type::upper ? Eigen::Upper : Eigen::Lower;
822  b.template triangularView<uplo>() = LL_x.matrixLLT().adjoint();
823  }
824  }
825  else [[unlikely]]
826  {
827  // If covariance is not positive definite, use the more robust LDLT decomposition.
828  auto LDL_x = nested_object(std::forward<A>(a)).ldlt();
829  if ((not LDL_x.isPositive() and not LDL_x.isNegative()) or LDL_x.info() != Eigen::Success) [[unlikely]]
830  {
831  if (LDL_x.isPositive() and LDL_x.isNegative()) // Covariance is zero, even though decomposition failed.
832  {
833  if constexpr(tri == triangle_type::lower)
834  b.template triangularView<Eigen::Lower>() = make_zero(nested_object(a));
835  else
836  b.template triangularView<Eigen::Upper>() = make_zero(nested_object(a));
837  }
838  else // Covariance is indefinite, so throw an exception.
839  {
840  throw (std::runtime_error("cholesky_factor of HermitianAdapter: covariance is indefinite"));
841  }
842  }
843  else if constexpr(tri == triangle_type::lower)
844  {
845  b.template triangularView<Eigen::Lower>() =
846  LDL_x.matrixL().toDenseMatrix() * LDL_x.vectorD().cwiseSqrt().asDiagonal();
847  }
848  else
849  {
850  b.template triangularView<Eigen::Upper>() =
851  LDL_x.vectorD().cwiseSqrt().asDiagonal() * LDL_x.matrixU().toDenseMatrix();
852  }
853  }
854  return TriangularAdapter<M, tri> {std::move(b)};
855  }
856  }
857 
858 
859  template<HermitianAdapterType significant_triangle, typename A, typename U, typename Alpha>
860  static decltype(auto)
861  rank_update_hermitian(A&& a, U&& u, const Alpha alpha)
862  {
863  if constexpr (OpenKalman::Eigen3::eigen_ArrayWrapper<A>)
864  {
865  return rank_update_hermitian<significant_triangle>(nested_object(std::forward<A>(a)), std::forward<U>(u), alpha);
866  }
867  else if constexpr (OpenKalman::Eigen3::eigen_array_general<A, true>)
868  {
869  return rank_update_hermitian<significant_triangle>(std::forward<A>(a).matrix(), std::forward<U>(u), alpha);
870  }
871  else
872  {
873  static_assert(writable<A>);
874  constexpr auto s = significant_triangle == HermitianAdapterType::lower ? Eigen::Lower : Eigen::Upper;
875  a.template selfadjointView<s>().template rankUpdate(std::forward<U>(u), alpha);
876  return std::forward<A>(a);
877  }
878  }
879 
880 
881  template<triangle_type triangle, typename A, typename U, typename Alpha>
882  static decltype(auto)
883  rank_update_triangular(A&& a, U&& u, const Alpha alpha)
884  {
885  if constexpr (OpenKalman::Eigen3::eigen_ArrayWrapper<A>)
886  {
887  return rank_update_triangular<triangle>(nested_object(std::forward<A>(a)), std::forward<U>(u), alpha);
888  }
889  else if constexpr (OpenKalman::Eigen3::eigen_array_general<A, true>)
890  {
891  return rank_update_triangular<triangle>(std::forward<A>(a).matrix(), std::forward<U>(u), alpha);
892  }
893  else
894  {
895  static_assert(writable<A>);
896  constexpr auto t = triangle == triangle_type::lower ? Eigen::Lower : Eigen::Upper;
897  using Scalar = scalar_type_of_t<A>;
898  for (std::size_t i = 0; i < get_index_dimension_of<1>(u); i++)
899  {
900  if (Eigen::internal::llt_inplace<Scalar, t>::rankUpdate(a, get_chip<1>(u, i), alpha) >= 0)
901  throw (std::runtime_error("rank_update_triangular: product is not positive definite"));
902  }
903  return std::forward<A>(a);
904  }
905  }
906 
907 
908  template<bool must_be_unique, bool must_be_exact, typename A, typename B>
909  static constexpr auto
910  solve(A&& a, B&& b)
911  {
912  using Scalar = scalar_type_of_t<A>;
913 
914  constexpr std::size_t a_rows = dynamic_dimension<A, 0> ? index_dimension_of_v<B, 0> : index_dimension_of_v<A, 0>;
915  constexpr std::size_t a_cols = index_dimension_of_v<A, 1>;
916  constexpr std::size_t b_cols = index_dimension_of_v<B, 1>;
917 
918  if constexpr (not Eigen3::eigen_matrix_general<A, true>)
919  {
920  auto&& n = OpenKalman::to_native_matrix(std::forward<A>(a));
921  static_assert(Eigen3::eigen_matrix_general<decltype(n), true>);
922  return solve<must_be_unique, must_be_exact>(std::forward<decltype(n)>(n), std::forward<B>(b));
923  }
924  else if constexpr (triangular_matrix<A>)
925  {
926  constexpr auto uplo = triangular_matrix<A, triangle_type::upper> ? Eigen::Upper : Eigen::Lower;
927  return make_self_contained<A, B>(
928  Eigen::Solve {std::forward<A>(a).template triangularView<uplo>(), std::forward<B>(b)});
929  }
930  else if constexpr (hermitian_matrix<A>)
931  {
932  constexpr auto uplo = hermitian_adapter_type_of_v<A> == triangle_type::upper ? Eigen::Upper : Eigen::Lower;
933  auto v {std::forward<A>(a).template selfadjointView<uplo>()};
934  auto llt {v.llt()};
935 
936  Eigen3::eigen_matrix_t<Scalar, a_cols, b_cols> ret;
937  if (llt.info() == Eigen::Success)
938  {
939  ret = Eigen::Solve {llt, std::forward<B>(b)};
940  }
941  else [[unlikely]]
942  {
943  // A is semidefinite. Use LDLT decomposition instead.
944  auto ldlt {v.ldlt()};
945  if ((not ldlt.isPositive() and not ldlt.isNegative()) or ldlt.info() != Eigen::Success)
946  {
947  throw (std::runtime_error("Eigen solve (hermitian case): A is indefinite"));
948  }
949  ret = Eigen::Solve {ldlt, std::forward<B>(b)};
950  }
951  return ret;
952  }
953  else
954  {
955  if constexpr (must_be_exact or must_be_unique or true)
956  {
957  auto a_cols_rt = get_index_dimension_of<1>(a);
958  Eigen::ColPivHouseholderQR<Eigen3::eigen_matrix_t<Scalar, a_rows, a_cols>> QR {std::forward<A>(a)};
959  if constexpr (must_be_unique)
960  {
961  if (QR.rank() < a_cols_rt) throw std::runtime_error {"solve function requests a "
962  "unique solution, but A is rank-deficient, so result X is not unique"};
963  }
964 
965  auto res = QR.solve(std::forward<B>(b));
966 
967  if constexpr (must_be_exact)
968  {
969  bool a_solution_exists = (a*res).isApprox(b, a_cols_rt * std::numeric_limits<scalar_type_of_t<A>>::epsilon());
970 
971  if (a_solution_exists)
972  return make_self_contained(std::move(res));
973  else
974  throw std::runtime_error {"solve function requests an exact solution, "
975  "but the solution is only an approximation"};
976  }
977  else
978  {
979  return make_self_contained(std::move(res));
980  }
981  }
982  else
983  {
984  Eigen::HouseholderQR<Eigen3::eigen_matrix_t<Scalar, a_rows, a_cols>> QR {std::forward<A>(a)};
985  return make_self_contained(QR.solve(std::forward<B>(b)));
986  }
987  }
988  }
989 
990  private:
991 
992  template<typename A>
993  static constexpr auto
994  QR_decomp_impl(A&& a)
995  {
996  using Scalar = scalar_type_of_t<A>;
997  constexpr auto rows = index_dimension_of_v<A, 0>;
998  constexpr auto cols = index_dimension_of_v<A, 1>;
999  using MatrixType = Eigen3::eigen_matrix_t<Scalar, rows, cols>;
1000  using ResultType = Eigen3::eigen_matrix_t<Scalar, cols, cols>;
1001 
1002  Eigen::HouseholderQR<MatrixType> QR {std::forward<A>(a)};
1003 
1004  if constexpr (dynamic_dimension<A, 1>)
1005  {
1006  auto rt_cols = get_index_dimension_of<1>(a);
1007 
1008  ResultType ret {rt_cols, rt_cols};
1009 
1010  if constexpr (dynamic_dimension<A, 0>)
1011  {
1012  auto rt_rows = get_index_dimension_of<0>(a);
1013 
1014  if (rt_rows < rt_cols)
1015  ret << QR.matrixQR().topRows(rt_rows),
1016  Eigen3::eigen_matrix_t<Scalar, stdex::dynamic_extent, stdex::dynamic_extent>::Zero(rt_cols - rt_rows, rt_cols);
1017  else
1018  ret = QR.matrixQR().topRows(rt_cols);
1019  }
1020  else
1021  {
1022  if (rows < rt_cols)
1023  ret << QR.matrixQR().template topRows<rows>(),
1024  Eigen3::eigen_matrix_t<Scalar, stdex::dynamic_extent, stdex::dynamic_extent>::Zero(rt_cols - rows, rt_cols);
1025  else
1026  ret = QR.matrixQR().topRows(rt_cols);
1027  }
1028 
1029  return ret;
1030  }
1031  else
1032  {
1033  ResultType ret;
1034 
1035  if constexpr (dynamic_dimension<A, 0>)
1036  {
1037  auto rt_rows = get_index_dimension_of<0>(a);
1038 
1039  if (rt_rows < cols)
1040  ret << QR.matrixQR().topRows(rt_rows),
1041  Eigen3::eigen_matrix_t<Scalar, stdex::dynamic_extent, stdex::dynamic_extent>::Zero(cols - rt_rows, cols);
1042  else
1043  ret = QR.matrixQR().template topRows<cols>();
1044  }
1045  else
1046  {
1047  if constexpr (rows < cols)
1048  ret << QR.matrixQR().template topRows<rows>(), Eigen3::eigen_matrix_t<Scalar, cols - rows, cols>::Zero();
1049  else
1050  ret = QR.matrixQR().template topRows<cols>();
1051  }
1052 
1053  return ret;
1054  }
1055  }
1056 
1057  public:
1058 
1059  template<typename A>
1060  static constexpr auto
1061  LQ_decomposition(A&& a)
1062  {
1063  return make_triangular_matrix<triangle_type::lower>(make_self_contained(adjoint(QR_decomp_impl(adjoint(std::forward<A>(a))))));
1064  }
1065 
1066 
1067  template<typename A>
1068  static constexpr auto
1069  QR_decomposition(A&& a)
1070  {
1071  return make_triangular_matrix<triangle_type::upper>(QR_decomp_impl(std::forward<A>(a)));
1072  }*/
1073 
1074  };
1075 
1076 
1077 }
1078 
1079 #endif
Definition: basics.hpp:41
Forward declarations for OpenKalman&#39;s Eigen Tensor module interface.
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
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 value
T is a fixed or dynamic value that is reducible to a number.
Definition: value.hpp:45
decltype(auto) constexpr apply(F &&f, T &&t)
A generalization of std::apply.
Definition: apply.hpp:49
constexpr auto get_dimension(const Arg &arg)
Get the vector dimension of coordinates::pattern Arg.
Definition: get_dimension.hpp:54
An interface to various routines from the linear algebra library associated with indexible object T...
Definition: library_interface.hpp:42
Definition: eigen-tensor-forward-declarations.hpp:30
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:65
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 bool dynamic
T is a value that is not fixed at compile time.
Definition: dynamic.hpp:30
Definition: trait_backports.hpp:64
constexpr auto make_constant(C c, stdex::extents< IndexType, Extents... > extents)
Make an indexible object in which every element is a constant value.
Definition: make_constant.hpp:39
decltype(auto) constexpr nested_object(Arg &&arg)
Retrieve a nested object of Arg, if it exists.
Definition: nested_object.hpp:35
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