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 namespace Eigen {
14 
31 namespace internal {
32 template<typename MatrixType, unsigned int UpLo>
33 struct traits<SelfAdjointView<MatrixType, UpLo> > : traits<MatrixType>
34 {
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 
49 template<typename _MatrixType, unsigned int UpLo> class SelfAdjointView
50  : public TriangularBase<SelfAdjointView<_MatrixType, UpLo> >
51 {
52  public:
53 
54  typedef _MatrixType MatrixType;
59 
62  typedef typename MatrixType::StorageIndex StorageIndex;
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 
75  {
76  EIGEN_STATIC_ASSERT(UpLo==Lower || UpLo==Upper,SELFADJOINTVIEW_ACCEPTS_UPPER_AND_LOWER_MODE_ONLY);
77  }
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  {
133  return Product<OtherDerived,SelfAdjointView>(lhs.derived(),rhs);
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  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 
198  inline const ConjugateReturnType conjugate() const
199  { return ConjugateReturnType(m_matrix.conjugate()); }
200 
204  template<bool Cond>
207  conjugateIf() const
208  {
210  return ReturnType(m_matrix.template conjugateIf<Cond>());
211  }
212 
216  inline const AdjointReturnType adjoint() const
217  { return AdjointReturnType(m_matrix.adjoint()); }
218 
223  {
225  typename MatrixType::TransposeReturnType tmp(m_matrix);
226  return TransposeReturnType(tmp);
227  }
228 
232  inline const ConstTransposeReturnType transpose() const
233  {
234  return ConstTransposeReturnType(m_matrix.transpose());
235  }
236 
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 
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>
286 {
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:
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 
306  typedef typename Base::Scalar Scalar;
308 
309 
311  : Base(dst, src, func, dstExpr)
312  {}
313 
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 
323  {
324  Base::assignCoeff(id,id);
325  }
326 
328  { eigen_internal_assert(false && "should never be called"); }
329 };
330 
331 } // end namespace internal
332 
333 /***************************************************************************
334 * Implementation of MatrixBase methods
335 ***************************************************************************/
336 
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
gtsam.examples.DogLegOptimizerExample.int
int
Definition: DogLegOptimizerExample.py:111
Eigen::internal::generic_dense_assignment_kernel::m_functor
const Functor & m_functor
Definition: AssignEvaluator.h:720
Eigen::SelfAdjointView::ConstSelfAdjointView
SelfAdjointView< typename internal::add_const< MatrixType >::type, UpLo > ConstSelfAdjointView
Definition: SelfAdjointView.h:64
Eigen::internal::traits< SelfAdjointView< MatrixType, UpLo > >::FullMatrixType
MatrixType::PlainObject FullMatrixType
Definition: SelfAdjointView.h:38
Eigen::SelfAdjointView::ConjugateReturnType
SelfAdjointView< const MatrixConjugateReturnType, UpLo > ConjugateReturnType
Definition: SelfAdjointView.h:195
EIGEN_DEVICE_FUNC
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:976
Eigen::internal::generic_dense_assignment_kernel
Definition: AssignEvaluator.h:618
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Eigen::SelfAdjointView::RealScalar
NumTraits< Scalar >::Real RealScalar
Definition: SelfAdjointView.h:256
Eigen::SelfAdjointView::diagonal
EIGEN_DEVICE_FUNC MatrixType::ConstDiagonalReturnType diagonal() const
Definition: SelfAdjointView.h:243
col
m col(1)
Eigen::internal::traits< SelfAdjointView< MatrixType, UpLo > >::ExpressionType
MatrixType ExpressionType
Definition: SelfAdjointView.h:37
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::DstEvaluatorType
Base::DstEvaluatorType DstEvaluatorType
Definition: SelfAdjointView.h:304
Eigen::SelfAdjointView::Base
TriangularBase< SelfAdjointView > Base
Definition: SelfAdjointView.h:55
gtsam.examples.DogLegOptimizerExample.type
type
Definition: DogLegOptimizerExample.py:111
Eigen::internal::traits< SelfAdjointView< MatrixType, UpLo > >::MatrixTypeNested
ref_selector< MatrixType >::non_const_type MatrixTypeNested
Definition: SelfAdjointView.h:35
alpha
RealScalar alpha
Definition: level1_cplx_impl.h:147
s
RealScalar s
Definition: level1_cplx_impl.h:126
Eigen::EigenBase< SelfAdjointView< const Derived, UpLo > >::Index
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:39
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Eigen::SelfAdjointView::NestedExpression
MatrixTypeNestedCleaned NestedExpression
Definition: SelfAdjointView.h:58
Eigen::SelfAdjointView::operator*
friend const EIGEN_DEVICE_FUNC SelfAdjointView< const EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(Scalar, MatrixType, product), UpLo > operator*(const Scalar &s, const SelfAdjointView &mat)
Definition: SelfAdjointView.h:138
Eigen::SelfAdjointView::conjugateIf
EIGEN_DEVICE_FUNC internal::conditional< Cond, ConjugateReturnType, ConstSelfAdjointView >::type conjugateIf() const
Definition: SelfAdjointView.h:207
Eigen::SelfAdjointView::ldlt
const LDLT< PlainObject, UpLo > ldlt() const
Definition: LDLT.h:670
Eigen::internal::is_lvalue
Definition: XprHelper.h:659
Eigen::SelfAdjointView::outerStride
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index outerStride() const EIGEN_NOEXCEPT
Definition: SelfAdjointView.h:84
Eigen::SelfAdjointView::llt
const LLT< PlainObject, UpLo > llt() const
Definition: LLT.h:551
Eigen::internal::evaluator_traits< SelfAdjointView< MatrixType, Mode > >::Kind
storage_kind_to_evaluator_kind< typename MatrixType::StorageKind >::Kind Kind
Definition: SelfAdjointView.h:287
EIGEN_CONSTEXPR
#define EIGEN_CONSTEXPR
Definition: Macros.h:787
Eigen::Upper
@ Upper
Definition: Constants.h:211
type
Definition: pytypes.h:1525
Eigen::MatrixBase::selfadjointView
EIGEN_DEVICE_FUNC SelfAdjointViewReturnType< UpLo >::Type selfadjointView()
Eigen::internal::evaluator_traits
Definition: CoreEvaluators.h:79
Eigen::SelfAdjointView::innerStride
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index innerStride() const EIGEN_NOEXCEPT
Definition: SelfAdjointView.h:86
Eigen::SelfAdjointView::PlainObject
MatrixType::PlainObject PlainObject
Definition: SelfAdjointView.h:71
Eigen::SelfAdjointView
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Definition: SelfAdjointView.h:49
Eigen::SelfAdjointView::StorageIndex
MatrixType::StorageIndex StorageIndex
Definition: SelfAdjointView.h:62
Eigen::SelfAdjointView::eigenvalues
EIGEN_DEVICE_FUNC EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Definition: MatrixBaseEigenvalues.h:88
Eigen::SelfAdjointView::TransposeMode
@ TransposeMode
Definition: SelfAdjointView.h:69
mat
MatrixXf mat
Definition: Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
Eigen::internal::copy_using_evaluator_traits
Definition: AssignEvaluator.h:28
Eigen::SelfAdjointView::adjoint
const EIGEN_DEVICE_FUNC AdjointReturnType adjoint() const
Definition: SelfAdjointView.h:216
Eigen::internal::generic_dense_assignment_kernel::SrcEvaluatorType
SrcEvaluatorTypeT SrcEvaluatorType
Definition: AssignEvaluator.h:628
eigen_internal_assert
#define eigen_internal_assert(x)
Definition: Macros.h:1043
Eigen::DirectAccessBit
const unsigned int DirectAccessBit
Definition: Constants.h:155
Eigen::internal::generic_dense_assignment_kernel::DstEvaluatorType
DstEvaluatorTypeT DstEvaluatorType
Definition: AssignEvaluator.h:627
Eigen::SelfAdjointView::Flags
@ Flags
Definition: SelfAdjointView.h:68
Eigen::PacketAccessBit
const unsigned int PacketAccessBit
Definition: Constants.h:94
Eigen::internal::generic_dense_assignment_kernel::Scalar
DstEvaluatorType::Scalar Scalar
Definition: AssignEvaluator.h:629
Eigen::SelfAdjointView::AdjointReturnType
SelfAdjointView< const typename MatrixType::AdjointReturnType, TransposeMode > AdjointReturnType
Definition: SelfAdjointView.h:213
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::assignOppositeCoeff
EIGEN_DEVICE_FUNC void assignOppositeCoeff(Index, Index)
Definition: SelfAdjointView.h:327
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::Base
generic_dense_assignment_kernel< DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version > Base
Definition: SelfAdjointView.h:296
Eigen::SelfAdjointView::_expression
const EIGEN_DEVICE_FUNC MatrixTypeNestedCleaned & _expression() const
Definition: SelfAdjointView.h:111
Eigen::SelfAdjointView::TransposeReturnType
SelfAdjointView< typename MatrixType::TransposeReturnType, TransposeMode > TransposeReturnType
Definition: SelfAdjointView.h:219
Eigen::internal::traits< SelfAdjointView< MatrixType, UpLo > >::MatrixTypeNestedCleaned
remove_all< MatrixTypeNested >::type MatrixTypeNestedCleaned
Definition: SelfAdjointView.h:36
Eigen::SelfAdjointView::Scalar
internal::traits< SelfAdjointView >::Scalar Scalar
The type of coefficients in this matrix.
Definition: SelfAdjointView.h:61
Eigen::SelfAdjointView::MatrixTypeNestedCleaned
internal::traits< SelfAdjointView >::MatrixTypeNestedCleaned MatrixTypeNestedCleaned
Definition: SelfAdjointView.h:57
Eigen::internal::true_type
Definition: Meta.h:96
Eigen::SelfAdjointView::operatorNorm
EIGEN_DEVICE_FUNC RealScalar operatorNorm() const
Computes the L2 operator norm.
Definition: MatrixBaseEigenvalues.h:151
Functor
Definition: NonLinearOptimization.cpp:117
Eigen::SelfAdjointView::Mode
@ Mode
Definition: SelfAdjointView.h:67
Eigen::LvalueBit
const unsigned int LvalueBit
Definition: Constants.h:144
Eigen::SelfAdjointView::transpose
const EIGEN_DEVICE_FUNC ConstTransposeReturnType transpose() const
Definition: SelfAdjointView.h:232
Eigen::TriangularBase
Base class for triangular part in a matrix.
Definition: TriangularMatrix.h:27
Eigen::SelfAdjointView::SelfAdjointView
EIGEN_DEVICE_FUNC SelfAdjointView(MatrixType &matrix)
Definition: SelfAdjointView.h:74
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::AssignmentTraits
Base::AssignmentTraits AssignmentTraits
Definition: SelfAdjointView.h:307
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::Scalar
Base::Scalar Scalar
Definition: SelfAdjointView.h:306
Eigen::internal::generic_dense_assignment_kernel::SrcXprType
SrcEvaluatorTypeT::XprType SrcXprType
Definition: AssignEvaluator.h:624
Eigen::Architecture::Type
Type
Definition: Constants.h:471
EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE
#define EIGEN_SCALAR_BINARYOP_EXPR_RETURN_TYPE(SCALAR, EXPR, OPNAME)
Definition: Macros.h:1351
Eigen::internal::generic_dense_assignment_kernel::m_src
const SrcEvaluatorType & m_src
Definition: AssignEvaluator.h:719
Eigen::SelfAdjointView::transpose
EIGEN_DEVICE_FUNC TransposeReturnType transpose()
Definition: SelfAdjointView.h:222
Eigen::Product
Expression of the product of two arbitrary matrices or vectors.
Definition: Product.h:71
Eigen::SelfAdjointView::cols
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: SelfAdjointView.h:82
Eigen::LDLT
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:59
Eigen::Triplet< double >
Eigen::SelfAdjointView::nestedExpression
const EIGEN_DEVICE_FUNC MatrixTypeNestedCleaned & nestedExpression() const
Definition: SelfAdjointView.h:114
Eigen::SelfAdjointView::nestedExpression
EIGEN_DEVICE_FUNC MatrixTypeNestedCleaned & nestedExpression()
Definition: SelfAdjointView.h:116
conj
AnnoyingScalar conj(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:104
Eigen::Lower
@ Lower
Definition: Constants.h:209
Eigen::internal::ref_selector
Definition: XprHelper.h:416
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::SrcEvaluatorType
Base::SrcEvaluatorType SrcEvaluatorType
Definition: SelfAdjointView.h:305
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::assignDiagonalCoeff
EIGEN_DEVICE_FUNC void assignDiagonalCoeff(Index id)
Definition: SelfAdjointView.h:322
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
Eigen::SelfAdjointView::MatrixTypeNested
internal::traits< SelfAdjointView >::MatrixTypeNested MatrixTypeNested
Definition: SelfAdjointView.h:56
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::assignCoeff
EIGEN_DEVICE_FUNC void assignCoeff(Index row, Index col)
Definition: SelfAdjointView.h:314
Eigen::LinearAccessBit
const unsigned int LinearAccessBit
Definition: Constants.h:130
Eigen::LLT
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:66
Eigen::SelfAdjointView::ConstTransposeReturnType
SelfAdjointView< const typename MatrixType::ConstTransposeReturnType, TransposeMode > ConstTransposeReturnType
Definition: SelfAdjointView.h:229
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::triangular_dense_assignment_kernel
EIGEN_DEVICE_FUNC triangular_dense_assignment_kernel(DstEvaluatorType &dst, const SrcEvaluatorType &src, const Functor &func, DstXprType &dstExpr)
Definition: SelfAdjointView.h:310
Eigen::internal::traits
Definition: ForwardDeclarations.h:17
Eigen::SelfAdjointView::conjugate
const EIGEN_DEVICE_FUNC ConjugateReturnType conjugate() const
Definition: SelfAdjointView.h:198
EIGEN_STATIC_ASSERT
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:127
row
m row(1)
Eigen::internal::conditional
Definition: Meta.h:109
Eigen::SelfAdjointShape
Definition: Constants.h:534
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
Eigen::internal::generic_dense_assignment_kernel::DstXprType
DstEvaluatorTypeT::XprType DstXprType
Definition: AssignEvaluator.h:623
product
void product(const MatrixType &m)
Definition: product.h:20
Eigen::SelfAdjointView::EigenvaluesReturnType
Matrix< RealScalar, internal::traits< MatrixType >::ColsAtCompileTime, 1 > EigenvaluesReturnType
Definition: SelfAdjointView.h:258
Eigen::SelfAdjointView::operator*
friend const EIGEN_DEVICE_FUNC Product< OtherDerived, SelfAdjointView > operator*(const MatrixBase< OtherDerived > &lhs, const SelfAdjointView &rhs)
Definition: SelfAdjointView.h:131
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
Eigen::internal::triangular_dense_assignment_kernel
Definition: TriangularMatrix.h:758
EIGEN_NOEXCEPT
#define EIGEN_NOEXCEPT
Definition: Macros.h:1418
Eigen::internal::generic_dense_assignment_kernel::assignCoeff
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:654
internal
Definition: BandTriangularSolver.h:13
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Eigen::SelfAdjointView::m_matrix
MatrixTypeNested m_matrix
Definition: SelfAdjointView.h:266
Eigen::SelfAdjointView::rows
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: SelfAdjointView.h:80
EIGEN_STATIC_ASSERT_LVALUE
#define EIGEN_STATIC_ASSERT_LVALUE(Derived)
Definition: StaticAssert.h:202
Eigen::internal::evaluator_traits< SelfAdjointView< MatrixType, Mode > >::Shape
SelfAdjointShape Shape
Definition: SelfAdjointView.h:288
Eigen::SelfAdjointView::coeff
EIGEN_DEVICE_FUNC Scalar coeff(Index row, Index col) const
Definition: SelfAdjointView.h:92
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::DstXprType
Base::DstXprType DstXprType
Definition: SelfAdjointView.h:297
func
Definition: benchGeometry.cpp:23
Eigen::SelfAdjointView::coeffRef
EIGEN_DEVICE_FUNC Scalar & coeffRef(Index row, Index col)
Definition: SelfAdjointView.h:102
Eigen::HereditaryBits
const unsigned int HereditaryBits
Definition: Constants.h:195
Eigen::TriangularView
Expression of a triangular part in a matrix.
Definition: TriangularMatrix.h:187
Eigen::internal::IndexBased
Definition: Constants.h:542
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:232
Eigen::SelfAdjointView::operator*
const EIGEN_DEVICE_FUNC Product< SelfAdjointView, OtherDerived > operator*(const MatrixBase< OtherDerived > &rhs) const
Definition: SelfAdjointView.h:122
Eigen::SelfAdjointView::triangularView
EIGEN_DEVICE_FUNC internal::conditional<(TriMode &(Upper|Lower))==(UpLo &(Upper|Lower)), TriangularView< MatrixType, TriMode >, TriangularView< typename MatrixType::AdjointReturnType, TriMode > >::type triangularView() const
Definition: SelfAdjointView.h:186
Eigen::internal::triangular_dense_assignment_kernel< UpLo, SelfAdjoint, SetOpposite, DstEvaluatorTypeT, SrcEvaluatorTypeT, Functor, Version >::SrcXprType
Base::SrcXprType SrcXprType
Definition: SelfAdjointView.h:298
Eigen::TriangularBase::check_coordinates_internal
void check_coordinates_internal(Index, Index) const
Definition: TriangularMatrix.h:146
Eigen::internal::generic_dense_assignment_kernel::m_dst
DstEvaluatorType & m_dst
Definition: AssignEvaluator.h:718
Eigen::SelfAdjoint
@ SelfAdjoint
Definition: Constants.h:225
Eigen::SelfAdjointView::MatrixConjugateReturnType
internal::remove_all< typename MatrixType::ConjugateReturnType >::type MatrixConjugateReturnType
Definition: SelfAdjointView.h:63
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Eigen::SelfAdjointView::MatrixType
_MatrixType MatrixType
Definition: SelfAdjointView.h:54
Eigen::SelfAdjointView::rankUpdate
EIGEN_DEVICE_FUNC SelfAdjointView & rankUpdate(const MatrixBase< DerivedU > &u, const MatrixBase< DerivedV > &v, const Scalar &alpha=Scalar(1))
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:04:10