Medial Code Documentation
Loading...
Searching...
No Matches
SelfAdjointView.h
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
13namespace Eigen {
14
31namespace internal {
32template<typename MatrixType, unsigned int UpLo>
33struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
34{
36 typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
37 typedef MatrixType ExpressionType;
38 typedef typename MatrixType::PlainObject FullMatrixType;
39 enum {
40 Mode = UpLo | SelfAdjoint,
41 FlagsLvalueBit = is_lvalue<MatrixType>::value ? LvalueBit : 0,
42 Flags = MatrixTypeNestedCleaned::Flags & (HereditaryBits|FlagsLvalueBit)
43 & (~(PacketAccessBit | DirectAccessBit | LinearAccessBit)) // FIXME these flags should be preserved
44 };
45};
46}
47
48
49template<typename _MatrixType, unsigned int UpLo> class SelfAdjointView
50 : public TriangularBase<SelfAdjointView<_MatrixType, UpLo> >
51{
52 public:
53
54 typedef _MatrixType MatrixType;
56 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNested MatrixTypeNested;
57 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
58 typedef MatrixTypeNestedCleaned NestedExpression;
59
62 typedef typename MatrixType::StorageIndex StorageIndex;
63 typedef typename internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType;
65
66 enum {
69 TransposeMode = ((int(Mode) & int(Upper)) ? Lower : 0) | ((int(Mode) & int(Lower)) ? Upper : 0)
70 };
71 typedef typename MatrixType::PlainObject PlainObject;
72
73 EIGEN_DEVICE_FUNC
74 explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix)
75 {
76 EIGEN_STATIC_ASSERT(UpLo==Lower || UpLo==Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY);
77 }
78
79 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
80 inline Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); }
81 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
82 inline Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
83 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
84 inline Index outerStride() const EIGEN_NOEXCEPT { return m_matrix.outerStride(); }
85 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR
86 inline Index innerStride() const EIGEN_NOEXCEPT { return m_matrix.innerStride(); }
87
91 EIGEN_DEVICE_FUNC
92 inline Scalar coeff(Index row, Index col) const
93 {
94 Base::check_coordinates_internal(row, col);
95 return m_matrix.coeff(row, col);
96 }
97
101 EIGEN_DEVICE_FUNC
102 inline Scalar& coeffRef(Index row, Index col)
103 {
104 EIGEN_STATIC_ASSERT_LVALUE(SelfAdjointView);
105 Base::check_coordinates_internal(row, col);
106 return m_matrix.coeffRef(row, col);
107 }
108
110 EIGEN_DEVICE_FUNC
111 const MatrixTypeNestedCleaned& _expression() const { return m_matrix; }
112
113 EIGEN_DEVICE_FUNC
114 const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
115 EIGEN_DEVICE_FUNC
116 MatrixTypeNestedCleaned& nestedExpression() { return m_matrix; }
117
119 template<typename OtherDerived>
120 EIGEN_DEVICE_FUNC
121 const Product<SelfAdjointView,OtherDerived>
123 {
124 return Product<SelfAdjointView,OtherDerived>(*this, rhs.derived());
125 }
126
128 template<typename OtherDerived> friend
129 EIGEN_DEVICE_FUNC
132 {
133 return Product<OtherDerived,SelfAdjointView>(lhs.derived(),rhs);
134 }
135
136 friend EIGEN_DEVICE_FUNC
137 const SelfAdjointView<const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar,MatrixType,product),UpLo>
138 operator*(const Scalar& s, const SelfAdjointView& mat)
139 {
140 return (s*mat.nestedExpression()).template selfadjointView<UpLo>();
141 }
142
153 template<typename DerivedU, typename DerivedV>
154 EIGEN_DEVICE_FUNC
156
167 template<typename DerivedU>
168 EIGEN_DEVICE_FUNC
170
181 template<unsigned int TriMode>
182 EIGEN_DEVICE_FUNC
183 typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
187 {
188 typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::ConstTransposeReturnType>::type tmp1(m_matrix);
189 typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)), MatrixType&, typename MatrixType::AdjointReturnType>::type tmp2(tmp1);
190 return typename internal::conditional<(TriMode&(Upper|Lower))==(UpLo&(Upper|Lower)),
193 }
194
197 EIGEN_DEVICE_FUNC
198 inline const ConjugateReturnType conjugate() const
199 { return ConjugateReturnType(m_matrix.conjugate()); }
200
204 template<bool Cond>
205 EIGEN_DEVICE_FUNC
208 {
210 return ReturnType(m_matrix.template conjugateIf<Cond>());
211 }
212
215 EIGEN_DEVICE_FUNC
216 inline const AdjointReturnType adjoint() const
217 { return AdjointReturnType(m_matrix.adjoint()); }
218
221 EIGEN_DEVICE_FUNC
223 {
224 EIGEN_STATIC_ASSERT_LVALUE(MatrixType)
225 typename MatrixType::TransposeReturnType tmp(m_matrix);
226 return TransposeReturnType(tmp);
227 }
228
231 EIGEN_DEVICE_FUNC
233 {
234 return ConstTransposeReturnType(m_matrix.transpose());
235 }
236
242 EIGEN_DEVICE_FUNC
243 typename MatrixType::ConstDiagonalReturnType diagonal() const
244 {
245 return typename MatrixType::ConstDiagonalReturnType(m_matrix);
246 }
247
249
250 const LLT<PlainObject, UpLo> llt() const;
251 const LDLT<PlainObject, UpLo> ldlt() const;
252
254
259
260 EIGEN_DEVICE_FUNC
262 EIGEN_DEVICE_FUNC
263 RealScalar operatorNorm() const;
264
265 protected:
266 MatrixTypeNested m_matrix;
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
279namespace 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)
284template<typename MatrixType, unsigned int Mode>
290
291template<int UpLo, int SetOpposite, typename DstEvaluatorTypeT, typename SrcEvaluatorTypeT, typename Functor, int Version>
293 : public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version>
294{
295protected:
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;
302public:
303
304 typedef typename Base::DstEvaluatorType DstEvaluatorType;
305 typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
306 typedef typename Base::Scalar Scalar;
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 {
316 eigen_internal_assert(row!=col);
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/***************************************************************************
334* Implementation of MatrixBase methods
335***************************************************************************/
336
338template<typename Derived>
339template<unsigned int UpLo>
340EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
345
355template<typename Derived>
356template<unsigned int UpLo>
357EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type
362
363} // end namespace Eigen
364
365#endif // EIGEN_SELFADJOINTMATRIX_H
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition SelfAdjointView.h:51
friend EIGEN_DEVICE_FUNC const Product< OtherDerived, SelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
Efficient vector/matrix times triangular matrix product.
Definition SelfAdjointView.h:131
EIGEN_DEVICE_FUNC internal::conditional<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typenameMatrixType::AdjointReturnType, TriMode > >::type triangularView() const
Definition SelfAdjointView.h:186
EIGEN_DEVICE_FUNC SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
Perform a symmetric rank 2 update of the selfadjoint matrix *this: .
EIGEN_DEVICE_FUNC const AdjointReturnType adjoint() const
Definition SelfAdjointView.h:216
EIGEN_DEVICE_FUNC SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
Perform a symmetric rank K update of the selfadjoint matrix *this: where u is a vector or matrix.
const LLT< PlainObject, UpLo > llt() const
\cholesky_module
Definition LLT.h:551
EIGEN_DEVICE_FUNC const Product< SelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
Efficient triangular matrix times vector/matrix product.
Definition SelfAdjointView.h:122
const LDLT< PlainObject, UpLo > ldlt() const
\cholesky_module
Definition LDLT.h:670
EIGEN_DEVICE_FUNC RealScalar operatorNorm() const
Computes the L2 operator norm.
Definition MatrixBaseEigenvalues.h:151
EIGEN_DEVICE_FUNC EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Definition MatrixBaseEigenvalues.h:88
Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
Return type of eigenvalues()
Definition SelfAdjointView.h:258
EIGEN_DEVICE_FUNC const ConstTransposeReturnType transpose() const
Definition SelfAdjointView.h:232
EIGEN_DEVICE_FUNC Scalar & coeffRef(Index row, Index col)
Definition SelfAdjointView.h:102
EIGEN_DEVICE_FUNC internal::conditional< Cond, ConjugateReturnType, ConstSelfAdjointView >::type conjugateIf() const
Definition SelfAdjointView.h:207
EIGEN_DEVICE_FUNC TransposeReturnType transpose()
Definition SelfAdjointView.h:222
EIGEN_DEVICE_FUNC MatrixType::ConstDiagonalReturnType diagonal() const
Definition SelfAdjointView.h:243
EIGEN_DEVICE_FUNC Scalar coeff(Index row, Index col) const
Definition SelfAdjointView.h:92
EIGEN_DEVICE_FUNC const ConjugateReturnType conjugate() const
Definition SelfAdjointView.h:198
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
Definition SelfAdjointView.h:61
NumTraits< Scalar >::Real RealScalar
Real part of Scalar.
Definition SelfAdjointView.h:256
Base class for triangular part in a matrix.
Definition TriangularMatrix.h:28
Definition AssignEvaluator.h:619
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void assignCoeff(Index row, Index col)
Assign src(row,col) to dst(row,col) through the assignment functor.
Definition AssignEvaluator.h:652
@ SelfAdjoint
Used in BandMatrix and SelfAdjointView to indicate that the matrix is self-adjoint.
Definition Constants.h:225
@ Lower
View matrix as a lower triangular matrix.
Definition Constants.h:209
@ Upper
View matrix as an upper triangular matrix.
Definition Constants.h:211
const unsigned int PacketAccessBit
Short version: means the expression might be vectorized.
Definition Constants.h:94
const unsigned int LinearAccessBit
Short version: means the expression can be seen as 1D vector.
Definition Constants.h:130
const unsigned int DirectAccessBit
Means that the underlying array of coefficients can be directly accessed as a plain strided array.
Definition Constants.h:155
const unsigned int LvalueBit
Means the expression has a coeffRef() method, i.e.
Definition Constants.h:144
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
Eigen::Index Index
The interface type of indices.
Definition EigenBase.h:39
Definition Constants.h:534
Definition Constants.h:542
Definition Meta.h:109
Definition CoreEvaluators.h:80
Definition XprHelper.h:660
Definition ForwardDeclarations.h:17
Definition Meta.h:96