mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-01-06 14:14:46 +08:00
Made AutoDiffJacobian more intuitive to use and updated for C++11
Changes: * Removed unnecessary types from the Functor by inferring from its types * Removed inputs() function reference, replaced with .rows() * Updated the forward constructor to use variadic templates * Added optional parameters to the Fuctor for passing parameters, control signals, etc * Has been tested with fixed size and dynamic matricies Ammendment by chtz: overload operator() for compatibility with not fully conforming compilers
This commit is contained in:
parent
4adeababf9
commit
6edd2e2851
@ -20,37 +20,60 @@ public:
|
|||||||
AutoDiffJacobian(const Functor& f) : Functor(f) {}
|
AutoDiffJacobian(const Functor& f) : Functor(f) {}
|
||||||
|
|
||||||
// forward constructors
|
// forward constructors
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
template<typename... T>
|
||||||
|
AutoDiffJacobian(const T& ...Values) : Functor(Values...) {}
|
||||||
|
#else
|
||||||
template<typename T0>
|
template<typename T0>
|
||||||
AutoDiffJacobian(const T0& a0) : Functor(a0) {}
|
AutoDiffJacobian(const T0& a0) : Functor(a0) {}
|
||||||
template<typename T0, typename T1>
|
template<typename T0, typename T1>
|
||||||
AutoDiffJacobian(const T0& a0, const T1& a1) : Functor(a0, a1) {}
|
AutoDiffJacobian(const T0& a0, const T1& a1) : Functor(a0, a1) {}
|
||||||
template<typename T0, typename T1, typename T2>
|
template<typename T0, typename T1, typename T2>
|
||||||
AutoDiffJacobian(const T0& a0, const T1& a1, const T2& a2) : Functor(a0, a1, a2) {}
|
AutoDiffJacobian(const T0& a0, const T1& a1, const T2& a2) : Functor(a0, a1, a2) {}
|
||||||
|
#endif
|
||||||
enum {
|
|
||||||
InputsAtCompileTime = Functor::InputsAtCompileTime,
|
|
||||||
ValuesAtCompileTime = Functor::ValuesAtCompileTime
|
|
||||||
};
|
|
||||||
|
|
||||||
typedef typename Functor::InputType InputType;
|
typedef typename Functor::InputType InputType;
|
||||||
typedef typename Functor::ValueType ValueType;
|
typedef typename Functor::ValueType ValueType;
|
||||||
typedef typename Functor::JacobianType JacobianType;
|
typedef typename ValueType::Scalar Scalar;
|
||||||
typedef typename JacobianType::Scalar Scalar;
|
|
||||||
|
enum {
|
||||||
|
InputsAtCompileTime = InputType::RowsAtCompileTime,
|
||||||
|
ValuesAtCompileTime = ValueType::RowsAtCompileTime
|
||||||
|
};
|
||||||
|
|
||||||
|
typedef Matrix<Scalar, ValuesAtCompileTime, InputsAtCompileTime> JacobianType;
|
||||||
typedef typename JacobianType::Index Index;
|
typedef typename JacobianType::Index Index;
|
||||||
|
|
||||||
typedef Matrix<Scalar,InputsAtCompileTime,1> DerivativeType;
|
typedef Matrix<Scalar, InputsAtCompileTime, 1> DerivativeType;
|
||||||
typedef AutoDiffScalar<DerivativeType> ActiveScalar;
|
typedef AutoDiffScalar<DerivativeType> ActiveScalar;
|
||||||
|
|
||||||
|
|
||||||
typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
|
typedef Matrix<ActiveScalar, InputsAtCompileTime, 1> ActiveInput;
|
||||||
typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
|
typedef Matrix<ActiveScalar, ValuesAtCompileTime, 1> ActiveValue;
|
||||||
|
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
// Some compilers don't accept variadic parameters after a default parameter,
|
||||||
|
// i.e., we can't just write _jac=0 but we need to overload operator():
|
||||||
|
EIGEN_STRONG_INLINE
|
||||||
|
void operator() (const InputType& x, ValueType* v) const
|
||||||
|
{
|
||||||
|
this->operator()(x, v, 0);
|
||||||
|
}
|
||||||
|
template<typename... ParamsType>
|
||||||
|
void operator() (const InputType& x, ValueType* v, JacobianType* _jac,
|
||||||
|
const ParamsType&... Params) const
|
||||||
|
#else
|
||||||
void operator() (const InputType& x, ValueType* v, JacobianType* _jac=0) const
|
void operator() (const InputType& x, ValueType* v, JacobianType* _jac=0) const
|
||||||
|
#endif
|
||||||
{
|
{
|
||||||
eigen_assert(v!=0);
|
eigen_assert(v!=0);
|
||||||
|
|
||||||
if (!_jac)
|
if (!_jac)
|
||||||
{
|
{
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
Functor::operator()(x, v, Params...);
|
||||||
|
#else
|
||||||
Functor::operator()(x, v);
|
Functor::operator()(x, v);
|
||||||
|
#endif
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -61,12 +84,16 @@ public:
|
|||||||
|
|
||||||
if(InputsAtCompileTime==Dynamic)
|
if(InputsAtCompileTime==Dynamic)
|
||||||
for (Index j=0; j<jac.rows(); j++)
|
for (Index j=0; j<jac.rows(); j++)
|
||||||
av[j].derivatives().resize(this->inputs());
|
av[j].derivatives().resize(x.rows());
|
||||||
|
|
||||||
for (Index i=0; i<jac.cols(); i++)
|
for (Index i=0; i<jac.cols(); i++)
|
||||||
ax[i].derivatives() = DerivativeType::Unit(this->inputs(),i);
|
ax[i].derivatives() = DerivativeType::Unit(x.rows(),i);
|
||||||
|
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
Functor::operator()(ax, &av, Params...);
|
||||||
|
#else
|
||||||
Functor::operator()(ax, &av);
|
Functor::operator()(ax, &av);
|
||||||
|
#endif
|
||||||
|
|
||||||
for (Index i=0; i<jac.rows(); i++)
|
for (Index i=0; i<jac.rows(); i++)
|
||||||
{
|
{
|
||||||
@ -74,8 +101,6 @@ public:
|
|||||||
jac.row(i) = av[i].derivatives();
|
jac.row(i) = av[i].derivatives();
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
protected:
|
|
||||||
|
|
||||||
};
|
};
|
||||||
|
|
||||||
}
|
}
|
||||||
|
@ -105,6 +105,89 @@ struct TestFunc1
|
|||||||
}
|
}
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
/* Test functor for the C++11 features. */
|
||||||
|
template <typename Scalar>
|
||||||
|
struct integratorFunctor
|
||||||
|
{
|
||||||
|
typedef Matrix<Scalar, 2, 1> InputType;
|
||||||
|
typedef Matrix<Scalar, 2, 1> ValueType;
|
||||||
|
|
||||||
|
/*
|
||||||
|
* Implementation starts here.
|
||||||
|
*/
|
||||||
|
integratorFunctor(const Scalar gain) : _gain(gain) {}
|
||||||
|
integratorFunctor(const integratorFunctor& f) : _gain(f._gain) {}
|
||||||
|
const Scalar _gain;
|
||||||
|
|
||||||
|
template <typename T1, typename T2>
|
||||||
|
void operator() (const T1 &input, T2 *output, const Scalar dt) const
|
||||||
|
{
|
||||||
|
T2 &o = *output;
|
||||||
|
|
||||||
|
/* Integrator to test the AD. */
|
||||||
|
o[0] = input[0] + input[1] * dt * _gain;
|
||||||
|
o[1] = input[1] * _gain;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* Only needed for the test */
|
||||||
|
template <typename T1, typename T2, typename T3>
|
||||||
|
void operator() (const T1 &input, T2 *output, T3 *jacobian, const Scalar dt) const
|
||||||
|
{
|
||||||
|
T2 &o = *output;
|
||||||
|
|
||||||
|
/* Integrator to test the AD. */
|
||||||
|
o[0] = input[0] + input[1] * dt * _gain;
|
||||||
|
o[1] = input[1] * _gain;
|
||||||
|
|
||||||
|
if (jacobian)
|
||||||
|
{
|
||||||
|
T3 &j = *jacobian;
|
||||||
|
|
||||||
|
j(0, 0) = 1;
|
||||||
|
j(0, 1) = dt * _gain;
|
||||||
|
j(1, 0) = 0;
|
||||||
|
j(1, 1) = _gain;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
};
|
||||||
|
|
||||||
|
template<typename Func> void forward_jacobian_cpp11(const Func& f)
|
||||||
|
{
|
||||||
|
typedef typename Func::ValueType::Scalar Scalar;
|
||||||
|
typedef typename Func::ValueType ValueType;
|
||||||
|
typedef typename Func::InputType InputType;
|
||||||
|
typedef typename AutoDiffJacobian<Func>::JacobianType JacobianType;
|
||||||
|
|
||||||
|
InputType x = InputType::Random(InputType::RowsAtCompileTime);
|
||||||
|
ValueType y, yref;
|
||||||
|
JacobianType j, jref;
|
||||||
|
|
||||||
|
const Scalar dt = internal::random<double>();
|
||||||
|
|
||||||
|
jref.setZero();
|
||||||
|
yref.setZero();
|
||||||
|
f(x, &yref, &jref, dt);
|
||||||
|
|
||||||
|
//std::cerr << "y, yref, jref: " << "\n";
|
||||||
|
//std::cerr << y.transpose() << "\n\n";
|
||||||
|
//std::cerr << yref << "\n\n";
|
||||||
|
//std::cerr << jref << "\n\n";
|
||||||
|
|
||||||
|
AutoDiffJacobian<Func> autoj(f);
|
||||||
|
autoj(x, &y, &j, dt);
|
||||||
|
|
||||||
|
//std::cerr << "y j (via autodiff): " << "\n";
|
||||||
|
//std::cerr << y.transpose() << "\n\n";
|
||||||
|
//std::cerr << j << "\n\n";
|
||||||
|
|
||||||
|
VERIFY_IS_APPROX(y, yref);
|
||||||
|
VERIFY_IS_APPROX(j, jref);
|
||||||
|
}
|
||||||
|
#endif
|
||||||
|
|
||||||
template<typename Func> void forward_jacobian(const Func& f)
|
template<typename Func> void forward_jacobian(const Func& f)
|
||||||
{
|
{
|
||||||
typename Func::InputType x = Func::InputType::Random(f.inputs());
|
typename Func::InputType x = Func::InputType::Random(f.inputs());
|
||||||
@ -128,7 +211,6 @@ template<typename Func> void forward_jacobian(const Func& f)
|
|||||||
VERIFY_IS_APPROX(j, jref);
|
VERIFY_IS_APPROX(j, jref);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// TODO also check actual derivatives!
|
// TODO also check actual derivatives!
|
||||||
template <int>
|
template <int>
|
||||||
void test_autodiff_scalar()
|
void test_autodiff_scalar()
|
||||||
@ -141,6 +223,7 @@ void test_autodiff_scalar()
|
|||||||
VERIFY_IS_APPROX(res.value(), foo(p.x(),p.y()));
|
VERIFY_IS_APPROX(res.value(), foo(p.x(),p.y()));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// TODO also check actual derivatives!
|
// TODO also check actual derivatives!
|
||||||
template <int>
|
template <int>
|
||||||
void test_autodiff_vector()
|
void test_autodiff_vector()
|
||||||
@ -164,6 +247,9 @@ void test_autodiff_jacobian()
|
|||||||
CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,2>()) ));
|
CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,2>()) ));
|
||||||
CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,3>()) ));
|
CALL_SUBTEST(( forward_jacobian(TestFunc1<double,3,3>()) ));
|
||||||
CALL_SUBTEST(( forward_jacobian(TestFunc1<double>(3,3)) ));
|
CALL_SUBTEST(( forward_jacobian(TestFunc1<double>(3,3)) ));
|
||||||
|
#if EIGEN_HAS_VARIADIC_TEMPLATES
|
||||||
|
CALL_SUBTEST(( forward_jacobian_cpp11(integratorFunctor<double>(10)) ));
|
||||||
|
#endif
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user