AutoDiffJacobian.h revision 2b8756b6f1de65d3f8bffab45be6c44ceb7411fc
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:
19  AutoDiffJacobian() : Functor() {}
20  AutoDiffJacobian(const Functor& f) : Functor(f) {}
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
44  typedef Matrix<Scalar, ValuesAtCompileTime, InputsAtCompileTime> JacobianType;
45  typedef typename JacobianType::Index Index;
46
47  typedef Matrix<Scalar, InputsAtCompileTime, 1> DerivativeType;
48  typedef AutoDiffScalar<DerivativeType> ActiveScalar;
49
50  typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
51  typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
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():
56  EIGEN_STRONG_INLINE
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
85    if(InputsAtCompileTime==Dynamic)
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
109