CompleteOrthogonalDecomposition.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2016 Rasmus Munk Larsen <rmlarsen@google.com>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
11 #define EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 template <typename MatrixType_, typename PermutationIndex_>
20 struct traits<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_>> : traits<MatrixType_> {
21  typedef MatrixXpr XprKind;
23  typedef PermutationIndex_ PermutationIndex;
24  enum { Flags = 0 };
25 };
26 
27 } // end namespace internal
28 
52 template <typename MatrixType_, typename PermutationIndex_>
54  : public SolverBase<CompleteOrthogonalDecomposition<MatrixType_, PermutationIndex_>> {
55  public:
56  typedef MatrixType_ MatrixType;
58 
59  template <typename Derived>
61  typedef PermutationIndex_ PermutationIndex;
63  enum {
64  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
65  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
66  };
74  typedef typename MatrixType::PlainObject PlainObject;
75 
76  public:
85 
93  : m_cpqr(rows, cols), m_zCoeffs((std::min)(rows, cols)), m_temp(cols) {}
94 
111  template <typename InputType>
113  : m_cpqr(matrix.rows(), matrix.cols()),
114  m_zCoeffs((std::min)(matrix.rows(), matrix.cols())),
115  m_temp(matrix.cols()) {
116  compute(matrix.derived());
117  }
118 
126  template <typename InputType>
129  computeInPlace();
130  }
131 
132 #ifdef EIGEN_PARSED_BY_DOXYGEN
142  template <typename Rhs>
144 #endif
145 
148 
151  MatrixType matrixZ() const {
152  MatrixType Z = MatrixType::Identity(m_cpqr.cols(), m_cpqr.cols());
153  applyZOnTheLeftInPlace<false>(Z);
154  return Z;
155  }
156 
160  const MatrixType& matrixQTZ() const { return m_cpqr.matrixQR(); }
161 
173  const MatrixType& matrixT() const { return m_cpqr.matrixQR(); }
174 
175  template <typename InputType>
177  // Compute the column pivoted QR factorization A P = Q R.
179  computeInPlace();
180  return *this;
181  }
182 
185 
199  typename MatrixType::Scalar determinant() const;
200 
214  typename MatrixType::RealScalar absDeterminant() const;
215 
229  typename MatrixType::RealScalar logAbsDeterminant() const;
230 
244  typename MatrixType::Scalar signDeterminant() const;
245 
253  inline Index rank() const { return m_cpqr.rank(); }
254 
262  inline Index dimensionOfKernel() const { return m_cpqr.dimensionOfKernel(); }
263 
271  inline bool isInjective() const { return m_cpqr.isInjective(); }
272 
280  inline bool isSurjective() const { return m_cpqr.isSurjective(); }
281 
289  inline bool isInvertible() const { return m_cpqr.isInvertible(); }
290 
297  eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
299  }
300 
301  inline Index rows() const { return m_cpqr.rows(); }
302  inline Index cols() const { return m_cpqr.cols(); }
303 
309  inline const HCoeffsType& hCoeffs() const { return m_cpqr.hCoeffs(); }
310 
316  const HCoeffsType& zCoeffs() const { return m_zCoeffs; }
317 
339  return *this;
340  }
341 
352  return *this;
353  }
354 
359  RealScalar threshold() const { return m_cpqr.threshold(); }
360 
368  inline Index nonzeroPivots() const { return m_cpqr.nonzeroPivots(); }
369 
373  inline RealScalar maxPivot() const { return m_cpqr.maxPivot(); }
374 
384  eigen_assert(m_cpqr.m_isInitialized && "Decomposition is not initialized.");
385  return Success;
386  }
387 
388 #ifndef EIGEN_PARSED_BY_DOXYGEN
389  template <typename RhsType, typename DstType>
390  void _solve_impl(const RhsType& rhs, DstType& dst) const;
391 
392  template <bool Conjugate, typename RhsType, typename DstType>
393  void _solve_impl_transposed(const RhsType& rhs, DstType& dst) const;
394 #endif
395 
396  protected:
398 
399  template <bool Transpose_, typename Rhs>
400  void _check_solve_assertion(const Rhs& b) const {
402  eigen_assert(m_cpqr.m_isInitialized && "CompleteOrthogonalDecomposition is not initialized.");
403  eigen_assert((Transpose_ ? derived().cols() : derived().rows()) == b.rows() &&
404  "CompleteOrthogonalDecomposition::solve(): invalid number of rows of the right hand side matrix b");
405  }
406 
407  void computeInPlace();
408 
413  template <bool Conjugate, typename Rhs>
414  void applyZOnTheLeftInPlace(Rhs& rhs) const;
415 
418  template <typename Rhs>
419  void applyZAdjointOnTheLeftInPlace(Rhs& rhs) const;
420 
424 };
425 
426 template <typename MatrixType, typename PermutationIndex>
428  return m_cpqr.determinant();
429 }
430 
431 template <typename MatrixType, typename PermutationIndex>
433  return m_cpqr.absDeterminant();
434 }
435 
436 template <typename MatrixType, typename PermutationIndex>
438  const {
439  return m_cpqr.logAbsDeterminant();
440 }
441 
442 template <typename MatrixType, typename PermutationIndex>
444  return m_cpqr.signDeterminant();
445 }
446 
454 template <typename MatrixType, typename PermutationIndex>
457 
458  const Index rank = m_cpqr.rank();
459  const Index cols = m_cpqr.cols();
460  const Index rows = m_cpqr.rows();
461  m_zCoeffs.resize((std::min)(rows, cols));
462  m_temp.resize(cols);
463 
464  if (rank < cols) {
465  // We have reduced the (permuted) matrix to the form
466  // [R11 R12]
467  // [ 0 R22]
468  // where R11 is r-by-r (r = rank) upper triangular, R12 is
469  // r-by-(n-r), and R22 is empty or the norm of R22 is negligible.
470  // We now compute the complete orthogonal decomposition by applying
471  // Householder transformations from the right to the upper trapezoidal
472  // matrix X = [R11 R12] to zero out R12 and obtain the factorization
473  // [R11 R12] = [T11 0] * Z, where T11 is r-by-r upper triangular and
474  // Z = Z(0) * Z(1) ... Z(r-1) is an n-by-n orthogonal matrix.
475  // We store the data representing Z in R12 and m_zCoeffs.
476  for (Index k = rank - 1; k >= 0; --k) {
477  if (k != rank - 1) {
478  // Given the API for Householder reflectors, it is more convenient if
479  // we swap the leading parts of columns k and r-1 (zero-based) to form
480  // the matrix X_k = [X(0:k, k), X(0:k, r:n)]
481  m_cpqr.m_qr.col(k).head(k + 1).swap(m_cpqr.m_qr.col(rank - 1).head(k + 1));
482  }
483  // Construct Householder reflector Z(k) to zero out the last row of X_k,
484  // i.e. choose Z(k) such that
485  // [X(k, k), X(k, r:n)] * Z(k) = [beta, 0, .., 0].
487  m_cpqr.m_qr.row(k).tail(cols - rank + 1).makeHouseholderInPlace(m_zCoeffs(k), beta);
488  m_cpqr.m_qr(k, rank - 1) = beta;
489  if (k > 0) {
490  // Apply Z(k) to the first k rows of X_k
491  m_cpqr.m_qr.topRightCorner(k, cols - rank + 1)
492  .applyHouseholderOnTheRight(m_cpqr.m_qr.row(k).tail(cols - rank).adjoint(), m_zCoeffs(k), &m_temp(0));
493  }
494  if (k != rank - 1) {
495  // Swap X(0:k,k) back to its proper location.
496  m_cpqr.m_qr.col(k).head(k + 1).swap(m_cpqr.m_qr.col(rank - 1).head(k + 1));
497  }
498  }
499  }
500 }
501 
502 template <typename MatrixType, typename PermutationIndex>
503 template <bool Conjugate, typename Rhs>
505  const Index cols = this->cols();
506  const Index nrhs = rhs.cols();
507  const Index rank = this->rank();
509  for (Index k = rank - 1; k >= 0; --k) {
510  if (k != rank - 1) {
511  rhs.row(k).swap(rhs.row(rank - 1));
512  }
513  rhs.middleRows(rank - 1, cols - rank + 1)
514  .applyHouseholderOnTheLeft(matrixQTZ().row(k).tail(cols - rank).transpose().template conjugateIf<!Conjugate>(),
515  zCoeffs().template conjugateIf<Conjugate>()(k), &temp(0));
516  if (k != rank - 1) {
517  rhs.row(k).swap(rhs.row(rank - 1));
518  }
519  }
520 }
521 
522 template <typename MatrixType, typename PermutationIndex>
523 template <typename Rhs>
525  const Index cols = this->cols();
526  const Index nrhs = rhs.cols();
527  const Index rank = this->rank();
529  for (Index k = 0; k < rank; ++k) {
530  if (k != rank - 1) {
531  rhs.row(k).swap(rhs.row(rank - 1));
532  }
533  rhs.middleRows(rank - 1, cols - rank + 1)
534  .applyHouseholderOnTheLeft(matrixQTZ().row(k).tail(cols - rank).adjoint(), zCoeffs()(k), &temp(0));
535  if (k != rank - 1) {
536  rhs.row(k).swap(rhs.row(rank - 1));
537  }
538  }
539 }
540 
541 #ifndef EIGEN_PARSED_BY_DOXYGEN
542 template <typename MatrixType_, typename PermutationIndex_>
543 template <typename RhsType, typename DstType>
545  DstType& dst) const {
546  const Index rank = this->rank();
547  if (rank == 0) {
548  dst.setZero();
549  return;
550  }
551 
552  // Compute c = Q^* * rhs
553  typename RhsType::PlainObject c(rhs);
554  c.applyOnTheLeft(matrixQ().setLength(rank).adjoint());
555 
556  // Solve T z = c(1:rank, :)
557  dst.topRows(rank) = matrixT().topLeftCorner(rank, rank).template triangularView<Upper>().solve(c.topRows(rank));
558 
559  const Index cols = this->cols();
560  if (rank < cols) {
561  // Compute y = Z^* * [ z ]
562  // [ 0 ]
563  dst.bottomRows(cols - rank).setZero();
564  applyZAdjointOnTheLeftInPlace(dst);
565  }
566 
567  // Undo permutation to get x = P^{-1} * y.
568  dst = colsPermutation() * dst;
569 }
570 
571 template <typename MatrixType_, typename PermutationIndex_>
572 template <bool Conjugate, typename RhsType, typename DstType>
574  DstType& dst) const {
575  const Index rank = this->rank();
576 
577  if (rank == 0) {
578  dst.setZero();
579  return;
580  }
581 
582  typename RhsType::PlainObject c(colsPermutation().transpose() * rhs);
583 
584  if (rank < cols()) {
585  applyZOnTheLeftInPlace<!Conjugate>(c);
586  }
587 
588  matrixT()
589  .topLeftCorner(rank, rank)
590  .template triangularView<Upper>()
591  .transpose()
592  .template conjugateIf<Conjugate>()
593  .solveInPlace(c.topRows(rank));
594 
595  dst.topRows(rank) = c.topRows(rank);
596  dst.bottomRows(rows() - rank).setZero();
597 
598  dst.applyOnTheLeft(householderQ().setLength(rank).template conjugateIf<!Conjugate>());
599 }
600 #endif
601 
602 namespace internal {
603 
604 template <typename MatrixType, typename PermutationIndex>
606  : traits<typename Transpose<typename MatrixType::PlainObject>::PlainObject> {
607  enum { Flags = 0 };
608 };
609 
610 template <typename DstXprType, typename MatrixType, typename PermutationIndex>
611 struct Assignment<DstXprType, Inverse<CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>>,
612  internal::assign_op<typename DstXprType::Scalar,
613  typename CompleteOrthogonalDecomposition<MatrixType, PermutationIndex>::Scalar>,
614  Dense2Dense> {
617  static void run(DstXprType& dst, const SrcXprType& src,
619  typedef Matrix<typename CodType::Scalar, CodType::RowsAtCompileTime, CodType::RowsAtCompileTime, 0,
620  CodType::MaxRowsAtCompileTime, CodType::MaxRowsAtCompileTime>
621  IdentityMatrixType;
622  dst = src.nestedExpression().solve(IdentityMatrixType::Identity(src.cols(), src.cols()));
623  }
624 };
625 
626 } // end namespace internal
627 
629 template <typename MatrixType, typename PermutationIndex>
632  return m_cpqr.householderQ();
633 }
634 
639 template <typename Derived>
640 template <typename PermutationIndex>
644 }
645 
646 } // end namespace Eigen
647 
648 #endif // EIGEN_COMPLETEORTHOGONALDECOMPOSITION_H
#define EIGEN_GENERIC_PUBLIC_INTERFACE(Derived)
Definition: Macros.h:1149
#define EIGEN_ONLY_USED_FOR_DEBUG(x)
Definition: Macros.h:922
#define eigen_assert(x)
Definition: Macros.h:910
m row(1)
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition: StaticAssert.h:74
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
void adjoint(const MatrixType &m)
Definition: adjoint.cpp:85
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
RealScalar maxPivot() const
Definition: ColPivHouseholderQR.h:398
Index rows() const
Definition: ColPivHouseholderQR.h:326
Index rank() const
Definition: ColPivHouseholderQR.h:261
bool isSurjective() const
Definition: ColPivHouseholderQR.h:300
Index nonzeroPivots() const
Definition: ColPivHouseholderQR.h:390
bool isInjective() const
Definition: ColPivHouseholderQR.h:288
ColPivHouseholderQR & setThreshold(const RealScalar &threshold)
Definition: ColPivHouseholderQR.h:352
Index cols() const
Definition: ColPivHouseholderQR.h:327
const HCoeffsType & hCoeffs() const
Definition: ColPivHouseholderQR.h:333
bool isInvertible() const
Definition: ColPivHouseholderQR.h:311
Index dimensionOfKernel() const
Definition: ColPivHouseholderQR.h:276
bool m_isInitialized
Definition: ColPivHouseholderQR.h:433
const PermutationType & colsPermutation() const
Definition: ColPivHouseholderQR.h:192
HouseholderSequenceType householderQ() const
Definition: ColPivHouseholderQR.h:656
ColPivHouseholderQR & compute(const EigenBase< InputType > &matrix)
RealScalar threshold() const
Definition: ColPivHouseholderQR.h:375
const MatrixType & matrixQR() const
Definition: ColPivHouseholderQR.h:169
Complete orthogonal decomposition (COD) of a matrix.
Definition: CompleteOrthogonalDecomposition.h:54
internal::plain_diag_type< MatrixType >::type HCoeffsType
Definition: CompleteOrthogonalDecomposition.h:67
MatrixType matrixZ() const
Definition: CompleteOrthogonalDecomposition.h:151
MatrixType_ MatrixType
Definition: CompleteOrthogonalDecomposition.h:56
void applyZAdjointOnTheLeftInPlace(Rhs &rhs) const
Definition: CompleteOrthogonalDecomposition.h:524
CompleteOrthogonalDecomposition()
Default Constructor.
Definition: CompleteOrthogonalDecomposition.h:84
bool isSurjective() const
Definition: CompleteOrthogonalDecomposition.h:280
MatrixType::Scalar signDeterminant() const
Definition: CompleteOrthogonalDecomposition.h:443
Index dimensionOfKernel() const
Definition: CompleteOrthogonalDecomposition.h:262
HouseholderSequenceType matrixQ(void) const
Definition: CompleteOrthogonalDecomposition.h:147
CompleteOrthogonalDecomposition & compute(const EigenBase< InputType > &matrix)
Definition: CompleteOrthogonalDecomposition.h:176
PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex > PermutationType
Definition: CompleteOrthogonalDecomposition.h:68
CompleteOrthogonalDecomposition(const EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition: CompleteOrthogonalDecomposition.h:112
CompleteOrthogonalDecomposition(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition: CompleteOrthogonalDecomposition.h:92
CompleteOrthogonalDecomposition & setThreshold(const RealScalar &threshold)
Definition: CompleteOrthogonalDecomposition.h:337
MatrixType::RealScalar absDeterminant() const
Definition: CompleteOrthogonalDecomposition.h:432
const HCoeffsType & zCoeffs() const
Definition: CompleteOrthogonalDecomposition.h:316
void _solve_impl(const RhsType &rhs, DstType &dst) const
Definition: CompleteOrthogonalDecomposition.h:544
Index cols() const
Definition: CompleteOrthogonalDecomposition.h:302
HouseholderSequenceType householderQ(void) const
Definition: CompleteOrthogonalDecomposition.h:631
MatrixType::RealScalar logAbsDeterminant() const
Definition: CompleteOrthogonalDecomposition.h:437
RealScalar threshold() const
Definition: CompleteOrthogonalDecomposition.h:359
ColPivHouseholderQR< MatrixType, PermutationIndex > m_cpqr
Definition: CompleteOrthogonalDecomposition.h:421
RowVectorType m_temp
Definition: CompleteOrthogonalDecomposition.h:423
ComputationInfo info() const
Reports whether the complete orthogonal decomposition was successful.
Definition: CompleteOrthogonalDecomposition.h:383
void _solve_impl_transposed(const RhsType &rhs, DstType &dst) const
Definition: CompleteOrthogonalDecomposition.h:573
const PermutationType & colsPermutation() const
Definition: CompleteOrthogonalDecomposition.h:184
internal::plain_row_type< MatrixType >::type RowVectorType
Definition: CompleteOrthogonalDecomposition.h:70
SolverBase< CompleteOrthogonalDecomposition > Base
Definition: CompleteOrthogonalDecomposition.h:57
internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
Definition: CompleteOrthogonalDecomposition.h:69
bool isInvertible() const
Definition: CompleteOrthogonalDecomposition.h:289
HouseholderSequence< MatrixType, internal::remove_all_t< typename HCoeffsType::ConjugateReturnType > > HouseholderSequenceType
Definition: CompleteOrthogonalDecomposition.h:73
const HCoeffsType & hCoeffs() const
Definition: CompleteOrthogonalDecomposition.h:309
void _check_solve_assertion(const Rhs &b) const
Definition: CompleteOrthogonalDecomposition.h:400
MatrixType::PlainObject PlainObject
Definition: CompleteOrthogonalDecomposition.h:74
internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
Definition: CompleteOrthogonalDecomposition.h:71
bool isInjective() const
Definition: CompleteOrthogonalDecomposition.h:271
Index rows() const
Definition: CompleteOrthogonalDecomposition.h:301
CompleteOrthogonalDecomposition & setThreshold(Default_t)
Definition: CompleteOrthogonalDecomposition.h:350
void applyZOnTheLeftInPlace(Rhs &rhs) const
Definition: CompleteOrthogonalDecomposition.h:504
void computeInPlace()
Definition: CompleteOrthogonalDecomposition.h:455
const MatrixType & matrixQTZ() const
Definition: CompleteOrthogonalDecomposition.h:160
Index nonzeroPivots() const
Definition: CompleteOrthogonalDecomposition.h:368
const Inverse< CompleteOrthogonalDecomposition > pseudoInverse() const
Definition: CompleteOrthogonalDecomposition.h:296
HCoeffsType m_zCoeffs
Definition: CompleteOrthogonalDecomposition.h:422
MatrixType::Scalar determinant() const
Definition: CompleteOrthogonalDecomposition.h:427
Index rank() const
Definition: CompleteOrthogonalDecomposition.h:253
RealScalar maxPivot() const
Definition: CompleteOrthogonalDecomposition.h:373
@ MaxColsAtCompileTime
Definition: CompleteOrthogonalDecomposition.h:65
@ MaxRowsAtCompileTime
Definition: CompleteOrthogonalDecomposition.h:64
const MatrixType & matrixT() const
Definition: CompleteOrthogonalDecomposition.h:173
CompleteOrthogonalDecomposition(EigenBase< InputType > &matrix)
Constructs a complete orthogonal decomposition from a given matrix.
Definition: CompleteOrthogonalDecomposition.h:127
PermutationIndex_ PermutationIndex
Definition: CompleteOrthogonalDecomposition.h:61
Sequence of Householder reflections acting on subspaces with decreasing size.
Definition: HouseholderSequence.h:117
Expression of the inverse of another expression.
Definition: Inverse.h:43
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: Inverse.h:55
EIGEN_DEVICE_FUNC const XprTypeNestedCleaned & nestedExpression() const
Definition: Inverse.h:57
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const CompleteOrthogonalDecomposition< PlainObject, PermutationIndex > completeOrthogonalDecomposition() const
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Pseudo expression representing a solving operation.
Definition: Solve.h:62
A base class for matrix decomposition and solvers.
Definition: SolverBase.h:72
constexpr EIGEN_DEVICE_FUNC CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ > & derived()
Definition: EigenBase.h:49
const Solve< CompleteOrthogonalDecomposition< MatrixType_, PermutationIndex_ >, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SolverBase.h:106
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
ComputationInfo
Definition: Constants.h:438
@ Success
Definition: Constants.h:440
#define Z
Definition: icosphere.cpp:21
Scalar beta
Definition: level2_cplx_impl.h:36
char char char int int * k
Definition: level2_impl.h:374
@ Rhs
Definition: TensorContractionMapper.h:20
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
Definition: AutoDiffScalar.h:494
Default_t
Definition: Constants.h:361
@ Default
Definition: Constants.h:361
Extend namespace for flags.
Definition: fsi_chan_precond_driver.cc:56
void transpose()
Definition: skew_symmetric_matrix3.cpp:135
int c
Definition: calibrate.py:100
Definition: Eigen_Colamd.h:49
internal::nested_eval< T, 1 >::type eval(const T &xpr)
Definition: sparse_permutations.cpp:47
Definition: EigenBase.h:33
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
Definition: Constants.h:534
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: Constants.h:525
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< typename DstXprType::Scalar, typename CodType::Scalar > &)
Definition: CompleteOrthogonalDecomposition.h:617
Definition: AssignEvaluator.h:773
Definition: AssignEvaluator.h:756
Template functor for scalar/packet assignment.
Definition: AssignmentFunctors.h:25
std::conditional_t< is_same< typename traits< ExpressionType >::XprKind, MatrixXpr >::value, MatrixDiagType, ArrayDiagType > type
Definition: XprHelper.h:797
std::conditional_t< is_same< typename traits< ExpressionType >::XprKind, MatrixXpr >::value, MatrixRowType, ArrayRowType > type
Definition: XprHelper.h:768
Definition: SolverBase.h:21
PermutationIndex_ PermutationIndex
Definition: CompleteOrthogonalDecomposition.h:23
SolverStorage StorageKind
Definition: CompleteOrthogonalDecomposition.h:22
Definition: ForwardDeclarations.h:21