EigenBase.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 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_EIGENBASE_H
12 #define EIGEN_EIGENBASE_H
13 
14 namespace Eigen {
15 
26 template<typename Derived> struct EigenBase
27 {
28 // typedef typename internal::plain_matrix_type<Derived>::type PlainObject;
29 
32 
34  Derived& derived() { return *static_cast<Derived*>(this); }
36  const Derived& derived() const { return *static_cast<const Derived*>(this); }
37 
38  inline Derived& const_cast_derived() const
39  { return *static_cast<Derived*>(const_cast<EigenBase*>(this)); }
40  inline const Derived& const_derived() const
41  { return *static_cast<const Derived*>(this); }
42 
44  inline Index rows() const { return derived().rows(); }
46  inline Index cols() const { return derived().cols(); }
49  inline Index size() const { return rows() * cols(); }
50 
52  template<typename Dest> inline void evalTo(Dest& dst) const
53  { derived().evalTo(dst); }
54 
56  template<typename Dest> inline void addTo(Dest& dst) const
57  {
58  // This is the default implementation,
59  // derived class can reimplement it in a more optimized way.
60  typename Dest::PlainObject res(rows(),cols());
61  evalTo(res);
62  dst += res;
63  }
64 
66  template<typename Dest> inline void subTo(Dest& dst) const
67  {
68  // This is the default implementation,
69  // derived class can reimplement it in a more optimized way.
70  typename Dest::PlainObject res(rows(),cols());
71  evalTo(res);
72  dst -= res;
73  }
74 
76  template<typename Dest> inline void applyThisOnTheRight(Dest& dst) const
77  {
78  // This is the default implementation,
79  // derived class can reimplement it in a more optimized way.
80  dst = dst * this->derived();
81  }
82 
84  template<typename Dest> inline void applyThisOnTheLeft(Dest& dst) const
85  {
86  // This is the default implementation,
87  // derived class can reimplement it in a more optimized way.
88  dst = this->derived() * dst;
89  }
90 
91 };
92 
93 /***************************************************************************
94 * Implementation of matrix base methods
95 ***************************************************************************/
96 
105 template<typename Derived>
106 template<typename OtherDerived>
108 {
109  other.derived().evalTo(derived());
110  return derived();
111 }
112 
113 template<typename Derived>
114 template<typename OtherDerived>
116 {
117  other.derived().addTo(derived());
118  return derived();
119 }
120 
121 template<typename Derived>
122 template<typename OtherDerived>
124 {
125  other.derived().subTo(derived());
126  return derived();
127 }
128 
133 template<typename Derived>
134 template<typename OtherDerived>
135 inline Derived&
137 {
138  other.derived().applyThisOnTheRight(derived());
139  return derived();
140 }
141 
144 template<typename Derived>
145 template<typename OtherDerived>
147 {
148  other.derived().applyThisOnTheRight(derived());
149 }
150 
152 template<typename Derived>
153 template<typename OtherDerived>
155 {
156  other.derived().applyThisOnTheLeft(derived());
157 }
158 
159 } // end namespace Eigen
160 
161 #endif // EIGEN_EIGENBASE_H
const Derived & derived() const
Definition: EigenBase.h:36
Derived & operator-=(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:123
void applyThisOnTheLeft(Dest &dst) const
Definition: EigenBase.h:84
Derived & const_cast_derived() const
Definition: EigenBase.h:38
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:154
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:146
Definition: LDLT.h:16
void applyThisOnTheRight(Dest &dst) const
Definition: EigenBase.h:76
Derived & derived()
Definition: EigenBase.h:34
const Derived & const_derived() const
Definition: EigenBase.h:40
Index rows() const
Definition: EigenBase.h:44
void addTo(Dest &dst) const
Definition: EigenBase.h:56
internal::traits< Derived >::Index Index
Definition: EigenBase.h:31
void evalTo(Dest &dst) const
Definition: EigenBase.h:52
Index cols() const
Definition: EigenBase.h:46
internal::traits< Derived >::StorageKind StorageKind
Definition: EigenBase.h:30
void subTo(Dest &dst) const
Definition: EigenBase.h:66
Derived & operator=(const DenseBase< OtherDerived > &other)
Derived & operator*=(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:136
Derived & operator+=(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:115
Index size() const
Definition: EigenBase.h:49


tuw_aruco
Author(s): Lukas Pfeifhofer
autogenerated on Mon Jun 10 2019 15:40:48