$darkmode
Eigen-unsupported  5.0.1-dev
AutoDiffJacobian.h
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 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 template <typename Functor>
19 class AutoDiffJacobian : public Functor {
20  public:
21  AutoDiffJacobian() : Functor() {}
22  AutoDiffJacobian(const Functor& f) : Functor(f) {}
23 
24  // forward constructors
25  template <typename... T>
26  AutoDiffJacobian(const T&... Values) : Functor(Values...) {}
27 
28  typedef typename Functor::InputType InputType;
29  typedef typename Functor::ValueType ValueType;
30  typedef typename ValueType::Scalar Scalar;
31 
32  enum { InputsAtCompileTime = InputType::RowsAtCompileTime, ValuesAtCompileTime = ValueType::RowsAtCompileTime };
33 
34  typedef Matrix<Scalar, ValuesAtCompileTime, InputsAtCompileTime> JacobianType;
35  typedef typename JacobianType::Index Index;
36 
37  typedef Matrix<Scalar, InputsAtCompileTime, 1> DerivativeType;
38  typedef AutoDiffScalar<DerivativeType> ActiveScalar;
39 
40  typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
41  typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
42 
43  // Some compilers don't accept variadic parameters after a default parameter,
44  // i.e., we can't just write _jac=0 but we need to overload operator():
45  EIGEN_STRONG_INLINE void operator()(const InputType& x, ValueType* v) const { this->operator()(x, v, 0); }
46  template <typename... ParamsType>
47  void operator()(const InputType& x, ValueType* v, JacobianType* _jac, const ParamsType&... Params) const {
48  eigen_assert(v != 0);
49 
50  if (!_jac) {
51  Functor::operator()(x, v, Params...);
52  return;
53  }
54 
55  JacobianType& jac = *_jac;
56 
57  ActiveInput ax = x.template cast<ActiveScalar>();
58  ActiveValue av(jac.rows());
59 
60  if (InputsAtCompileTime == Dynamic)
61  for (Index j = 0; j < jac.rows(); j++) av[j].derivatives().resize(x.rows());
62 
63  for (Index i = 0; i < jac.cols(); i++) ax[i].derivatives() = DerivativeType::Unit(x.rows(), i);
64 
65  Functor::operator()(ax, &av, Params...);
66 
67  for (Index i = 0; i < jac.rows(); i++) {
68  (*v)[i] = av[i].value();
69  jac.row(i) = av[i].derivatives();
70  }
71  }
72 };
73 
74 } // namespace Eigen
75 
76 #endif // EIGEN_AUTODIFF_JACOBIAN_H
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const int Dynamic