SparseTranspose.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) 2008-2015 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_SPARSETRANSPOSE_H
11 #define EIGEN_SPARSETRANSPOSE_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18  template<typename MatrixType,int CompressedAccess=int(MatrixType::Flags&CompressedAccessBit)>
19  class SparseTransposeImpl
20  : public SparseMatrixBase<Transpose<MatrixType> >
21  {};
22 
23  template<typename MatrixType>
24  class SparseTransposeImpl<MatrixType,CompressedAccessBit>
25  : public SparseCompressedBase<Transpose<MatrixType> >
26  {
27  typedef SparseCompressedBase<Transpose<MatrixType> > Base;
28  public:
29  using Base::derived;
30  typedef typename Base::Scalar Scalar;
31  typedef typename Base::StorageIndex StorageIndex;
32 
33  inline Index nonZeros() const { return derived().nestedExpression().nonZeros(); }
34 
35  inline const Scalar* valuePtr() const { return derived().nestedExpression().valuePtr(); }
36  inline const StorageIndex* innerIndexPtr() const { return derived().nestedExpression().innerIndexPtr(); }
37  inline const StorageIndex* outerIndexPtr() const { return derived().nestedExpression().outerIndexPtr(); }
38  inline const StorageIndex* innerNonZeroPtr() const { return derived().nestedExpression().innerNonZeroPtr(); }
39 
40  inline Scalar* valuePtr() { return derived().nestedExpression().valuePtr(); }
41  inline StorageIndex* innerIndexPtr() { return derived().nestedExpression().innerIndexPtr(); }
42  inline StorageIndex* outerIndexPtr() { return derived().nestedExpression().outerIndexPtr(); }
43  inline StorageIndex* innerNonZeroPtr() { return derived().nestedExpression().innerNonZeroPtr(); }
44  };
45 }
46 
47 template<typename MatrixType> class TransposeImpl<MatrixType,Sparse>
48  : public internal::SparseTransposeImpl<MatrixType>
49 {
50  protected:
51  typedef internal::SparseTransposeImpl<MatrixType> Base;
52 };
53 
54 namespace internal {
55 
56 template<typename ArgType>
57 struct unary_evaluator<Transpose<ArgType>, IteratorBased>
58  : public evaluator_base<Transpose<ArgType> >
59 {
60  typedef typename evaluator<ArgType>::InnerIterator EvalIterator;
61  public:
62  typedef Transpose<ArgType> XprType;
63 
64  inline Index nonZerosEstimate() const {
65  return m_argImpl.nonZerosEstimate();
66  }
67 
68  class InnerIterator : public EvalIterator
69  {
70  public:
71  EIGEN_STRONG_INLINE InnerIterator(const unary_evaluator& unaryOp, Index outer)
72  : EvalIterator(unaryOp.m_argImpl,outer)
73  {}
74 
75  Index row() const { return EvalIterator::col(); }
76  Index col() const { return EvalIterator::row(); }
77  };
78 
79  enum {
80  CoeffReadCost = evaluator<ArgType>::CoeffReadCost,
81  Flags = XprType::Flags
82  };
83 
84  explicit unary_evaluator(const XprType& op) :m_argImpl(op.nestedExpression()) {}
85 
86  protected:
87  evaluator<ArgType> m_argImpl;
88 };
89 
90 } // end namespace internal
91 
92 } // end namespace Eigen
93 
94 #endif // EIGEN_SPARSETRANSPOSE_H
RowXpr row(Index i)
This is the const version of row(). *‍/.
ColXpr col(Index i)
This is the const version of col().
Matrix< float, 1, Dynamic > MatrixType
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
internal::traits< Transpose< MatrixType > >::StorageIndex StorageIndex
internal::traits< Transpose< MatrixType > >::Scalar Scalar
internal::SparseTransposeImpl< MatrixType > Base
Expression of the transpose of a matrix.
Definition: Transpose.h:56
const unsigned int CompressedAccessBit
Definition: Constants.h:193
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
Transpose< MatrixType > & derived()
Definition: EigenBase.h:48