10#ifndef TLAPACK_EIGEN_HH
11#define TLAPACK_EIGEN_HH
28 template <
class Derived>
29 std::true_type is_eigen_dense_f(
const Eigen::DenseBase<Derived>*);
30 std::false_type is_eigen_dense_f(
const void*);
32 template <
class Derived>
33 std::true_type is_eigen_matrix_f(
const Eigen::MatrixBase<Derived>*);
34 std::false_type is_eigen_matrix_f(
const void*);
36 template <
class XprType,
int BlockRows,
int BlockCols,
bool InnerPanel>
37 std::true_type is_eigen_block_f(
38 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>*);
39 std::false_type is_eigen_block_f(
const void*);
45 decltype(is_eigen_dense_f(std::declval<T*>()))::value;
51 decltype(is_eigen_matrix_f(std::declval<T*>()))::value;
57 decltype(is_eigen_block_f(std::declval<T*>()))::value;
79 template <
class matrix_t>
83 (matrix_t::InnerStrideAtCompileTime == 1 ||
84 matrix_t::OuterStrideAtCompileTime == 1),
86 static constexpr Layout value =
87 (matrix_t::IsVectorAtCompileTime)
89 : ((matrix_t::IsRowMajor) ? Layout::RowMajor
93 template <
class matrix_t>
97 using type = Eigen::Matrix<real_type<typename matrix_t::Scalar>,
98 matrix_t::RowsAtCompileTime,
99 matrix_t::ColsAtCompileTime,
100 matrix_t::IsRowMajor ? Eigen::RowMajor
102 matrix_t::MaxRowsAtCompileTime,
103 matrix_t::MaxColsAtCompileTime>;
106 template <
class matrix_t>
110 using type = Eigen::Matrix<complex_type<typename matrix_t::Scalar>,
111 matrix_t::RowsAtCompileTime,
112 matrix_t::ColsAtCompileTime,
113 matrix_t::IsRowMajor ? Eigen::RowMajor
115 matrix_t::MaxRowsAtCompileTime,
116 matrix_t::MaxColsAtCompileTime>;
120 template <
typename U>
123 static constexpr int Rows_ =
124 (U::RowsAtCompileTime == 1) ? 1 : Eigen::Dynamic;
125 static constexpr int Cols_ =
126 (U::ColsAtCompileTime == 1) ? 1 : Eigen::Dynamic;
127 static constexpr int Options_ =
128 (U::IsRowMajor) ? Eigen::RowMajor : Eigen::ColMajor;
130 template <
typename T>
133 Eigen::Index n = 1)
const
137 return Eigen::Matrix<T, Rows_, Cols_, Options_>(m, n);
142 template <
typename U,
int m,
int n>
148 static constexpr int Options_ =
149 (U::IsRowMajor) ? Eigen::RowMajor : Eigen::ColMajor;
150 static_assert(m >= 0 && n >= -1);
152 template <
typename T>
155 constexpr int Cols_ = (n == -1) ? 1 : n;
156 return Eigen::Matrix<T, m, Cols_, Options_, m, Cols_>();
172#if __cplusplus >= 201703L
180 const Eigen::Matrix<T, Rows, Cols, Options, MaxRows, MaxCols>&
x)
noexcept
184template <
class Derived>
185constexpr auto size(
const Eigen::EigenBase<Derived>& x)
noexcept
191constexpr auto nrows(
const Eigen::EigenBase<T>& x)
noexcept
197constexpr auto ncols(
const Eigen::EigenBase<T>& x)
noexcept
218#define isSlice(SliceSpec) !std::is_convertible<SliceSpec, Eigen::Index>::value
226 typename std::enable_if<isSlice(SliceSpecRow) && isSlice(SliceSpecCol) &&
227 traits::is_eigen_type<T> &&
228 !traits::internal::is_eigen_block<T>,
230constexpr auto slice(T& A, SliceSpecRow&& rows, SliceSpecCol&& cols)
noexcept
232 return A.block(rows.first, cols.first, rows.second - rows.first,
233 cols.second - cols.first);
237 typename SliceSpecCol,
238 typename std::enable_if<traits::is_eigen_type<T> &&
239 !traits::internal::is_eigen_block<T>,
241constexpr auto slice(T& A, Eigen::Index rowIdx, SliceSpecCol&& cols)
noexcept
243 return A.row(rowIdx).segment(cols.first, cols.second - cols.first);
247 typename SliceSpecRow,
248 typename std::enable_if<traits::is_eigen_type<T> &&
249 !traits::internal::is_eigen_block<T>,
251constexpr auto slice(T& A, SliceSpecRow&& rows, Eigen::Index colIdx)
noexcept
253 return A.col(colIdx).segment(rows.first, rows.second - rows.first);
258 typename std::enable_if<traits::is_eigen_type<T> &&
259 !traits::internal::is_eigen_block<T>,
261constexpr auto slice(T& x, SliceSpec&& range)
noexcept
263 return x.segment(range.first, range.second - range.first);
268 typename std::enable_if<traits::is_eigen_type<T> &&
269 !traits::internal::is_eigen_block<T>,
271constexpr auto rows(T& A, SliceSpec&& rows)
noexcept
273 return A.middleRows(rows.first, rows.second - rows.first);
277 typename std::enable_if<traits::is_eigen_type<T> &&
278 !traits::internal::is_eigen_block<T>,
280constexpr auto row(T& A, Eigen::Index rowIdx)
noexcept
282 return A.row(rowIdx);
287 typename std::enable_if<traits::is_eigen_type<T> &&
288 !traits::internal::is_eigen_block<T>,
290constexpr auto cols(T& A, SliceSpec&& cols)
noexcept
292 return A.middleCols(cols.first, cols.second - cols.first);
296 typename std::enable_if<traits::is_eigen_type<T> &&
297 !traits::internal::is_eigen_block<T>,
299constexpr auto col(T& A, Eigen::Index colIdx)
noexcept
301 return A.col(colIdx);
313 typename std::enable_if<isSlice(SliceSpecRow) && isSlice(SliceSpecCol),
315constexpr auto slice(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
317 SliceSpecCol&& cols)
noexcept
319 assert(rows.second <= A.rows());
320 assert(cols.second <= A.cols());
322 return Eigen::Block<XprType, Eigen::Dynamic, Eigen::Dynamic, InnerPanel>(
323 A.nestedExpression(), A.startRow() + rows.first,
324 A.startCol() + cols.first, rows.second - rows.first,
325 cols.second - cols.first);
328template <
class XprType,
332 typename SliceSpecCol>
333constexpr auto slice(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
335 SliceSpecCol&& cols)
noexcept
337 assert(rowIdx < A.rows());
338 assert(cols.second <= A.cols());
340 return Eigen::Block<XprType, 1, Eigen::Dynamic, InnerPanel>(
341 A.nestedExpression(), A.startRow() + rowIdx, A.startCol() + cols.first,
342 1, cols.second - cols.first);
345template <
class XprType,
349 typename SliceSpecRow>
350constexpr auto slice(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
352 Eigen::Index colIdx)
noexcept
354 assert(rows.second <= A.rows());
355 assert(colIdx < A.cols());
357 return Eigen::Block<XprType, Eigen::Dynamic, 1, InnerPanel>(
358 A.nestedExpression(), A.startRow() + rows.first, A.startCol() + colIdx,
359 rows.second - rows.first, 1);
362template <
class XprType,
int BlockRows,
bool InnerPanel,
typename SliceSpec>
363constexpr auto slice(Eigen::Block<XprType, BlockRows, 1, InnerPanel>& x,
364 SliceSpec&& range)
noexcept
366 assert(range.second <= x.size());
368 return Eigen::Block<XprType, Eigen::Dynamic, 1, InnerPanel>(
369 x.nestedExpression(), x.startRow() + range.first, x.startCol(),
370 range.second - range.first, 1);
373template <
class XprType,
int BlockCols,
bool InnerPanel,
typename SliceSpec>
374constexpr auto slice(Eigen::Block<XprType, 1, BlockCols, InnerPanel>& x,
375 SliceSpec&& range)
noexcept
377 assert(range.second <= x.size());
379 return Eigen::Block<XprType, 1, Eigen::Dynamic, InnerPanel>(
380 x.nestedExpression(), x.startRow(), x.startCol() + range.first, 1,
381 range.second - range.first);
384template <
class XprType,
389constexpr auto rows(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
390 SliceSpec&& rows)
noexcept
392 assert(rows.second <= A.rows());
394 return Eigen::Block<XprType, Eigen::Dynamic, BlockCols, InnerPanel>(
395 A.nestedExpression(), A.startRow() + rows.first, A.startCol(),
396 rows.second - rows.first, A.cols());
399template <
class XprType,
int BlockRows,
int BlockCols,
bool InnerPanel>
400constexpr auto row(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
401 Eigen::Index rowIdx)
noexcept
403 assert(rowIdx < A.rows());
405 return Eigen::Block<XprType, 1, BlockCols, InnerPanel>(
406 A.nestedExpression(), A.startRow() + rowIdx, A.startCol(), 1, A.cols());
409template <
class XprType,
414constexpr auto cols(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
415 SliceSpec&& cols)
noexcept
417 assert(cols.second <= A.cols());
419 return Eigen::Block<XprType, BlockRows, Eigen::Dynamic, InnerPanel>(
420 A.nestedExpression(), A.startRow(), A.startCol() + cols.first, A.rows(),
421 cols.second - cols.first);
424template <
class XprType,
int BlockRows,
int BlockCols,
bool InnerPanel>
425constexpr auto col(Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
426 Eigen::Index colIdx)
noexcept
428 assert(colIdx < A.cols());
430 return Eigen::Block<XprType, BlockRows, 1, InnerPanel>(
431 A.nestedExpression(), A.startRow(), A.startCol() + colIdx, A.rows(), 1);
441 typename std::enable_if<isSlice(SliceSpecRow) && isSlice(SliceSpecCol),
444 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
446 SliceSpecCol&& cols)
noexcept
448 assert(rows.second <= A.rows());
449 assert(cols.second <= A.cols());
451 return Eigen::Block<
const XprType, Eigen::Dynamic, Eigen::Dynamic,
453 A.nestedExpression(), A.startRow() + rows.first,
454 A.startCol() + cols.first, rows.second - rows.first,
455 cols.second - cols.first);
458template <
class XprType,
462 typename SliceSpecCol>
464 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
466 SliceSpecCol&& cols)
noexcept
468 assert(rowIdx < A.rows());
469 assert(cols.second <= A.cols());
471 return Eigen::Block<const XprType, 1, Eigen::Dynamic, InnerPanel>(
472 A.nestedExpression(), A.startRow() + rowIdx, A.startCol() + cols.first,
473 1, cols.second - cols.first);
476template <
class XprType,
480 typename SliceSpecRow>
482 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
484 Eigen::Index colIdx)
noexcept
486 assert(rows.second <= A.rows());
487 assert(colIdx < A.cols());
489 return Eigen::Block<const XprType, Eigen::Dynamic, 1, InnerPanel>(
490 A.nestedExpression(), A.startRow() + rows.first, A.startCol() + colIdx,
491 rows.second - rows.first, 1);
494template <
class XprType,
int BlockRows,
bool InnerPanel,
typename SliceSpec>
495constexpr auto slice(
const Eigen::Block<XprType, BlockRows, 1, InnerPanel>& x,
496 SliceSpec&& range)
noexcept
498 assert(range.second <= x.size());
500 return Eigen::Block<const XprType, Eigen::Dynamic, 1, InnerPanel>(
501 x.nestedExpression(), x.startRow() + range.first, x.startCol(),
502 range.second - range.first, 1);
505template <
class XprType,
int BlockCols,
bool InnerPanel,
typename SliceSpec>
506constexpr auto slice(
const Eigen::Block<XprType, 1, BlockCols, InnerPanel>& x,
507 SliceSpec&& range)
noexcept
509 assert(range.second <= x.size());
511 return Eigen::Block<const XprType, 1, Eigen::Dynamic, InnerPanel>(
512 x.nestedExpression(), x.startRow(), x.startCol() + range.first, 1,
513 range.second - range.first);
516template <
class XprType,
522 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
523 SliceSpec&& rows)
noexcept
525 assert(rows.second <= A.rows());
527 return Eigen::Block<const XprType, Eigen::Dynamic, BlockCols, InnerPanel>(
528 A.nestedExpression(), A.startRow() + rows.first, A.startCol(),
529 rows.second - rows.first, A.cols());
532template <
class XprType,
int BlockRows,
int BlockCols,
bool InnerPanel>
534 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
535 Eigen::Index rowIdx)
noexcept
537 assert(rowIdx < A.rows());
539 return Eigen::Block<const XprType, 1, BlockCols, InnerPanel>(
540 A.nestedExpression(), A.startRow() + rowIdx, A.startCol(), 1, A.cols());
543template <
class XprType,
549 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
550 SliceSpec&& cols)
noexcept
552 assert(cols.second <= A.cols());
554 return Eigen::Block<const XprType, BlockRows, Eigen::Dynamic, InnerPanel>(
555 A.nestedExpression(), A.startRow(), A.startCol() + cols.first, A.rows(),
556 cols.second - cols.first);
559template <
class XprType,
int BlockRows,
int BlockCols,
bool InnerPanel>
561 const Eigen::Block<XprType, BlockRows, BlockCols, InnerPanel>& A,
562 Eigen::Index colIdx)
noexcept
564 assert(colIdx < A.cols());
566 return Eigen::Block<const XprType, BlockRows, 1, InnerPanel>(
567 A.nestedExpression(), A.startRow(), A.startCol() + colIdx, A.rows(), 1);
574 typename std::enable_if<traits::internal::is_eigen_matrix<T>,
582template <
class matrix_t,
583 typename std::enable_if<(traits::is_eigen_type<matrix_t> &&
584 !matrix_t::IsVectorAtCompileTime),
586constexpr auto transpose_view(matrix_t& A)
noexcept
588 using T =
typename matrix_t::Scalar;
589 using Stride = Eigen::OuterStride<>;
590 assert(A.innerStride() == 1);
592 constexpr int Rows_ = matrix_t::ColsAtCompileTime;
593 constexpr int Cols_ = matrix_t::RowsAtCompileTime;
595 using transpose_t = Eigen::Matrix<
597 (matrix_t::IsRowMajor) ? Eigen::ColMajor : Eigen::RowMajor,
598 matrix_t::MaxColsAtCompileTime, matrix_t::MaxRowsAtCompileTime>;
600 using map_t = Eigen::Map<transpose_t, Eigen::Unaligned, Stride>;
602 return map_t((T*)A.data(), A.cols(), A.rows(), A.outerStride());
606template <
class matrix_t,
607 typename std::enable_if<(traits::is_eigen_type<matrix_t> &&
608 !matrix_t::IsVectorAtCompileTime),
610auto reshape(matrix_t& A, Eigen::Index m, Eigen::Index n)
612 using idx_t = Eigen::Index;
613 using T =
typename matrix_t::Scalar;
615 Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic,
616 matrix_t::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor>;
617 using map_t = Eigen::Map<newm_t, Eigen::Unaligned, Eigen::OuterStride<>>;
620 const idx_t size = A.size();
621 const idx_t new_size = m * n;
622 const idx_t stride = A.outerStride();
623 const bool is_contiguous =
625 (matrix_t::IsRowMajor ? (stride == A.cols() || A.rows() <= 1)
626 : (stride == A.rows() || A.cols() <= 1));
630 throw std::domain_error(
"New size is larger than current size");
631 if (A.innerStride() != 1)
632 throw std::domain_error(
633 "Reshaping is not available for matrices with inner stride "
634 "different than 1.");
637 const idx_t s = size - new_size;
638 if constexpr (matrix_t::IsRowMajor)
639 return std::make_pair(map_t((T*)A.data(), m, n, n),
640 map_t((T*)A.data() + new_size, 1, s, s));
642 return std::make_pair(map_t((T*)A.data(), m, n, m),
643 map_t((T*)A.data() + new_size, s, 1, s));
646 if (m == A.rows() || n == 0) {
647 return std::make_pair(
648 map_t((T*)A.data(), m, n, stride),
649 map_t((T*)A.data() + (matrix_t::IsRowMajor ? n : n * stride),
650 A.rows(), A.cols() - n, stride));
652 else if (n == A.cols() || m == 0) {
653 return std::make_pair(
654 map_t((T*)A.data(), m, n, stride),
655 map_t((T*)A.data() + (matrix_t::IsRowMajor ? m * stride : m),
656 A.rows() - m, A.cols(), stride));
659 throw std::domain_error(
660 "Cannot reshape to non-contiguous matrix if both the number of "
661 "rows and columns are different.");
665template <
class vector_t,
666 typename std::enable_if<(traits::is_eigen_type<vector_t> &&
667 vector_t::IsVectorAtCompileTime),
669auto reshape(vector_t& v, Eigen::Index m, Eigen::Index n)
671 using idx_t = Eigen::Index;
672 constexpr idx_t Rows_ = vector_t::IsRowMajor ? 1 : Eigen::Dynamic;
673 constexpr idx_t Cols_ = vector_t::IsRowMajor ? Eigen::Dynamic : 1;
675 using T =
typename vector_t::Scalar;
676 using newm_t = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic,
677 (vector_t::IsRowMajor) ? Eigen::RowMajor
679 using newv_t = Eigen::Matrix<T, Rows_, Cols_,
680 (vector_t::IsRowMajor) ? Eigen::RowMajor
682 using map0_t = Eigen::Map<newm_t, Eigen::Unaligned, Eigen::OuterStride<>>;
683 using map1_t = Eigen::Map<newv_t, Eigen::Unaligned, Eigen::InnerStride<>>;
686 const idx_t size = v.size();
687 const idx_t new_size = m * n;
688 const idx_t s = size - new_size;
689 const idx_t stride = v.innerStride();
690 const bool is_contiguous = (size <= 1 || stride == 1);
694 throw std::domain_error(
"New size is larger than current size");
695 if (!is_contiguous && (newm_t::IsRowMajor ? (n > 1) : (m > 1)))
696 throw std::domain_error(
697 "New sizes are not compatible with the current vector.");
700 return std::make_pair(
701 map0_t((T*)v.data(), m, n, (newm_t::IsRowMajor) ? n : m),
702 map1_t((T*)v.data() + new_size, s, Eigen::InnerStride<>(1)));
705 return std::make_pair(map0_t((T*)v.data(), m, n, stride),
706 map1_t((T*)v.data() + new_size * stride, s,
707 Eigen::InnerStride<>(stride)));
712template <
class matrix_t,
713 typename std::enable_if<(traits::is_eigen_type<matrix_t> &&
714 !matrix_t::IsVectorAtCompileTime),
716auto reshape(matrix_t& A, Eigen::Index n)
718 using idx_t = Eigen::Index;
719 constexpr idx_t Rows_ = matrix_t::IsRowMajor ? 1 : Eigen::Dynamic;
720 constexpr idx_t Cols_ = matrix_t::IsRowMajor ? Eigen::Dynamic : 1;
722 using T =
typename matrix_t::Scalar;
724 Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic,
725 matrix_t::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor>;
727 Eigen::Matrix<T, Rows_, Cols_,
728 matrix_t::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor>;
729 using map0_t = Eigen::Map<newv_t, Eigen::Unaligned, Eigen::InnerStride<>>;
730 using map1_t = Eigen::Map<newm_t, Eigen::Unaligned, Eigen::OuterStride<>>;
733 const idx_t size = A.size();
734 const idx_t stride = A.outerStride();
735 const bool is_contiguous =
737 (matrix_t::IsRowMajor ? (stride == A.cols() || A.rows() <= 1)
738 : (stride == A.rows() || A.cols() <= 1));
742 throw std::domain_error(
"New size is larger than current size");
743 if (A.innerStride() != 1)
744 throw std::domain_error(
745 "Reshaping is not available for matrices with inner stride "
746 "different than 1.");
749 const idx_t s = size - n;
750 return std::make_pair(map0_t((T*)A.data(), n, Eigen::InnerStride<>(1)),
751 (matrix_t::IsRowMajor)
752 ? map1_t((T*)A.data() + n, 1, s, s)
753 : map1_t((T*)A.data() + n, s, 1, s));
757 return std::make_pair(
758 map0_t((T*)A.data(), 0, Eigen::InnerStride<>(1)),
759 map1_t((T*)A.data(), A.rows(), A.cols(), stride));
761 else if (n == A.rows()) {
762 if constexpr (matrix_t::IsRowMajor)
763 return std::make_pair(
764 map0_t((T*)A.data(), n, Eigen::InnerStride<>(stride)),
765 map1_t((T*)A.data() + 1, A.rows(), A.cols() - 1, stride));
767 return std::make_pair(
768 map0_t((T*)A.data(), n, Eigen::InnerStride<>(1)),
769 map1_t((T*)A.data() + stride, A.rows(), A.cols() - 1,
772 else if (n == A.cols()) {
773 if constexpr (matrix_t::IsRowMajor)
774 return std::make_pair(
775 map0_t((T*)A.data(), n, Eigen::InnerStride<>(1)),
776 map1_t((T*)A.data() + stride, A.rows() - 1, A.cols(),
779 return std::make_pair(
780 map0_t((T*)A.data(), n, Eigen::InnerStride<>(stride)),
781 map1_t((T*)A.data() + 1, A.rows() - 1, A.cols(), stride));
784 throw std::domain_error(
785 "Cannot reshape to non-contiguous matrix into a vector if the "
786 "new size is different from the number of rows and columns.");
790template <
class vector_t,
791 typename std::enable_if<(traits::is_eigen_type<vector_t> &&
792 vector_t::IsVectorAtCompileTime),
794auto reshape(vector_t& v, Eigen::Index n)
796 using idx_t = Eigen::Index;
797 constexpr idx_t Rows_ = vector_t::IsRowMajor ? 1 : Eigen::Dynamic;
798 constexpr idx_t Cols_ = vector_t::IsRowMajor ? Eigen::Dynamic : 1;
800 using T =
typename vector_t::Scalar;
802 Eigen::Matrix<T, Rows_, Cols_,
803 vector_t::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor>;
804 using map_t = Eigen::Map<newv_t, Eigen::Unaligned, Eigen::InnerStride<>>;
807 const idx_t stride = v.innerStride();
811 throw std::domain_error(
"New size is larger than current size");
813 return std::make_pair(map_t((T*)v.data(), n, Eigen::InnerStride<>(stride)),
814 map_t((T*)v.data() + n * stride, v.size() - n,
815 Eigen::InnerStride<>(stride)));
823 template <
typename T>
824 constexpr bool cast_to_eigen_type = is_eigen_type<T> || is_stdvector_type<T>
825#ifdef TLAPACK_LEGACYARRAY_HH
828#ifdef TLAPACK_MDSPAN_HH
835 template <
typename matrixA_t,
typename matrixB_t>
839 typename std::enable_if<
848 static constexpr int L =
849 ((LA == Layout::RowMajor) && (LB == Layout::RowMajor))
853 using type = Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, L>;
858 template <
typename vecA_t,
typename vecB_t>
863 ((is_eigen_type<vecA_t> ||
864 is_eigen_type<vecB_t>)&&cast_to_eigen_type<vecA_t> &&
865 cast_to_eigen_type<vecB_t>),
868 using type = Eigen::Matrix<T, Eigen::Dynamic, 1>;
871#if !defined(TLAPACK_MDSPAN_HH) && !defined(TLAPACK_LEGACYARRAY_HH)
872 template <
class vecA_t,
class vecB_t>
876 std::enable_if_t<traits::is_stdvector_type<vecA_t> &&
877 traits::is_stdvector_type<vecB_t>,
880 using type = Eigen::Matrix<T, Eigen::Dynamic, 1, Eigen::ColMajor>;
883 template <
class vecA_t,
class vecB_t>
887 std::enable_if_t<traits::is_stdvector_type<vecA_t> &&
888 traits::is_stdvector_type<vecB_t>,
891 using type = Eigen::Matrix<T, Eigen::Dynamic, 1, Eigen::ColMajor>;
900template <
class T,
int Rows,
int Cols,
int Options,
int MaxRows,
int MaxCols>
901constexpr auto legacy_matrix(
902 const Eigen::Matrix<T, Rows, Cols, Options, MaxRows, MaxCols>&
A)
noexcept
904 using matrix_t = Eigen::Matrix<T, Rows, Cols, Options, MaxRows, MaxCols>;
905 using idx_t = Eigen::Index;
907 if constexpr (matrix_t::IsVectorAtCompileTime)
909 (T*)
A.data(),
A.innerStride()};
911 constexpr Layout L = layout<matrix_t>;
912 return legacy::Matrix<T, idx_t>{L, A.rows(), A.cols(), (T*)A.data(),
917template <
class Derived>
918constexpr auto legacy_matrix(
919 const Eigen::MapBase<Derived, Eigen::ReadOnlyAccessors>& A)
noexcept
921 using T =
typename Derived::Scalar;
922 using idx_t = Eigen::Index;
924 if constexpr (Derived::IsVectorAtCompileTime) {
925 assert(A.outerStride() == 1 ||
926 (A.innerStride() == 1 &&
927 traits::internal::is_eigen_block<Derived>));
928 return legacy::Matrix<T, idx_t>{Layout::ColMajor, 1, A.size(),
929 (T*)A.data(), A.innerStride()};
932 assert(A.innerStride() == 1);
933 constexpr Layout L = layout<Derived>;
934 return legacy::Matrix<T, idx_t>{L, A.rows(), A.cols(), (T*)A.data(),
939template <
class T,
int Rows,
int Cols,
int Options,
int MaxRows,
int MaxCols>
940constexpr auto legacy_vector(
941 const Eigen::Matrix<T, Rows, Cols, Options, MaxRows, MaxCols>& A)
noexcept
943 using matrix_t = Eigen::Matrix<T, Rows, Cols, Options, MaxRows, MaxCols>;
944 using idx_t = Eigen::Index;
946 if constexpr (matrix_t::IsVectorAtCompileTime)
947 return legacy::Vector<T, idx_t>{A.size(), (T*)A.data(),
950 assert(A.rows() == 1 || A.cols() == 1);
951 return legacy::Vector<T, idx_t>{A.size(), (T*)A.data(),
956template <
class Derived>
957constexpr auto legacy_vector(
958 const Eigen::MapBase<Derived, Eigen::ReadOnlyAccessors>& A)
noexcept
960 using T =
typename Derived::Scalar;
961 using idx_t = Eigen::Index;
963 if constexpr (Derived::IsVectorAtCompileTime) {
964 assert(A.outerStride() == 1 ||
965 (A.innerStride() == 1 &&
966 traits::internal::is_eigen_block<Derived>));
967 return legacy::Vector<T, idx_t>{A.size(), (T*)A.data(),
971 assert(A.innerStride() == 1);
972 assert(A.rows() == 1 || A.cols() == 1);
973 return legacy::Vector<T, idx_t>{A.size(), (T*)A.data(),
Layout
Definition types.hpp:29
@ ColMajor
Column-major layout.
constexpr bool is_eigen_type
True if T is derived from Eigen::DenseBase.
Definition eigen.hpp:69
constexpr bool is_eigen_matrix
True if T is derived from Eigen::EigenMatrix<T>
Definition eigen.hpp:50
constexpr auto diag(T &A, int diagIdx=0) noexcept
Get the Diagonal of an Eigen Matrix.
Definition eigen.hpp:576
constexpr bool is_eigen_block
True if T is derived from Eigen::EigenBlock.
Definition eigen.hpp:56
constexpr bool is_eigen_dense
True if T is derived from Eigen::EigenDense<T>
Definition eigen.hpp:44
typename traits::real_type_traits< Types..., int >::type real_type
The common real type of the list of types.
Definition scalar_type_traits.hpp:113
Describes a row- or column-major matrix.
Definition types.hpp:448
Functor for data creation.
Definition arrayTraits.hpp:89
constexpr auto operator()(std::vector< T > &v, idx_t m, idx_t n=1) const
Creates a m-by-n matrix with entries of type T.
Definition arrayTraits.hpp:105
Functor for data creation with static size.
Definition arrayTraits.hpp:141
constexpr auto operator()(T *v) const
Creates a m-by-n matrix or, if n == -1, a vector of size m.
Definition arrayTraits.hpp:157
Complex type traits for the list of types Types.
Definition scalar_type_traits.hpp:145
Trait to determine the layout of a given data structure.
Definition arrayTraits.hpp:75
Matrix type deduction.
Definition arrayTraits.hpp:176
Real type traits for the list of types Types.
Definition scalar_type_traits.hpp:71
Vector type deduction.
Definition arrayTraits.hpp:203