SelfAdjointView.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) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
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_SELFADJOINTMATRIX_H
11 #define EIGEN_SELFADJOINTMATRIX_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
33 namespace internal {
34 template<typename MatrixType, unsigned int UpLo>
35 struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
36 {
37  typedef typename ref_selector<MatrixType>::non_const_type MatrixTypeNested;
38  typedef remove_all_t<MatrixTypeNested> MatrixTypeNestedCleaned;
39  typedef MatrixType ExpressionType;
40  typedef typename MatrixType::PlainObject FullMatrixType;
41  enum {
42  Mode = UpLo | SelfAdjoint,
43  FlagsLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
44  Flags = MatrixTypeNestedCleaned::Flags & (HereditaryBits|FlagsLvalueBit)
45  & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)) // FIXME these flags should be preserved
46  };
47 };
48 }
49 
50 
51 template<typename MatrixType_, unsigned int UpLo> class SelfAdjointView
52  : public TriangularBase<SelfAdjointView<MatrixType_, UpLo> >
53 {
54  public:
55  EIGEN_STATIC_ASSERT(UpLo==Lower || UpLo==Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY)
56 
57  typedef MatrixType_ MatrixType;
62 
64  typedef typename internal::traits<SelfAdjointView>::Scalar Scalar;
67  typedef SelfAdjointView<std::add_const_t<MatrixType>, UpLo> ConstSelfAdjointView;
68 
69  enum {
70  Mode = internal::traits<SelfAdjointView>::Mode,
71  Flags = internal::traits<SelfAdjointView>::Flags,
72  TransposeMode = ((int(Mode) & int(Upper)) ? Lower : 0) | ((int(Mode) & int(Lower)) ? Upper : 0)
73  };
75 
77  explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix) { }
78 
80  inline Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); }
82  inline Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
84  inline Index outerStride() const EIGEN_NOEXCEPT { return m_matrix.outerStride(); }
86  inline Index innerStride() const EIGEN_NOEXCEPT { return m_matrix.innerStride(); }
87 
92  inline Scalar coeff(Index row, Index col) const
93  {
95  return m_matrix.coeff(row, col);
96  }
97 
103  {
106  return m_matrix.coeffRef(row, col);
107  }
108 
111  const MatrixTypeNestedCleaned& _expression() const { return m_matrix; }
112 
117 
119  template<typename OtherDerived>
123  {
124  return Product<SelfAdjointView,OtherDerived>(*this, rhs.derived());
125  }
126 
128  template<typename OtherDerived> friend
132  {
134  }
135 
136  friend EIGEN_DEVICE_FUNC
139  {
140  return (s*mat.nestedExpression()).template selfadjointView<UpLo>();
141  }
142 
153  template<typename DerivedU, typename DerivedV>
156 
167  template<typename DerivedU>
170 
181  template<unsigned int TriMode>
183  std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
187  {
188  std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::ConstTransposeReturnType> tmp1(m_matrix);
189  std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::AdjointReturnType> tmp2(tmp1);
190  return std::conditional_t<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
193  }
194 
198  inline const ConjugateReturnType conjugate() const
199  { return ConjugateReturnType(m_matrix.conjugate()); }
200 
204  template<bool Cond>
206  inline std::conditional_t<Cond,ConjugateReturnType,ConstSelfAdjointView>
207  conjugateIf() const
208  {
209  typedef std::conditional_t<Cond,ConjugateReturnType,ConstSelfAdjointView> ReturnType;
210  return ReturnType(m_matrix.template conjugateIf<Cond>());
211  }
212 
216  inline const AdjointReturnType adjoint() const
217  { return AdjointReturnType(m_matrix.adjoint()); }
218 
221  template<class Dummy=int>
223  inline TransposeReturnType transpose(std::enable_if_t<Eigen::internal::is_lvalue<MatrixType>::value, Dummy*> = nullptr)
224  {
226  return TransposeReturnType(tmp);
227  }
228 
232  inline const ConstTransposeReturnType transpose() const
233  {
234  return ConstTransposeReturnType(m_matrix.transpose());
235  }
236 
244  {
246  }
247 
249 
250  const LLT<PlainObject, UpLo> llt() const;
251  const LDLT<PlainObject, UpLo> ldlt() const;
252 
254 
259 
263  RealScalar operatorNorm() const;
264 
265  protected:
267 };
268 
269 
270 // template<typename OtherDerived, typename MatrixType, unsigned int UpLo>
271 // internal::selfadjoint_matrix_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >
272 // operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView<MatrixType,UpLo>& rhs)
273 // {
274 // return internal::matrix_selfadjoint_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >(lhs.derived(),rhs);
275 // }
276 
277 // selfadjoint to dense matrix
278 
279 namespace internal {
280 
281 // TODO currently a selfadjoint expression has the form SelfAdjointView<.,.>
282 // in the future selfadjoint-ness should be defined by the expression traits
283 // such that Transpose<SelfAdjointView<.,.> > is valid. (currently TriangularBase::transpose() is overloaded to make it work)
284 template<typename MatrixType, unsigned int Mode>
285 struct evaluator_traits<SelfAdjointView<MatrixType,Mode> >
286 {
287  typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
288  typedef SelfAdjointShape Shape;
289 };
290 
291 template<int UpLo, int SetOpposite, typename DstEvaluatorTypeT, typename SrcEvaluatorTypeT, typename Functor, int Version>
292 class triangular_dense_assignment_kernel<UpLo,SelfAdjoint,SetOpposite,DstEvaluatorTypeT,SrcEvaluatorTypeT,Functor,Version>
293  : public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version>
294 {
295 protected:
296  typedef generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> Base;
297  typedef typename Base::DstXprType DstXprType;
298  typedef typename Base::SrcXprType SrcXprType;
299  using Base::m_dst;
300  using Base::m_src;
301  using Base::m_functor;
302 public:
303 
304  typedef typename Base::DstEvaluatorType DstEvaluatorType;
305  typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
306  typedef typename Base::Scalar Scalar;
307  typedef typename Base::AssignmentTraits AssignmentTraits;
308 
309 
310  EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType &dst, const SrcEvaluatorType &src, const Functor &func, DstXprType& dstExpr)
311  : Base(dst, src, func, dstExpr)
312  {}
313 
314  EIGEN_DEVICE_FUNC void assignCoeff(Index row, Index col)
315  {
317  Scalar tmp = m_src.coeff(row,col);
318  m_functor.assignCoeff(m_dst.coeffRef(row,col), tmp);
319  m_functor.assignCoeff(m_dst.coeffRef(col,row), numext::conj(tmp));
320  }
321 
322  EIGEN_DEVICE_FUNC void assignDiagonalCoeff(Index id)
323  {
324  Base::assignCoeff(id,id);
325  }
326 
327  EIGEN_DEVICE_FUNC void assignOppositeCoeff(Index, Index)
328  { eigen_internal_assert(false && "should never be called"); }
329 };
330 
331 } // end namespace internal
332 
333 
338 template<typename Derived>
339 template<unsigned int UpLo>
340 EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
342 {
343  return typename ConstSelfAdjointViewReturnType<UpLo>::Type(derived());
344 }
345 
355 template<typename Derived>
356 template<unsigned int UpLo>
359 {
360  return typename SelfAdjointViewReturnType<UpLo>::Type(derived());
361 }
362 
363 } // end namespace Eigen
364 
365 #endif // EIGEN_SELFADJOINTMATRIX_H
Array< int, Dynamic, 1 > v
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().
#define eigen_internal_assert(x)
Definition: Macros.h:908
#define EIGEN_NOEXCEPT
Definition: Macros.h:1260
#define EIGEN_CONSTEXPR
Definition: Macros.h:747
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(SCALAR, EXPR, OPNAME)
Definition: Macros.h:1207
#define EIGEN_STATIC_ASSERT_LVALUE(Derived)
Definition: StaticAssert.h:96
#define EIGEN_STATIC_ASSERT(X, MSG)
Definition: StaticAssert.h:26
Matrix< float, 1, Dynamic > MatrixType
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:67
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:63
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
SelfAdjointViewReturnType< UpLo >::Type selfadjointView()
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
Base::PlainObject PlainObject
Definition: Matrix.h:194
Expression of the product of two arbitrary matrices or vectors.
Definition: Product.h:77
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
std::conditional_t< Cond, ConjugateReturnType, ConstSelfAdjointView > conjugateIf() const
EIGEN_CONSTEXPR Index innerStride() const EIGEN_NOEXCEPT
RealScalar operatorNorm() const
Computes the L2 operator norm.
SelfAdjointView< const typename MatrixType::ConstTransposeReturnType, TransposeMode > ConstTransposeReturnType
const MatrixTypeNestedCleaned & nestedExpression() const
TransposeReturnType transpose(std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >=nullptr)
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
const MatrixTypeNestedCleaned & _expression() const
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:548
MatrixType::PlainObject PlainObject
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
const AdjointReturnType adjoint() const
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:667
internal::traits< SelfAdjointView >::MatrixTypeNested MatrixTypeNested
SelfAdjointView(MatrixType &matrix)
const Product< SelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
MatrixType::StorageIndex StorageIndex
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
SelfAdjointView< typename MatrixType::TransposeReturnType, TransposeMode > TransposeReturnType
friend const Product< OtherDerived, SelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
SelfAdjointView< const MatrixConjugateReturnType, UpLo > ConjugateReturnType
MatrixType::ConstDiagonalReturnType diagonal() const
const ConstTransposeReturnType transpose() const
Scalar coeff(Index row, Index col) const
internal::remove_all_t< typename MatrixType::ConjugateReturnType > MatrixConjugateReturnType
EIGEN_CONSTEXPR Index outerStride() const EIGEN_NOEXCEPT
std::conditional_t<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typename MatrixType::AdjointReturnType, TriMode > > triangularView() const
NumTraits< Scalar >::Real RealScalar
Scalar & coeffRef(Index row, Index col)
MatrixTypeNestedCleaned NestedExpression
const ConjugateReturnType conjugate() const
SelfAdjointView< const typename MatrixType::AdjointReturnType, TransposeMode > AdjointReturnType
MatrixTypeNestedCleaned & nestedExpression()
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
internal::traits< SelfAdjointView >::MatrixTypeNestedCleaned MatrixTypeNestedCleaned
MatrixTypeNested m_matrix
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
friend const SelfAdjointView< const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar, MatrixType, product), UpLo > operator*(const Scalar &s, const SelfAdjointView &mat)
Expression of the transpose of a matrix.
Definition: Transpose.h:56
Base class for triangular part in a matrix.
void check_coordinates_internal(Index, Index) const
Expression of a triangular part in a matrix.
@ SelfAdjoint
Definition: Constants.h:227
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
const unsigned int PacketAccessBit
Definition: Constants.h:96
const unsigned int LinearAccessBit
Definition: Constants.h:132
const unsigned int DirectAccessBit
Definition: Constants.h:157
const unsigned int LvalueBit
Definition: Constants.h:146
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
: InteropHeaders
Definition: Core:139
const unsigned int HereditaryBits
Definition: Constants.h:197
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
Definition: BFloat16.h:222
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231