AutoDiffVector.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_AUTODIFF_VECTOR_H
11 #define EIGEN_AUTODIFF_VECTOR_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 /* \class AutoDiffScalar
18  * \brief A scalar type replacement with automatic differentation capability
19  *
20  * \param DerType the vector type used to store/represent the derivatives (e.g. Vector3f)
21  *
22  * This class represents a scalar value while tracking its respective derivatives.
23  *
24  * It supports the following list of global math function:
25  * - std::abs, std::sqrt, std::pow, std::exp, std::log, std::sin, std::cos,
26  * - internal::abs, internal::sqrt, numext::pow, internal::exp, internal::log, internal::sin, internal::cos,
27  * - internal::conj, internal::real, internal::imag, numext::abs2.
28  *
29  * AutoDiffScalar can be used as the scalar type of an Eigen::Matrix object. However,
30  * in that case, the expression template mechanism only occurs at the top Matrix level,
31  * while derivatives are computed right away.
32  *
33  */
34 template<typename ValueType, typename JacobianType>
36 {
37  public:
38  //typedef typename internal::traits<ValueType>::Scalar Scalar;
39  typedef typename internal::traits<ValueType>::Scalar BaseScalar;
43  typedef typename JacobianType::Index Index;
44 
45  inline AutoDiffVector() {}
46 
47  inline AutoDiffVector(const ValueType& values)
48  : m_values(values)
49  {
50  m_jacobian.setZero();
51  }
52 
53 
55  const CoeffType operator[] (Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
56 
58  const CoeffType operator() (Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
59 
61  const CoeffType coeffRef(Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
62 
63  Index size() const { return m_values.size(); }
64 
65  // FIXME here we could return an expression of the sum
66  Scalar sum() const { /*std::cerr << "sum \n\n";*/ /*std::cerr << m_jacobian.rowwise().sum() << "\n\n";*/ return Scalar(m_values.sum(), m_jacobian.rowwise().sum()); }
67 
68 
69  inline AutoDiffVector(const ValueType& values, const JacobianType& jac)
70  : m_values(values), m_jacobian(jac)
71  {}
72 
73  template<typename OtherValueType, typename OtherJacobianType>
75  : m_values(other.values()), m_jacobian(other.jacobian())
76  {}
77 
78  inline AutoDiffVector(const AutoDiffVector& other)
79  : m_values(other.values()), m_jacobian(other.jacobian())
80  {}
81 
82  template<typename OtherValueType, typename OtherJacobianType>
84  {
85  m_values = other.values();
86  m_jacobian = other.jacobian();
87  return *this;
88  }
89 
90  inline AutoDiffVector& operator=(const AutoDiffVector& other)
91  {
92  m_values = other.values();
93  m_jacobian = other.jacobian();
94  return *this;
95  }
96 
97  inline const ValueType& values() const { return m_values; }
98  inline ValueType& values() { return m_values; }
99 
100  inline const JacobianType& jacobian() const { return m_jacobian; }
101  inline JacobianType& jacobian() { return m_jacobian; }
102 
103  template<typename OtherValueType,typename OtherJacobianType>
104  inline const AutoDiffVector<
105  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,ValueType,OtherValueType>::Type,
106  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,JacobianType,OtherJacobianType>::Type >
108  {
109  return AutoDiffVector<
110  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,ValueType,OtherValueType>::Type,
111  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,JacobianType,OtherJacobianType>::Type >(
112  m_values + other.values(),
113  m_jacobian + other.jacobian());
114  }
115 
116  template<typename OtherValueType, typename OtherJacobianType>
117  inline AutoDiffVector&
119  {
120  m_values += other.values();
121  m_jacobian += other.jacobian();
122  return *this;
123  }
124 
125  template<typename OtherValueType,typename OtherJacobianType>
126  inline const AutoDiffVector<
127  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,ValueType,OtherValueType>::Type,
128  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,JacobianType,OtherJacobianType>::Type >
130  {
131  return AutoDiffVector<
132  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,ValueType,OtherValueType>::Type,
133  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,JacobianType,OtherJacobianType>::Type >(
134  m_values - other.values(),
135  m_jacobian - other.jacobian());
136  }
137 
138  template<typename OtherValueType, typename OtherJacobianType>
139  inline AutoDiffVector&
141  {
142  m_values -= other.values();
143  m_jacobian -= other.jacobian();
144  return *this;
145  }
146 
147  inline const AutoDiffVector<
148  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, ValueType>::Type,
149  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type >
150  operator-() const
151  {
152  return AutoDiffVector<
153  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, ValueType>::Type,
154  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type >(
155  -m_values,
156  -m_jacobian);
157  }
158 
159  inline const AutoDiffVector<
160  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
161  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>
162  operator*(const BaseScalar& other) const
163  {
164  return AutoDiffVector<
165  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
166  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >(
167  m_values * other,
168  m_jacobian * other);
169  }
170 
171  friend inline const AutoDiffVector<
172  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
173  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >
174  operator*(const Scalar& other, const AutoDiffVector& v)
175  {
176  return AutoDiffVector<
177  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
178  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >(
179  v.values() * other,
180  v.jacobian() * other);
181  }
182 
183 // template<typename OtherValueType,typename OtherJacobianType>
184 // inline const AutoDiffVector<
185 // CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
186 // CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
187 // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
188 // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >
189 // operator*(const AutoDiffVector<OtherValueType,OtherJacobianType>& other) const
190 // {
191 // return AutoDiffVector<
192 // CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
193 // CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
194 // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
195 // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >(
196 // m_values.cwise() * other.values(),
197 // (m_jacobian * other.values()) + (m_values * other.jacobian()));
198 // }
199 
200  inline AutoDiffVector& operator*=(const Scalar& other)
201  {
202  m_values *= other;
203  m_jacobian *= other;
204  return *this;
205  }
206 
207  template<typename OtherValueType,typename OtherJacobianType>
209  {
210  *this = *this * other;
211  return *this;
212  }
213 
214  protected:
215  ValueType m_values;
216  JacobianType m_jacobian;
217 
218 };
219 
220 }
221 
222 #endif // EIGEN_AUTODIFF_VECTOR_H
Array< int, Dynamic, 1 > v
int i
A scalar type replacement with automatic differentiation capability.
AutoDiffScalar< Matrix< BaseScalar, JacobianType::RowsAtCompileTime, 1 > > ActiveScalar
const AutoDiffVector< typename MakeCwiseUnaryOp< internal::scalar_multiple_op< Scalar >, ValueType >::Type, typename MakeCwiseUnaryOp< internal::scalar_multiple_op< Scalar >, JacobianType >::Type > operator*(const BaseScalar &other) const
AutoDiffVector & operator=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
JacobianType & jacobian()
AutoDiffVector & operator+=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
AutoDiffVector & operator=(const AutoDiffVector &other)
AutoDiffScalar< typename JacobianType::ColXpr > CoeffType
CoeffType coeffRef(Index i)
AutoDiffVector(const ValueType &values)
CoeffType operator[](Index i)
JacobianType::Index Index
const AutoDiffVector< typename MakeCwiseUnaryOp< internal::scalar_opposite_op< Scalar >, ValueType >::Type, typename MakeCwiseUnaryOp< internal::scalar_opposite_op< Scalar >, JacobianType >::Type > operator-() const
AutoDiffVector(const ValueType &values, const JacobianType &jac)
internal::traits< ValueType >::Scalar BaseScalar
const CoeffType coeffRef(Index i) const
AutoDiffVector(const AutoDiffVector &other)
AutoDiffVector & operator-=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
const AutoDiffVector< typename MakeCwiseBinaryOp< internal::scalar_difference_op< Scalar >, ValueType, OtherValueType >::Type, typename MakeCwiseBinaryOp< internal::scalar_difference_op< Scalar >, JacobianType, OtherJacobianType >::Type > operator-(const AutoDiffVector< OtherValueType, OtherJacobianType > &other) const
const ValueType & values() const
AutoDiffVector & operator*=(const Scalar &other)
AutoDiffVector & operator*=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
AutoDiffVector(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
friend const AutoDiffVector< typename MakeCwiseUnaryOp< internal::scalar_multiple_op< Scalar >, ValueType >::Type, typename MakeCwiseUnaryOp< internal::scalar_multiple_op< Scalar >, JacobianType >::Type > operator*(const Scalar &other, const AutoDiffVector &v)
const JacobianType & jacobian() const
const AutoDiffVector< typename MakeCwiseBinaryOp< internal::scalar_sum_op< BaseScalar >, ValueType, OtherValueType >::Type, typename MakeCwiseBinaryOp< internal::scalar_sum_op< BaseScalar >, JacobianType, OtherJacobianType >::Type > operator+(const AutoDiffVector< OtherValueType, OtherJacobianType > &other) const
CoeffType operator()(Index i)
ValueType & values()
Scalar sum() const
Type
: TensorContractionSycl.h, provides various tensor contraction kernel for SYCL backend