10#ifndef EIGEN_SELFADJOINTMATRIX_H
11#define EIGEN_SELFADJOINTMATRIX_H
13#include "./InternalHeaderCheck.h"
34template<
typename MatrixType,
unsigned int UpLo>
35struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
37 typedef typename ref_selector<MatrixType>::non_const_type MatrixTypeNested;
38 typedef typename remove_all<MatrixTypeNested>::type MatrixTypeNestedCleaned;
39 typedef MatrixType ExpressionType;
40 typedef typename MatrixType::PlainObject FullMatrixType;
43 FlagsLvalueBit = is_lvalue<MatrixType>::value ?
LvalueBit : 0,
44 Flags = MatrixTypeNestedCleaned::Flags & (HereditaryBits|FlagsLvalueBit)
55 EIGEN_STATIC_ASSERT(UpLo==
Lower || UpLo==
Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY)
57 typedef MatrixType_ MatrixType;
59 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNested MatrixTypeNested;
60 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
61 typedef MatrixTypeNestedCleaned NestedExpression;
64 typedef typename internal::traits<SelfAdjointView>::Scalar
Scalar;
65 typedef typename MatrixType::StorageIndex StorageIndex;
66 typedef typename internal::remove_all<typename MatrixType::ConjugateReturnType>::type MatrixConjugateReturnType;
70 Mode = internal::traits<SelfAdjointView>::Mode,
71 Flags = internal::traits<SelfAdjointView>::Flags,
74 typedef typename MatrixType::PlainObject PlainObject;
77 explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix) { }
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(); }
94 Base::check_coordinates_internal(row, col);
95 return m_matrix.coeff(row, col);
105 Base::check_coordinates_internal(row, col);
106 return m_matrix.coeffRef(row, col);
111 const MatrixTypeNestedCleaned& _expression()
const {
return m_matrix; }
114 const MatrixTypeNestedCleaned& nestedExpression()
const {
return m_matrix; }
116 MatrixTypeNestedCleaned& nestedExpression() {
return m_matrix; }
119 template<
typename OtherDerived>
121 const Product<SelfAdjointView,OtherDerived>
128 template<
typename OtherDerived>
friend
136 friend EIGEN_DEVICE_FUNC
140 return (s*mat.nestedExpression()).template selfadjointView<UpLo>();
153 template<
typename DerivedU,
typename DerivedV>
167 template<
typename DerivedU>
181 template<
unsigned int TriMode>
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);
206 inline typename internal::conditional<Cond,ConjugateReturnType,ConstSelfAdjointView>::type
209 typedef typename internal::conditional<Cond,ConjugateReturnType,ConstSelfAdjointView>::type ReturnType;
210 return ReturnType(m_matrix.template conjugateIf<Cond>());
221 template<
class Dummy=
int>
225 typename MatrixType::TransposeReturnType tmp(m_matrix);
243 typename MatrixType::ConstDiagonalReturnType
diagonal()
const
245 return typename MatrixType::ConstDiagonalReturnType(m_matrix);
266 MatrixTypeNested m_matrix;
284template<
typename MatrixType,
unsigned int Mode>
287 typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
288 typedef SelfAdjointShape Shape;
291template<
int UpLo,
int SetOpposite,
typename DstEvaluatorTypeT,
typename SrcEvaluatorTypeT,
typename Functor,
int Version>
292class triangular_dense_assignment_kernel<UpLo,
SelfAdjoint,SetOpposite,DstEvaluatorTypeT,SrcEvaluatorTypeT,Functor,Version>
293 :
public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version>
296 typedef generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> Base;
297 typedef typename Base::DstXprType DstXprType;
298 typedef typename Base::SrcXprType SrcXprType;
301 using Base::m_functor;
304 typedef typename Base::DstEvaluatorType DstEvaluatorType;
305 typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
306 typedef typename Base::Scalar Scalar;
307 typedef typename Base::AssignmentTraits AssignmentTraits;
310 EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType &dst,
const SrcEvaluatorType &src,
const Functor &func, DstXprType& dstExpr)
311 : Base(dst, src, func, dstExpr)
314 EIGEN_DEVICE_FUNC
void assignCoeff(
Index row,
Index col)
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));
322 EIGEN_DEVICE_FUNC
void assignDiagonalCoeff(
Index id)
324 Base::assignCoeff(
id,
id);
327 EIGEN_DEVICE_FUNC
void assignOppositeCoeff(
Index,
Index)
328 { eigen_internal_assert(
false &&
"should never be called"); }
338template<
typename Derived>
339template<
unsigned int UpLo>
340EIGEN_DEVICE_FUNC
typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
355template<
typename Derived>
356template<
unsigned int UpLo>
Derived & derived()
Definition: EigenBase.h:48
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
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
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.
Definition: SelfAdjointView.h:53
Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
Definition: SelfAdjointView.h:258
RealScalar operatorNorm() const
Computes the L2 operator norm.
Definition: MatrixBaseEigenvalues.h:153
internal::conditional< Cond, ConjugateReturnType, ConstSelfAdjointView >::type conjugateIf() const
Definition: SelfAdjointView.h:207
TransposeReturnType transpose(typename internal::enable_if< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >::type=nullptr)
Definition: SelfAdjointView.h:223
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:548
const AdjointReturnType adjoint() const
Definition: SelfAdjointView.h:216
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:667
internal::conditional<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typenameMatrixType::AdjointReturnType, TriMode > >::type triangularView() const
Definition: SelfAdjointView.h:186
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
Definition: SelfAdjointView.h:64
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
MatrixType::ConstDiagonalReturnType diagonal() const
Definition: SelfAdjointView.h:243
const ConstTransposeReturnType transpose() const
Definition: SelfAdjointView.h:232
Scalar coeff(Index row, Index col) const
Definition: SelfAdjointView.h:92
NumTraits< Scalar >::Real RealScalar
Definition: SelfAdjointView.h:256
const ConjugateReturnType conjugate() const
Definition: SelfAdjointView.h:198
Scalar & coeffRef(Index row, Index col)
Definition: SelfAdjointView.h:102
const Product< SelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
Definition: SelfAdjointView.h:122
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Definition: MatrixBaseEigenvalues.h:90
friend const Product< OtherDerived, SelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
Definition: SelfAdjointView.h:131
Base class for triangular part in a matrix.
Definition: TriangularMatrix.h:30
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:191
@ 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
Namespace containing all symbols from the Eigen library.
Definition: B01_Experimental.dox:1
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
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:235