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 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 /* \class AutoDiffScalar
19  * \brief A scalar type replacement with automatic differentation capability
20  *
21  * \param DerType the vector type used to store/represent the derivatives (e.g. Vector3f)
22  *
23  * This class represents a scalar value while tracking its respective derivatives.
24  *
25  * It supports the following list of global math function:
26  * - std::abs, std::sqrt, std::pow, std::exp, std::log, std::sin, std::cos,
27  * - internal::abs, internal::sqrt, numext::pow, internal::exp, internal::log, internal::sin, internal::cos,
28  * - internal::conj, internal::real, internal::imag, numext::abs2.
29  *
30  * AutoDiffScalar can be used as the scalar type of an Eigen::Matrix object. However,
31  * in that case, the expression template mechanism only occurs at the top Matrix level,
32  * while derivatives are computed right away.
33  *
34  */
35 template <typename ValueType, typename JacobianType>
37  public:
38  // typedef typename internal::traits<ValueType>::Scalar Scalar;
43  typedef typename JacobianType::Index Index;
44 
45  inline AutoDiffVector() {}
46 
47  inline AutoDiffVector(const ValueType& values) : m_values(values) { m_jacobian.setZero(); }
48 
50  const CoeffType operator[](Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
51 
53  const CoeffType operator()(Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
54 
56  const CoeffType coeffRef(Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
57 
58  Index size() const { return m_values.size(); }
59 
60  // FIXME here we could return an expression of the sum
61  Scalar sum() const { /*std::cerr << "sum \n\n";*/ /*std::cerr << m_jacobian.rowwise().sum() << "\n\n";*/
62  return Scalar(m_values.sum(), m_jacobian.rowwise().sum());
63  }
64 
65  inline AutoDiffVector(const ValueType& values, const JacobianType& jac) : m_values(values), m_jacobian(jac) {}
66 
67  template <typename OtherValueType, typename OtherJacobianType>
69  : m_values(other.values()), m_jacobian(other.jacobian()) {}
70 
71  inline AutoDiffVector(const AutoDiffVector& other) : m_values(other.values()), m_jacobian(other.jacobian()) {}
72 
73  template <typename OtherValueType, typename OtherJacobianType>
75  m_values = other.values();
76  m_jacobian = other.jacobian();
77  return *this;
78  }
79 
80  inline AutoDiffVector& operator=(const AutoDiffVector& other) {
81  m_values = other.values();
82  m_jacobian = other.jacobian();
83  return *this;
84  }
85 
86  inline const ValueType& values() const { return m_values; }
87  inline ValueType& values() { return m_values; }
88 
89  inline const JacobianType& jacobian() const { return m_jacobian; }
90  inline JacobianType& jacobian() { return m_jacobian; }
91 
92  template <typename OtherValueType, typename OtherJacobianType>
93  inline const AutoDiffVector<
94  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>, ValueType, OtherValueType>::Type,
95  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>, JacobianType, OtherJacobianType>::Type>
97  return AutoDiffVector<
98  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>, ValueType, OtherValueType>::Type,
99  typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>, JacobianType, OtherJacobianType>::Type>(
100  m_values + other.values(), m_jacobian + other.jacobian());
101  }
102 
103  template <typename OtherValueType, typename OtherJacobianType>
105  m_values += other.values();
106  m_jacobian += other.jacobian();
107  return *this;
108  }
109 
110  template <typename OtherValueType, typename OtherJacobianType>
111  inline const AutoDiffVector<
112  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>, ValueType, OtherValueType>::Type,
113  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>, JacobianType, OtherJacobianType>::Type>
115  return AutoDiffVector<
116  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>, ValueType, OtherValueType>::Type,
117  typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>, JacobianType, OtherJacobianType>::Type>(
118  m_values - other.values(), m_jacobian - other.jacobian());
119  }
120 
121  template <typename OtherValueType, typename OtherJacobianType>
123  m_values -= other.values();
124  m_jacobian -= other.jacobian();
125  return *this;
126  }
127 
129  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type>
130  operator-() const {
132  typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type>(
133  -m_values, -m_jacobian);
134  }
135 
137  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>
138  operator*(const BaseScalar& other) const {
140  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>(
141  m_values * other, m_jacobian * other);
142  }
143 
144  friend inline const AutoDiffVector<
145  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
146  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>
147  operator*(const Scalar& other, const AutoDiffVector& v) {
149  typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>(
150  v.values() * other, v.jacobian() * other);
151  }
152 
153  // template<typename OtherValueType,typename OtherJacobianType>
154  // inline const AutoDiffVector<
155  // CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
156  // CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
157  // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
158  // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >
159  // operator*(const AutoDiffVector<OtherValueType,OtherJacobianType>& other) const
160  // {
161  // return AutoDiffVector<
162  // CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
163  // CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
164  // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
165  // CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >(
166  // m_values.cwise() * other.values(),
167  // (m_jacobian * other.values()) + (m_values * other.jacobian()));
168  // }
169 
170  inline AutoDiffVector& operator*=(const Scalar& other) {
171  m_values *= other;
172  m_jacobian *= other;
173  return *this;
174  }
175 
176  template <typename OtherValueType, typename OtherJacobianType>
178  *this = *this * other;
179  return *this;
180  }
181 
182  protected:
183  ValueType m_values;
184  JacobianType m_jacobian;
185 };
186 
187 } // namespace Eigen
188 
189 #endif // EIGEN_AUTODIFF_VECTOR_H
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
A scalar type replacement with automatic differentiation capability.
Definition: AutoDiffScalar.h:99
Definition: AutoDiffVector.h:36
AutoDiffVector & operator=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
Definition: AutoDiffVector.h:74
JacobianType & jacobian()
Definition: AutoDiffVector.h:90
Index size() const
Definition: AutoDiffVector.h:58
AutoDiffVector & operator+=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
Definition: AutoDiffVector.h:104
ValueType m_values
Definition: AutoDiffVector.h:183
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)
Definition: AutoDiffVector.h:147
const AutoDiffVector< typename MakeCwiseUnaryOp< internal::scalar_opposite_op< Scalar >, ValueType >::Type, typename MakeCwiseUnaryOp< internal::scalar_opposite_op< Scalar >, JacobianType >::Type > operator-() const
Definition: AutoDiffVector.h:130
AutoDiffVector & operator=(const AutoDiffVector &other)
Definition: AutoDiffVector.h:80
AutoDiffScalar< typename JacobianType::ColXpr > CoeffType
Definition: AutoDiffVector.h:42
CoeffType coeffRef(Index i)
Definition: AutoDiffVector.h:55
AutoDiffVector(const ValueType &values)
Definition: AutoDiffVector.h:47
CoeffType operator[](Index i)
Definition: AutoDiffVector.h:49
JacobianType::Index Index
Definition: AutoDiffVector.h:43
const CoeffType operator[](Index i) const
Definition: AutoDiffVector.h:50
AutoDiffVector()
Definition: AutoDiffVector.h:45
AutoDiffVector(const ValueType &values, const JacobianType &jac)
Definition: AutoDiffVector.h:65
internal::traits< ValueType >::Scalar BaseScalar
Definition: AutoDiffVector.h:39
const CoeffType coeffRef(Index i) const
Definition: AutoDiffVector.h:56
const CoeffType operator()(Index i) const
Definition: AutoDiffVector.h:53
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
Definition: AutoDiffVector.h:96
AutoDiffVector(const AutoDiffVector &other)
Definition: AutoDiffVector.h:71
AutoDiffVector & operator-=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
Definition: AutoDiffVector.h:122
AutoDiffScalar< Matrix< BaseScalar, JacobianType::RowsAtCompileTime, 1 > > ActiveScalar
Definition: AutoDiffVector.h:40
const ValueType & values() const
Definition: AutoDiffVector.h:86
AutoDiffVector & operator*=(const Scalar &other)
Definition: AutoDiffVector.h:170
AutoDiffVector & operator*=(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
Definition: AutoDiffVector.h:177
AutoDiffVector(const AutoDiffVector< OtherValueType, OtherJacobianType > &other)
Definition: AutoDiffVector.h:68
const JacobianType & jacobian() const
Definition: AutoDiffVector.h:89
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
Definition: AutoDiffVector.h:114
CoeffType operator()(Index i)
Definition: AutoDiffVector.h:52
ValueType & values()
Definition: AutoDiffVector.h:87
ActiveScalar Scalar
Definition: AutoDiffVector.h:41
JacobianType m_jacobian
Definition: AutoDiffVector.h:184
Scalar sum() const
Definition: AutoDiffVector.h:61
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
Definition: AutoDiffVector.h:138
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
Type
Type of JSON value.
Definition: rapidjson.h:513
Definition: ForwardDeclarations.h:21