Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
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 <[email protected]>
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// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
34namespace internal {
35template <typename MatrixType, unsigned int UpLo>
36struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType> {
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} // namespace internal
49
50template <typename MatrixType_, unsigned int UpLo>
51class SelfAdjointView : public TriangularBase<SelfAdjointView<MatrixType_, UpLo> > {
52 public:
53 EIGEN_STATIC_ASSERT(UpLo == Lower || UpLo == Upper, SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY)
54
55 typedef MatrixType_ MatrixType;
57 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNested MatrixTypeNested;
58 typedef typename internal::traits<SelfAdjointView>::MatrixTypeNestedCleaned MatrixTypeNestedCleaned;
59 typedef MatrixTypeNestedCleaned NestedExpression;
60
62 typedef typename internal::traits<SelfAdjointView>::Scalar Scalar;
63 typedef typename MatrixType::StorageIndex StorageIndex;
64 typedef internal::remove_all_t<typename MatrixType::ConjugateReturnType> MatrixConjugateReturnType;
66
67 enum {
68 Mode = internal::traits<SelfAdjointView>::Mode,
69 Flags = internal::traits<SelfAdjointView>::Flags,
70 TransposeMode = ((int(Mode) & int(Upper)) ? Lower : 0) | ((int(Mode) & int(Lower)) ? Upper : 0)
71 };
72 typedef typename MatrixType::PlainObject PlainObject;
73
74 EIGEN_DEVICE_FUNC explicit inline SelfAdjointView(MatrixType& matrix) : m_matrix(matrix) {}
75
76 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline Index rows() const EIGEN_NOEXCEPT { return m_matrix.rows(); }
77 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline Index cols() const EIGEN_NOEXCEPT { return m_matrix.cols(); }
78 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline Index outerStride() const EIGEN_NOEXCEPT { return m_matrix.outerStride(); }
79 EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR inline Index innerStride() const EIGEN_NOEXCEPT { return m_matrix.innerStride(); }
80
84 EIGEN_DEVICE_FUNC inline Scalar coeff(Index row, Index col) const {
85 Base::check_coordinates_internal(row, col);
86 return m_matrix.coeff(row, col);
87 }
88
92 EIGEN_DEVICE_FUNC inline Scalar& coeffRef(Index row, Index col) {
93 EIGEN_STATIC_ASSERT_LVALUE(SelfAdjointView);
94 Base::check_coordinates_internal(row, col);
95 return m_matrix.coeffRef(row, col);
96 }
97
99 EIGEN_DEVICE_FUNC const MatrixTypeNestedCleaned& _expression() const { return m_matrix; }
100
101 EIGEN_DEVICE_FUNC const MatrixTypeNestedCleaned& nestedExpression() const { return m_matrix; }
102 EIGEN_DEVICE_FUNC MatrixTypeNestedCleaned& nestedExpression() { return m_matrix; }
103
105 template <typename OtherDerived>
108 }
109
111 template <typename OtherDerived>
113 const SelfAdjointView& rhs) {
115 }
116
117 friend EIGEN_DEVICE_FUNC const
118 SelfAdjointView<const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar, MatrixType, product), UpLo>
119 operator*(const Scalar& s, const SelfAdjointView& mat) {
120 return (s * mat.nestedExpression()).template selfadjointView<UpLo>();
121 }
122
133 template <typename DerivedU, typename DerivedV>
135 const Scalar& alpha = Scalar(1));
136
147 template <typename DerivedU>
148 EIGEN_DEVICE_FUNC SelfAdjointView& rankUpdate(const MatrixBase<DerivedU>& u, const Scalar& alpha = Scalar(1));
149
161 template <unsigned int TriMode>
162 EIGEN_DEVICE_FUNC
163 std::conditional_t<(TriMode & (Upper | Lower)) == (UpLo & (Upper | Lower)), TriangularView<MatrixType, TriMode>,
166 std::conditional_t<(TriMode & (Upper | Lower)) == (UpLo & (Upper | Lower)), MatrixType&,
167 typename MatrixType::ConstTransposeReturnType>
168 tmp1(m_matrix);
169 std::conditional_t<(TriMode & (Upper | Lower)) == (UpLo & (Upper | Lower)), MatrixType&,
170 typename MatrixType::AdjointReturnType>
171 tmp2(tmp1);
172 return std::conditional_t<(TriMode & (Upper | Lower)) == (UpLo & (Upper | Lower)),
175 }
176
179 EIGEN_DEVICE_FUNC inline const ConjugateReturnType conjugate() const {
180 return ConjugateReturnType(m_matrix.conjugate());
181 }
182
186 template <bool Cond>
187 EIGEN_DEVICE_FUNC inline std::conditional_t<Cond, ConjugateReturnType, ConstSelfAdjointView> conjugateIf() const {
188 typedef std::conditional_t<Cond, ConjugateReturnType, ConstSelfAdjointView> ReturnType;
189 return ReturnType(m_matrix.template conjugateIf<Cond>());
190 }
191
194 EIGEN_DEVICE_FUNC inline const AdjointReturnType adjoint() const { return AdjointReturnType(m_matrix.adjoint()); }
195
198 template <class Dummy = int>
199 EIGEN_DEVICE_FUNC inline TransposeReturnType transpose(
200 std::enable_if_t<Eigen::internal::is_lvalue<MatrixType>::value, Dummy*> = nullptr) {
201 typename MatrixType::TransposeReturnType tmp(m_matrix);
202 return TransposeReturnType(tmp);
203 }
204
207 EIGEN_DEVICE_FUNC inline const ConstTransposeReturnType transpose() const {
208 return ConstTransposeReturnType(m_matrix.transpose());
209 }
210
216 EIGEN_DEVICE_FUNC typename MatrixType::ConstDiagonalReturnType diagonal() const {
217 return typename MatrixType::ConstDiagonalReturnType(m_matrix);
218 }
219
221
222 const LLT<PlainObject, UpLo> llt() const;
223 const LDLT<PlainObject, UpLo> ldlt() const;
224
226
228 typedef typename NumTraits<Scalar>::Real RealScalar;
231
232 EIGEN_DEVICE_FUNC EigenvaluesReturnType eigenvalues() const;
233 EIGEN_DEVICE_FUNC RealScalar operatorNorm() const;
234
235 protected:
236 MatrixTypeNested m_matrix;
237};
238
239// template<typename OtherDerived, typename MatrixType, unsigned int UpLo>
240// internal::selfadjoint_matrix_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo> >
241// operator*(const MatrixBase<OtherDerived>& lhs, const SelfAdjointView<MatrixType,UpLo>& rhs)
242// {
243// return internal::matrix_selfadjoint_product_returntype<OtherDerived,SelfAdjointView<MatrixType,UpLo>
244// >(lhs.derived(),rhs);
245// }
246
247// selfadjoint to dense matrix
248
249namespace internal {
250
251// TODO currently a selfadjoint expression has the form SelfAdjointView<.,.>
252// in the future selfadjoint-ness should be defined by the expression traits
253// such that Transpose<SelfAdjointView<.,.> > is valid. (currently TriangularBase::transpose() is overloaded to
254// make it work)
255template <typename MatrixType, unsigned int Mode>
256struct evaluator_traits<SelfAdjointView<MatrixType, Mode> > {
257 typedef typename storage_kind_to_evaluator_kind<typename MatrixType::StorageKind>::Kind Kind;
258 typedef SelfAdjointShape Shape;
259};
260
261template <int UpLo, int SetOpposite, typename DstEvaluatorTypeT, typename SrcEvaluatorTypeT, typename Functor,
262 int Version>
263class triangular_dense_assignment_kernel<UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor,
264 Version>
265 : public generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> {
266 protected:
267 typedef generic_dense_assignment_kernel<DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version> Base;
268 typedef typename Base::DstXprType DstXprType;
269 typedef typename Base::SrcXprType SrcXprType;
270 using Base::m_dst;
271 using Base::m_functor;
272 using Base::m_src;
273
274 public:
275 typedef typename Base::DstEvaluatorType DstEvaluatorType;
276 typedef typename Base::SrcEvaluatorType SrcEvaluatorType;
277 typedef typename Base::Scalar Scalar;
278 typedef typename Base::AssignmentTraits AssignmentTraits;
279
280 EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType& dst, const SrcEvaluatorType& src,
281 const Functor& func, DstXprType& dstExpr)
282 : Base(dst, src, func, dstExpr) {}
283
284 EIGEN_DEVICE_FUNC void assignCoeff(Index row, Index col) {
285 eigen_internal_assert(row != col);
286 Scalar tmp = m_src.coeff(row, col);
287 m_functor.assignCoeff(m_dst.coeffRef(row, col), tmp);
288 m_functor.assignCoeff(m_dst.coeffRef(col, row), numext::conj(tmp));
289 }
290
291 EIGEN_DEVICE_FUNC void assignDiagonalCoeff(Index id) { Base::assignCoeff(id, id); }
292
293 EIGEN_DEVICE_FUNC void assignOppositeCoeff(Index, Index) { eigen_internal_assert(false && "should never be called"); }
294};
295
296} // end namespace internal
297
298/***************************************************************************
299 * Implementation of MatrixBase methods
300 ***************************************************************************/
301
303template <typename Derived>
304template <unsigned int UpLo>
305EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::template ConstSelfAdjointViewReturnType<UpLo>::Type
309
320template <typename Derived>
321template <unsigned int UpLo>
322EIGEN_DEVICE_FUNC typename MatrixBase<Derived>::template SelfAdjointViewReturnType<UpLo>::Type
326
327} // end namespace Eigen
328
329#endif // EIGEN_SELFADJOINTMATRIX_H
Derived & derived()
Definition EigenBase.h:49
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:186
Expression of the product of two arbitrary matrices or vectors.
Definition Product.h:202
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition SelfAdjointView.h:51
Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
Definition SelfAdjointView.h:230
RealScalar operatorNorm() const
Computes the L2 operator norm.
Definition MatrixBaseEigenvalues.h:136
TransposeReturnType transpose(std::enable_if_t< Eigen::internal::is_lvalue< MatrixType >::value, Dummy * >=nullptr)
Definition SelfAdjointView.h:199
std::conditional_t<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typename MatrixType::AdjointReturnType, TriMode > > triangularView() const
Definition SelfAdjointView.h:165
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
const LLT< PlainObject, UpLo > llt() const
Definition LLT.h:507
const AdjointReturnType adjoint() const
Definition SelfAdjointView.h:194
const LDLT< PlainObject, UpLo > ldlt() const
Definition LDLT.h:634
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
Definition SelfAdjointView.h:62
SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const Scalar &alpha=Scalar(1))
MatrixType::ConstDiagonalReturnType diagonal() const
Definition SelfAdjointView.h:216
const ConstTransposeReturnType transpose() const
Definition SelfAdjointView.h:207
Scalar coeff(Index row, Index col) const
Definition SelfAdjointView.h:84
NumTraits< Scalar >::Real RealScalar
Definition SelfAdjointView.h:228
const ConjugateReturnType conjugate() const
Definition SelfAdjointView.h:179
Scalar & coeffRef(Index row, Index col)
Definition SelfAdjointView.h:92
const Product< SelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
Definition SelfAdjointView.h:106
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Definition MatrixBaseEigenvalues.h:83
friend const Product< OtherDerived, SelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
Definition SelfAdjointView.h:112
std::conditional_t< Cond, ConjugateReturnType, ConstSelfAdjointView > conjugateIf() const
Definition SelfAdjointView.h:187
Base class for triangular part in a matrix.
Definition TriangularMatrix.h:32
Expression of a triangular part in a matrix.
Definition TriangularMatrix.h:167
@ SelfAdjoint
Definition Constants.h:227
@ Lower
Definition Constants.h:211
@ Upper
Definition Constants.h:213
const unsigned int LvalueBit
Definition Constants.h:148
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83