Eigen  3.3.5
LDLT.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
6 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7 // Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
8 //
9 // This Source Code Form is subject to the terms of the Mozilla
10 // Public License v. 2.0. If a copy of the MPL was not distributed
11 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12 
13 #ifndef EIGEN_LDLT_H
14 #define EIGEN_LDLT_H
15 
16 namespace Eigen {
17 
18 namespace internal {
19  template<typename MatrixType, int UpLo> struct LDLT_Traits;
20 
21  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
22  enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
23 }
24 
50 template<typename _MatrixType, int _UpLo> class LDLT
51 {
52  public:
53  typedef _MatrixType MatrixType;
54  enum {
55  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
58  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
59  UpLo = _UpLo
60  };
61  typedef typename MatrixType::Scalar Scalar;
62  typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
63  typedef Eigen::Index Index;
64  typedef typename MatrixType::StorageIndex StorageIndex;
66 
69 
70  typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
71 
77  LDLT()
78  : m_matrix(),
79  m_transpositions(),
80  m_sign(internal::ZeroSign),
81  m_isInitialized(false)
82  {}
83 
90  explicit LDLT(Index size)
91  : m_matrix(size, size),
92  m_transpositions(size),
93  m_temporary(size),
94  m_sign(internal::ZeroSign),
95  m_isInitialized(false)
96  {}
97 
104  template<typename InputType>
105  explicit LDLT(const EigenBase<InputType>& matrix)
106  : m_matrix(matrix.rows(), matrix.cols()),
107  m_transpositions(matrix.rows()),
108  m_temporary(matrix.rows()),
109  m_sign(internal::ZeroSign),
110  m_isInitialized(false)
111  {
112  compute(matrix.derived());
113  }
114 
121  template<typename InputType>
122  explicit LDLT(EigenBase<InputType>& matrix)
123  : m_matrix(matrix.derived()),
124  m_transpositions(matrix.rows()),
125  m_temporary(matrix.rows()),
126  m_sign(internal::ZeroSign),
127  m_isInitialized(false)
128  {
129  compute(matrix.derived());
130  }
131 
135  void setZero()
136  {
137  m_isInitialized = false;
138  }
139 
141  inline typename Traits::MatrixU matrixU() const
142  {
143  eigen_assert(m_isInitialized && "LDLT is not initialized.");
144  return Traits::getU(m_matrix);
145  }
146 
148  inline typename Traits::MatrixL matrixL() const
149  {
150  eigen_assert(m_isInitialized && "LDLT is not initialized.");
151  return Traits::getL(m_matrix);
152  }
153 
156  inline const TranspositionType& transpositionsP() const
157  {
158  eigen_assert(m_isInitialized && "LDLT is not initialized.");
159  return m_transpositions;
160  }
161 
164  {
165  eigen_assert(m_isInitialized && "LDLT is not initialized.");
166  return m_matrix.diagonal();
167  }
168 
170  inline bool isPositive() const
171  {
172  eigen_assert(m_isInitialized && "LDLT is not initialized.");
173  return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
174  }
175 
177  inline bool isNegative(void) const
178  {
179  eigen_assert(m_isInitialized && "LDLT is not initialized.");
180  return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
181  }
182 
198  template<typename Rhs>
199  inline const Solve<LDLT, Rhs>
200  solve(const MatrixBase<Rhs>& b) const
201  {
202  eigen_assert(m_isInitialized && "LDLT is not initialized.");
203  eigen_assert(m_matrix.rows()==b.rows()
204  && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
205  return Solve<LDLT, Rhs>(*this, b.derived());
206  }
207 
208  template<typename Derived>
209  bool solveInPlace(MatrixBase<Derived> &bAndX) const;
210 
211  template<typename InputType>
212  LDLT& compute(const EigenBase<InputType>& matrix);
213 
217  RealScalar rcond() const
218  {
219  eigen_assert(m_isInitialized && "LDLT is not initialized.");
220  return internal::rcond_estimate_helper(m_l1_norm, *this);
221  }
222 
223  template <typename Derived>
224  LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
225 
230  inline const MatrixType& matrixLDLT() const
231  {
232  eigen_assert(m_isInitialized && "LDLT is not initialized.");
233  return m_matrix;
234  }
235 
236  MatrixType reconstructedMatrix() const;
237 
243  const LDLT& adjoint() const { return *this; };
244 
245  inline Index rows() const { return m_matrix.rows(); }
246  inline Index cols() const { return m_matrix.cols(); }
247 
254  {
255  eigen_assert(m_isInitialized && "LDLT is not initialized.");
256  return m_info;
257  }
258 
259  #ifndef EIGEN_PARSED_BY_DOXYGEN
260  template<typename RhsType, typename DstType>
261  EIGEN_DEVICE_FUNC
262  void _solve_impl(const RhsType &rhs, DstType &dst) const;
263  #endif
264 
265  protected:
266 
267  static void check_template_parameters()
268  {
269  EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
270  }
271 
278  MatrixType m_matrix;
279  RealScalar m_l1_norm;
280  TranspositionType m_transpositions;
281  TmpMatrixType m_temporary;
282  internal::SignMatrix m_sign;
283  bool m_isInitialized;
284  ComputationInfo m_info;
285 };
286 
287 namespace internal {
288 
289 template<int UpLo> struct ldlt_inplace;
290 
291 template<> struct ldlt_inplace<Lower>
292 {
293  template<typename MatrixType, typename TranspositionType, typename Workspace>
294  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
295  {
296  using std::abs;
297  typedef typename MatrixType::Scalar Scalar;
298  typedef typename MatrixType::RealScalar RealScalar;
299  typedef typename TranspositionType::StorageIndex IndexType;
300  eigen_assert(mat.rows()==mat.cols());
301  const Index size = mat.rows();
302  bool found_zero_pivot = false;
303  bool ret = true;
304 
305  if (size <= 1)
306  {
307  transpositions.setIdentity();
308  if (numext::real(mat.coeff(0,0)) > static_cast<RealScalar>(0) ) sign = PositiveSemiDef;
309  else if (numext::real(mat.coeff(0,0)) < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
310  else sign = ZeroSign;
311  return true;
312  }
313 
314  for (Index k = 0; k < size; ++k)
315  {
316  // Find largest diagonal element
317  Index index_of_biggest_in_corner;
318  mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
319  index_of_biggest_in_corner += k;
320 
321  transpositions.coeffRef(k) = IndexType(index_of_biggest_in_corner);
322  if(k != index_of_biggest_in_corner)
323  {
324  // apply the transposition while taking care to consider only
325  // the lower triangular part
326  Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
327  mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
328  mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
329  std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
330  for(Index i=k+1;i<index_of_biggest_in_corner;++i)
331  {
332  Scalar tmp = mat.coeffRef(i,k);
333  mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
334  mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
335  }
336  if(NumTraits<Scalar>::IsComplex)
337  mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
338  }
339 
340  // partition the matrix:
341  // A00 | - | -
342  // lu = A10 | A11 | -
343  // A20 | A21 | A22
344  Index rs = size - k - 1;
345  Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
346  Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
347  Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
348 
349  if(k>0)
350  {
351  temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
352  mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
353  if(rs>0)
354  A21.noalias() -= A20 * temp.head(k);
355  }
356 
357  // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
358  // was smaller than the cutoff value. However, since LDLT is not rank-revealing
359  // we should only make sure that we do not introduce INF or NaN values.
360  // Remark that LAPACK also uses 0 as the cutoff value.
361  RealScalar realAkk = numext::real(mat.coeffRef(k,k));
362  bool pivot_is_valid = (abs(realAkk) > RealScalar(0));
363 
364  if(k==0 && !pivot_is_valid)
365  {
366  // The entire diagonal is zero, there is nothing more to do
367  // except filling the transpositions, and checking whether the matrix is zero.
368  sign = ZeroSign;
369  for(Index j = 0; j<size; ++j)
370  {
371  transpositions.coeffRef(j) = IndexType(j);
372  ret = ret && (mat.col(j).tail(size-j-1).array()==Scalar(0)).all();
373  }
374  return ret;
375  }
376 
377  if((rs>0) && pivot_is_valid)
378  A21 /= realAkk;
379  else if(rs>0)
380  ret = ret && (A21.array()==Scalar(0)).all();
381 
382  if(found_zero_pivot && pivot_is_valid) ret = false; // factorization failed
383  else if(!pivot_is_valid) found_zero_pivot = true;
384 
385  if (sign == PositiveSemiDef) {
386  if (realAkk < static_cast<RealScalar>(0)) sign = Indefinite;
387  } else if (sign == NegativeSemiDef) {
388  if (realAkk > static_cast<RealScalar>(0)) sign = Indefinite;
389  } else if (sign == ZeroSign) {
390  if (realAkk > static_cast<RealScalar>(0)) sign = PositiveSemiDef;
391  else if (realAkk < static_cast<RealScalar>(0)) sign = NegativeSemiDef;
392  }
393  }
394 
395  return ret;
396  }
397 
398  // Reference for the algorithm: Davis and Hager, "Multiple Rank
399  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
400  // Trivial rearrangements of their computations (Timothy E. Holy)
401  // allow their algorithm to work for rank-1 updates even if the
402  // original matrix is not of full rank.
403  // Here only rank-1 updates are implemented, to reduce the
404  // requirement for intermediate storage and improve accuracy
405  template<typename MatrixType, typename WDerived>
406  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
407  {
408  using numext::isfinite;
409  typedef typename MatrixType::Scalar Scalar;
410  typedef typename MatrixType::RealScalar RealScalar;
411 
412  const Index size = mat.rows();
413  eigen_assert(mat.cols() == size && w.size()==size);
414 
415  RealScalar alpha = 1;
416 
417  // Apply the update
418  for (Index j = 0; j < size; j++)
419  {
420  // Check for termination due to an original decomposition of low-rank
421  if (!(isfinite)(alpha))
422  break;
423 
424  // Update the diagonal terms
425  RealScalar dj = numext::real(mat.coeff(j,j));
426  Scalar wj = w.coeff(j);
427  RealScalar swj2 = sigma*numext::abs2(wj);
428  RealScalar gamma = dj*alpha + swj2;
429 
430  mat.coeffRef(j,j) += swj2/alpha;
431  alpha += swj2/dj;
432 
433 
434  // Update the terms of L
435  Index rs = size-j-1;
436  w.tail(rs) -= wj * mat.col(j).tail(rs);
437  if(gamma != 0)
438  mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
439  }
440  return true;
441  }
442 
443  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
444  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
445  {
446  // Apply the permutation to the input w
447  tmp = transpositions * w;
448 
449  return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
450  }
451 };
452 
453 template<> struct ldlt_inplace<Upper>
454 {
455  template<typename MatrixType, typename TranspositionType, typename Workspace>
456  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
457  {
458  Transpose<MatrixType> matt(mat);
459  return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
460  }
461 
462  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
463  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
464  {
465  Transpose<MatrixType> matt(mat);
466  return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
467  }
468 };
469 
470 template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
471 {
472  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
473  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
474  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m); }
475  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m.adjoint()); }
476 };
477 
478 template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
479 {
480  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
481  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
482  static inline MatrixL getL(const MatrixType& m) { return MatrixL(m.adjoint()); }
483  static inline MatrixU getU(const MatrixType& m) { return MatrixU(m); }
484 };
485 
486 } // end namespace internal
487 
490 template<typename MatrixType, int _UpLo>
491 template<typename InputType>
493 {
494  check_template_parameters();
495 
496  eigen_assert(a.rows()==a.cols());
497  const Index size = a.rows();
498 
499  m_matrix = a.derived();
500 
501  // Compute matrix L1 norm = max abs column sum.
502  m_l1_norm = RealScalar(0);
503  // TODO move this code to SelfAdjointView
504  for (Index col = 0; col < size; ++col) {
505  RealScalar abs_col_sum;
506  if (_UpLo == Lower)
507  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
508  else
509  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
510  if (abs_col_sum > m_l1_norm)
511  m_l1_norm = abs_col_sum;
512  }
513 
514  m_transpositions.resize(size);
515  m_isInitialized = false;
516  m_temporary.resize(size);
517  m_sign = internal::ZeroSign;
518 
519  m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
520 
521  m_isInitialized = true;
522  return *this;
523 }
524 
530 template<typename MatrixType, int _UpLo>
531 template<typename Derived>
532 LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType,_UpLo>::RealScalar& sigma)
533 {
534  typedef typename TranspositionType::StorageIndex IndexType;
535  const Index size = w.rows();
536  if (m_isInitialized)
537  {
538  eigen_assert(m_matrix.rows()==size);
539  }
540  else
541  {
542  m_matrix.resize(size,size);
543  m_matrix.setZero();
544  m_transpositions.resize(size);
545  for (Index i = 0; i < size; i++)
546  m_transpositions.coeffRef(i) = IndexType(i);
547  m_temporary.resize(size);
548  m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
549  m_isInitialized = true;
550  }
551 
552  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
553 
554  return *this;
555 }
556 
557 #ifndef EIGEN_PARSED_BY_DOXYGEN
558 template<typename _MatrixType, int _UpLo>
559 template<typename RhsType, typename DstType>
560 void LDLT<_MatrixType,_UpLo>::_solve_impl(const RhsType &rhs, DstType &dst) const
561 {
562  eigen_assert(rhs.rows() == rows());
563  // dst = P b
564  dst = m_transpositions * rhs;
565 
566  // dst = L^-1 (P b)
567  matrixL().solveInPlace(dst);
568 
569  // dst = D^-1 (L^-1 P b)
570  // more precisely, use pseudo-inverse of D (see bug 241)
571  using std::abs;
572  const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
573  // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
574  // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
575  // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
576  // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
577  // diagonal element is not well justified and leads to numerical issues in some cases.
578  // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
579  // Using numeric_limits::min() gives us more robustness to denormals.
580  RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
581 
582  for (Index i = 0; i < vecD.size(); ++i)
583  {
584  if(abs(vecD(i)) > tolerance)
585  dst.row(i) /= vecD(i);
586  else
587  dst.row(i).setZero();
588  }
589 
590  // dst = L^-T (D^-1 L^-1 P b)
591  matrixU().solveInPlace(dst);
592 
593  // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
594  dst = m_transpositions.transpose() * dst;
595 }
596 #endif
597 
611 template<typename MatrixType,int _UpLo>
612 template<typename Derived>
613 bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
614 {
615  eigen_assert(m_isInitialized && "LDLT is not initialized.");
616  eigen_assert(m_matrix.rows() == bAndX.rows());
617 
618  bAndX = this->solve(bAndX);
619 
620  return true;
621 }
622 
626 template<typename MatrixType, int _UpLo>
628 {
629  eigen_assert(m_isInitialized && "LDLT is not initialized.");
630  const Index size = m_matrix.rows();
631  MatrixType res(size,size);
632 
633  // P
634  res.setIdentity();
635  res = transpositionsP() * res;
636  // L^* P
637  res = matrixU() * res;
638  // D(L^*P)
639  res = vectorD().real().asDiagonal() * res;
640  // L(DL^*P)
641  res = matrixL() * res;
642  // P^T (LDL^*P)
643  res = transpositionsP().transpose() * res;
644 
645  return res;
646 }
647 
652 template<typename MatrixType, unsigned int UpLo>
655 {
656  return LDLT<PlainObject,UpLo>(m_matrix);
657 }
658 
663 template<typename Derived>
666 {
667  return LDLT<PlainObject>(derived());
668 }
669 
670 } // end namespace Eigen
671 
672 #endif // EIGEN_LDLT_H
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: LDLT.h:253
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:50
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:654
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isfinite_op< typename Derived::Scalar >, const Derived > isfinite(const Eigen::ArrayBase< Derived > &x)
const TranspositionType & transpositionsP() const
Definition: LDLT.h:156
Traits::MatrixU matrixU() const
Definition: LDLT.h:141
RealScalar rcond() const
Definition: LDLT.h:217
LDLT(Index size)
Default Constructor with memory preallocation.
Definition: LDLT.h:90
Namespace containing all symbols from the Eigen library.
Definition: Core:306
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:150
Derived & derived()
Definition: EigenBase.h:45
void setZero()
Definition: LDLT.h:135
Traits::MatrixL matrixL() const
Definition: LDLT.h:148
bool isNegative(void) const
Definition: LDLT.h:177
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:163
Definition: EigenBase.h:29
LDLT()
Default Constructor.
Definition: LDLT.h:77
Definition: Constants.h:204
Index rows() const
Definition: EigenBase.h:59
Derived & derived()
Definition: EigenBase.h:45
LDLT(const EigenBase< InputType > &matrix)
Constructor with decomposition.
Definition: LDLT.h:105
Definition: Constants.h:434
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: LDLT.h:200
const LDLT & adjoint() const
Definition: LDLT.h:243
Definition: Constants.h:206
Definition: Constants.h:432
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:665
bool isPositive() const
Definition: LDLT.h:170
Eigen::Index Index
Definition: LDLT.h:63
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Definition: Eigen_Colamd.h:50
const MatrixType & matrixLDLT() const
Definition: LDLT.h:230
MatrixType reconstructedMatrix() const
Definition: LDLT.h:627
Index rows() const
Definition: EigenBase.h:59
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
Pseudo expression representing a solving operation.
Definition: Solve.h:62
Index cols() const
Definition: EigenBase.h:62
LDLT(EigenBase< InputType > &matrix)
Constructs a LDLT factorization from a given matrix.
Definition: LDLT.h:122
ComputationInfo
Definition: Constants.h:430
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48