TR-mbed 1.0
Loading...
Searching...
No Matches
AutoDiffJacobian.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_JACOBIAN_H
11#define EIGEN_AUTODIFF_JACOBIAN_H
12
13namespace Eigen
14{
15
16template<typename Functor> class AutoDiffJacobian : public Functor
17{
18public:
21
22 // forward constructors
23#if EIGEN_HAS_VARIADIC_TEMPLATES
24 template<typename... T>
25 AutoDiffJacobian(const T& ...Values) : Functor(Values...) {}
26#else
27 template<typename T0>
28 AutoDiffJacobian(const T0& a0) : Functor(a0) {}
29 template<typename T0, typename T1>
30 AutoDiffJacobian(const T0& a0, const T1& a1) : Functor(a0, a1) {}
31 template<typename T0, typename T1, typename T2>
32 AutoDiffJacobian(const T0& a0, const T1& a1, const T2& a2) : Functor(a0, a1, a2) {}
33#endif
34
35 typedef typename Functor::InputType InputType;
36 typedef typename Functor::ValueType ValueType;
37 typedef typename ValueType::Scalar Scalar;
38
39 enum {
40 InputsAtCompileTime = InputType::RowsAtCompileTime,
41 ValuesAtCompileTime = ValueType::RowsAtCompileTime
42 };
43
45 typedef typename JacobianType::Index Index;
46
49
52
53#if EIGEN_HAS_VARIADIC_TEMPLATES
54 // Some compilers don't accept variadic parameters after a default parameter,
55 // i.e., we can't just write _jac=0 but we need to overload operator():
57 void operator() (const InputType& x, ValueType* v) const
58 {
59 this->operator()(x, v, 0);
60 }
61 template<typename... ParamsType>
62 void operator() (const InputType& x, ValueType* v, JacobianType* _jac,
63 const ParamsType&... Params) const
64#else
65 void operator() (const InputType& x, ValueType* v, JacobianType* _jac=0) const
66#endif
67 {
68 eigen_assert(v!=0);
69
70 if (!_jac)
71 {
72#if EIGEN_HAS_VARIADIC_TEMPLATES
73 Functor::operator()(x, v, Params...);
74#else
75 Functor::operator()(x, v);
76#endif
77 return;
78 }
79
80 JacobianType& jac = *_jac;
81
82 ActiveInput ax = x.template cast<ActiveScalar>();
83 ActiveValue av(jac.rows());
84
86 for (Index j=0; j<jac.rows(); j++)
87 av[j].derivatives().resize(x.rows());
88
89 for (Index i=0; i<jac.cols(); i++)
90 ax[i].derivatives() = DerivativeType::Unit(x.rows(),i);
91
92#if EIGEN_HAS_VARIADIC_TEMPLATES
93 Functor::operator()(ax, &av, Params...);
94#else
95 Functor::operator()(ax, &av);
96#endif
97
98 for (Index i=0; i<jac.rows(); i++)
99 {
100 (*v)[i] = av[i].value();
101 jac.row(i) = av[i].derivatives();
102 }
103 }
104};
105
106}
107
108#endif // EIGEN_AUTODIFF_JACOBIAN_H
Array< int, Dynamic, 1 > v
Definition Array_initializer_list_vector_cxx11.cpp:1
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define eigen_assert(x)
Definition Macros.h:1037
#define EIGEN_STRONG_INLINE
Definition Macros.h:917
Definition AutoDiffJacobian.h:17
Functor::InputType InputType
Definition AutoDiffJacobian.h:35
AutoDiffJacobian(const T0 &a0, const T1 &a1)
Definition AutoDiffJacobian.h:30
AutoDiffJacobian()
Definition AutoDiffJacobian.h:19
ValueType::Scalar Scalar
Definition AutoDiffJacobian.h:37
Functor::ValueType ValueType
Definition AutoDiffJacobian.h:36
void operator()(const InputType &x, ValueType *v, JacobianType *_jac=0) const
Definition AutoDiffJacobian.h:65
AutoDiffJacobian(const T0 &a0)
Definition AutoDiffJacobian.h:28
AutoDiffScalar< DerivativeType > ActiveScalar
Definition AutoDiffJacobian.h:48
Matrix< Scalar, InputsAtCompileTime, 1 > DerivativeType
Definition AutoDiffJacobian.h:47
AutoDiffJacobian(const T0 &a0, const T1 &a1, const T2 &a2)
Definition AutoDiffJacobian.h:32
Matrix< Scalar, ValuesAtCompileTime, InputsAtCompileTime > JacobianType
Definition AutoDiffJacobian.h:44
@ ValuesAtCompileTime
Definition AutoDiffJacobian.h:41
@ InputsAtCompileTime
Definition AutoDiffJacobian.h:40
Matrix< ActiveScalar, ValuesAtCompileTime, 1 > ActiveValue
Definition AutoDiffJacobian.h:51
JacobianType::Index Index
Definition AutoDiffJacobian.h:45
Matrix< ActiveScalar, InputsAtCompileTime, 1 > ActiveInput
Definition AutoDiffJacobian.h:50
AutoDiffJacobian(const Functor &f)
Definition AutoDiffJacobian.h:20
A scalar type replacement with automatic differentiation capability.
Definition AutoDiffScalar.h:71
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition PlainObjectBase.h:145
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:271
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition PlainObjectBase.h:143
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
const int Dynamic
Definition Constants.h:22
std::vector< float > Values
Definition sparse_setter.cpp:45
Definition NonLinearOptimization.cpp:118
Matrix< Scalar, InputsAtCompileTime, 1 > InputType
Definition NonLinearOptimization.cpp:124
Matrix< Scalar, ValuesAtCompileTime, 1 > ValueType
Definition NonLinearOptimization.cpp:125
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2