Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:02

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_AUTODIFF_JACOBIAN_H
0011 #define EIGEN_AUTODIFF_JACOBIAN_H
0012 
0013 namespace Eigen
0014 {
0015 
0016 template<typename Functor> class AutoDiffJacobian : public Functor
0017 {
0018 public:
0019   AutoDiffJacobian() : Functor() {}
0020   AutoDiffJacobian(const Functor& f) : Functor(f) {}
0021 
0022   // forward constructors
0023 #if EIGEN_HAS_VARIADIC_TEMPLATES
0024   template<typename... T>
0025   AutoDiffJacobian(const T& ...Values) : Functor(Values...) {}
0026 #else
0027   template<typename T0>
0028   AutoDiffJacobian(const T0& a0) : Functor(a0) {}
0029   template<typename T0, typename T1>
0030   AutoDiffJacobian(const T0& a0, const T1& a1) : Functor(a0, a1) {}
0031   template<typename T0, typename T1, typename T2>
0032   AutoDiffJacobian(const T0& a0, const T1& a1, const T2& a2) : Functor(a0, a1, a2) {}
0033 #endif
0034 
0035   typedef typename Functor::InputType InputType;
0036   typedef typename Functor::ValueType ValueType;
0037   typedef typename ValueType::Scalar Scalar;
0038 
0039   enum {
0040     InputsAtCompileTime = InputType::RowsAtCompileTime,
0041     ValuesAtCompileTime = ValueType::RowsAtCompileTime
0042   };
0043 
0044   typedef Matrix<Scalar, ValuesAtCompileTime, InputsAtCompileTime> JacobianType;
0045   typedef typename JacobianType::Index Index;
0046 
0047   typedef Matrix<Scalar, InputsAtCompileTime, 1> DerivativeType;
0048   typedef AutoDiffScalar<DerivativeType> ActiveScalar;
0049 
0050   typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
0051   typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
0052 
0053 #if EIGEN_HAS_VARIADIC_TEMPLATES
0054   // Some compilers don't accept variadic parameters after a default parameter,
0055   // i.e., we can't just write _jac=0 but we need to overload operator():
0056   EIGEN_STRONG_INLINE
0057   void operator() (const InputType& x, ValueType* v) const
0058   {
0059       this->operator()(x, v, 0);
0060   }
0061   template<typename... ParamsType>
0062   void operator() (const InputType& x, ValueType* v, JacobianType* _jac,
0063                    const ParamsType&... Params) const
0064 #else
0065   void operator() (const InputType& x, ValueType* v, JacobianType* _jac=0) const
0066 #endif
0067   {
0068     eigen_assert(v!=0);
0069 
0070     if (!_jac)
0071     {
0072 #if EIGEN_HAS_VARIADIC_TEMPLATES
0073       Functor::operator()(x, v, Params...);
0074 #else
0075       Functor::operator()(x, v);
0076 #endif
0077       return;
0078     }
0079 
0080     JacobianType& jac = *_jac;
0081 
0082     ActiveInput ax = x.template cast<ActiveScalar>();
0083     ActiveValue av(jac.rows());
0084 
0085     if(InputsAtCompileTime==Dynamic)
0086       for (Index j=0; j<jac.rows(); j++)
0087         av[j].derivatives().resize(x.rows());
0088 
0089     for (Index i=0; i<jac.cols(); i++)
0090       ax[i].derivatives() = DerivativeType::Unit(x.rows(),i);
0091 
0092 #if EIGEN_HAS_VARIADIC_TEMPLATES
0093     Functor::operator()(ax, &av, Params...);
0094 #else
0095     Functor::operator()(ax, &av);
0096 #endif
0097 
0098     for (Index i=0; i<jac.rows(); i++)
0099     {
0100       (*v)[i] = av[i].value();
0101       jac.row(i) = av[i].derivatives();
0102     }
0103   }
0104 };
0105 
0106 }
0107 
0108 #endif // EIGEN_AUTODIFF_JACOBIAN_H