mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-12 03:39:01 +08:00
add atan2 support in AutoDiff and remove superfluous std:: specializations
This commit is contained in:
parent
063042bca3
commit
f1d98aad1b
@ -503,8 +503,6 @@ struct scalar_product_traits<AutoDiffScalar<DerType>,T>
|
|||||||
|
|
||||||
} // end namespace internal
|
} // end namespace internal
|
||||||
|
|
||||||
} // end namespace Eigen
|
|
||||||
|
|
||||||
#define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
|
#define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
|
||||||
template<typename DerType> \
|
template<typename DerType> \
|
||||||
inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > \
|
inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > \
|
||||||
@ -515,58 +513,6 @@ struct scalar_product_traits<AutoDiffScalar<DerType>,T>
|
|||||||
CODE; \
|
CODE; \
|
||||||
}
|
}
|
||||||
|
|
||||||
namespace std
|
|
||||||
{
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs,
|
|
||||||
return ReturnType(std::abs(x.value()), x.derivatives() * (sign(x.value())));)
|
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sqrt,
|
|
||||||
Scalar sqrtx = std::sqrt(x.value());
|
|
||||||
return ReturnType(sqrtx,x.derivatives() * (Scalar(0.5) / sqrtx));)
|
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(cos,
|
|
||||||
return ReturnType(std::cos(x.value()), x.derivatives() * (-std::sin(x.value())));)
|
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sin,
|
|
||||||
return ReturnType(std::sin(x.value()),x.derivatives() * std::cos(x.value()));)
|
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(exp,
|
|
||||||
Scalar expx = std::exp(x.value());
|
|
||||||
return ReturnType(expx,x.derivatives() * expx);)
|
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
|
|
||||||
return ReturnType(std::log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
|
|
||||||
|
|
||||||
template<typename DerType>
|
|
||||||
inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<DerType>::Scalar>, const DerType> >
|
|
||||||
pow(const Eigen::AutoDiffScalar<DerType>& x, typename Eigen::internal::traits<DerType>::Scalar y)
|
|
||||||
{
|
|
||||||
using namespace Eigen;
|
|
||||||
typedef typename Eigen::internal::traits<DerType>::Scalar Scalar;
|
|
||||||
return AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const DerType> >(
|
|
||||||
std::pow(x.value(),y),
|
|
||||||
x.derivatives() * (y * std::pow(x.value(),y-1)));
|
|
||||||
}
|
|
||||||
|
|
||||||
}
|
|
||||||
|
|
||||||
#undef EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY
|
|
||||||
#define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
|
|
||||||
template<typename DerType> \
|
|
||||||
struct FUNC##_impl<Eigen::AutoDiffScalar<DerType> > \
|
|
||||||
{ \
|
|
||||||
static inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > \
|
|
||||||
run(const Eigen::AutoDiffScalar<DerType>& x) { \
|
|
||||||
using namespace Eigen; \
|
|
||||||
typedef typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar Scalar; \
|
|
||||||
typedef AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > ReturnType; \
|
|
||||||
CODE; \
|
|
||||||
} };
|
|
||||||
|
|
||||||
namespace Eigen {
|
|
||||||
|
|
||||||
namespace internal {
|
|
||||||
|
|
||||||
template<typename DerType>
|
template<typename DerType>
|
||||||
inline const AutoDiffScalar<DerType>& conj(const AutoDiffScalar<DerType>& x) { return x; }
|
inline const AutoDiffScalar<DerType>& conj(const AutoDiffScalar<DerType>& x) { return x; }
|
||||||
template<typename DerType>
|
template<typename DerType>
|
||||||
@ -575,34 +521,70 @@ template<typename DerType>
|
|||||||
inline typename DerType::Scalar imag(const AutoDiffScalar<DerType>&) { return 0.; }
|
inline typename DerType::Scalar imag(const AutoDiffScalar<DerType>&) { return 0.; }
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs,
|
||||||
|
using std::abs;
|
||||||
return ReturnType(abs(x.value()), x.derivatives() * (sign(x.value())));)
|
return ReturnType(abs(x.value()), x.derivatives() * (sign(x.value())));)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs2,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs2,
|
||||||
|
using internal::abs2;
|
||||||
return ReturnType(abs2(x.value()), x.derivatives() * (Scalar(2)*x.value()));)
|
return ReturnType(abs2(x.value()), x.derivatives() * (Scalar(2)*x.value()));)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sqrt,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sqrt,
|
||||||
|
using std::sqrt;
|
||||||
Scalar sqrtx = sqrt(x.value());
|
Scalar sqrtx = sqrt(x.value());
|
||||||
return ReturnType(sqrtx,x.derivatives() * (Scalar(0.5) / sqrtx));)
|
return ReturnType(sqrtx,x.derivatives() * (Scalar(0.5) / sqrtx));)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(cos,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(cos,
|
||||||
|
using std::cos;
|
||||||
|
using std::sin;
|
||||||
return ReturnType(cos(x.value()), x.derivatives() * (-sin(x.value())));)
|
return ReturnType(cos(x.value()), x.derivatives() * (-sin(x.value())));)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sin,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sin,
|
||||||
|
using std::sin;
|
||||||
|
using std::cos;
|
||||||
return ReturnType(sin(x.value()),x.derivatives() * cos(x.value()));)
|
return ReturnType(sin(x.value()),x.derivatives() * cos(x.value()));)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(exp,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(exp,
|
||||||
|
using std::exp;
|
||||||
Scalar expx = exp(x.value());
|
Scalar expx = exp(x.value());
|
||||||
return ReturnType(expx,x.derivatives() * expx);)
|
return ReturnType(expx,x.derivatives() * expx);)
|
||||||
|
|
||||||
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
|
EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
|
||||||
|
using std::log;
|
||||||
return ReturnType(log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
|
return ReturnType(log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
|
||||||
|
|
||||||
template<typename DerType>
|
template<typename DerType>
|
||||||
inline const AutoDiffScalar<CwiseUnaryOp<scalar_multiple_op<typename traits<DerType>::Scalar>, DerType> >
|
inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<DerType>::Scalar>, const DerType> >
|
||||||
pow(const AutoDiffScalar<DerType>& x, typename traits<DerType>::Scalar y)
|
pow(const Eigen::AutoDiffScalar<DerType>& x, typename Eigen::internal::traits<DerType>::Scalar y)
|
||||||
{ return std::pow(x,y);}
|
{
|
||||||
|
using namespace Eigen;
|
||||||
|
typedef typename Eigen::internal::traits<DerType>::Scalar Scalar;
|
||||||
|
return AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const DerType> >(
|
||||||
|
std::pow(x.value(),y),
|
||||||
|
x.derivatives() * (y * std::pow(x.value(),y-1)));
|
||||||
|
}
|
||||||
|
|
||||||
} // end namespace internal
|
|
||||||
|
template<typename DerTypeA,typename DerTypeB>
|
||||||
|
inline const AutoDiffScalar<Matrix<typename internal::traits<DerTypeA>::Scalar,Dynamic,1> >
|
||||||
|
atan2(const AutoDiffScalar<DerTypeA>& a, const AutoDiffScalar<DerTypeB>& b)
|
||||||
|
{
|
||||||
|
using std::atan2;
|
||||||
|
typedef typename internal::traits<DerTypeA>::Scalar Scalar;
|
||||||
|
typedef AutoDiffScalar<Matrix<Scalar,Dynamic,1> > PlainADS;
|
||||||
|
PlainADS ret;
|
||||||
|
ret.value() = atan2(a.value(), b.value());
|
||||||
|
|
||||||
|
Scalar tmp2 = a.value() * a.value();
|
||||||
|
Scalar tmp3 = b.value() * b.value();
|
||||||
|
Scalar tmp4 = tmp3/(tmp2+tmp3);
|
||||||
|
|
||||||
|
if (tmp4!=0)
|
||||||
|
ret.derivatives() = (a.derivatives() * b.value() - a.value() * b.derivatives()) * (tmp2+tmp3);
|
||||||
|
else
|
||||||
|
ret.derivatives().setZero();
|
||||||
|
|
||||||
|
return ret;
|
||||||
|
}
|
||||||
|
|
||||||
#undef EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY
|
#undef EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY
|
||||||
|
|
||||||
|
@ -28,13 +28,21 @@
|
|||||||
template<typename Scalar>
|
template<typename Scalar>
|
||||||
EIGEN_DONT_INLINE Scalar foo(const Scalar& x, const Scalar& y)
|
EIGEN_DONT_INLINE Scalar foo(const Scalar& x, const Scalar& y)
|
||||||
{
|
{
|
||||||
|
using namespace std;
|
||||||
// return x+std::sin(y);
|
// return x+std::sin(y);
|
||||||
EIGEN_ASM_COMMENT("mybegin");
|
EIGEN_ASM_COMMENT("mybegin");
|
||||||
return static_cast<Scalar>(x*2 - std::pow(x,2) + 2*std::sqrt(y*y) - 4 * std::sin(x) + 2 * std::cos(y) - std::exp(-0.5*x*x));
|
return static_cast<Scalar>(x*2 - pow(x,2) + 2*sqrt(y*y) - 4 * sin(x) + 2 * cos(y) - exp(-0.5*x*x));
|
||||||
//return x+2*y*x;//x*2 -std::pow(x,2);//(2*y/x);// - y*2;
|
//return x+2*y*x;//x*2 -std::pow(x,2);//(2*y/x);// - y*2;
|
||||||
EIGEN_ASM_COMMENT("myend");
|
EIGEN_ASM_COMMENT("myend");
|
||||||
}
|
}
|
||||||
|
|
||||||
|
template<typename Vector>
|
||||||
|
EIGEN_DONT_INLINE typename Vector::Scalar foo(const Vector& p)
|
||||||
|
{
|
||||||
|
typedef typename Vector::Scalar Scalar;
|
||||||
|
return (p-Vector(Scalar(-1),Scalar(1.))).norm();
|
||||||
|
}
|
||||||
|
|
||||||
template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
|
template<typename _Scalar, int NX=Dynamic, int NY=Dynamic>
|
||||||
struct TestFunc1
|
struct TestFunc1
|
||||||
{
|
{
|
||||||
@ -140,9 +148,23 @@ void test_autodiff_scalar()
|
|||||||
typedef AutoDiffScalar<Vector2f> AD;
|
typedef AutoDiffScalar<Vector2f> AD;
|
||||||
AD ax(1,Vector2f::UnitX());
|
AD ax(1,Vector2f::UnitX());
|
||||||
AD ay(2,Vector2f::UnitY());
|
AD ay(2,Vector2f::UnitY());
|
||||||
foo<AD>(ax,ay);
|
AD res = foo<AD>(ax,ay);
|
||||||
std::cerr << foo<AD>(ax,ay).value() << " <> "
|
std::cerr << res.value() << " <> "
|
||||||
<< foo<AD>(ax,ay).derivatives().transpose() << "\n\n";
|
<< res.derivatives().transpose() << "\n\n";
|
||||||
|
}
|
||||||
|
|
||||||
|
void test_autodiff_vector()
|
||||||
|
{
|
||||||
|
std::cerr << foo<Vector2f>(Vector2f(1,2)) << "\n";
|
||||||
|
typedef AutoDiffScalar<Vector2f> AD;
|
||||||
|
typedef Matrix<AD,2,1> VectorAD;
|
||||||
|
VectorAD p(AD(1),AD(-1));
|
||||||
|
p.x().derivatives() = Vector2f::UnitX();
|
||||||
|
p.y().derivatives() = Vector2f::UnitY();
|
||||||
|
|
||||||
|
AD res = foo<VectorAD>(p);
|
||||||
|
std::cerr << res.value() << " <> "
|
||||||
|
<< res.derivatives().transpose() << "\n\n";
|
||||||
}
|
}
|
||||||
|
|
||||||
void test_autodiff_jacobian()
|
void test_autodiff_jacobian()
|
||||||
@ -159,6 +181,7 @@ void test_autodiff_jacobian()
|
|||||||
void test_autodiff()
|
void test_autodiff()
|
||||||
{
|
{
|
||||||
test_autodiff_scalar();
|
test_autodiff_scalar();
|
||||||
test_autodiff_jacobian();
|
test_autodiff_vector();
|
||||||
|
// test_autodiff_jacobian();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
Loading…
x
Reference in New Issue
Block a user