10 #ifndef EIGEN_AUTODIFF_SCALAR_H
11 #define EIGEN_AUTODIFF_SCALAR_H
17 template<
typename A,
typename B>
18 struct make_coherent_impl {
19 static void run(A&, B&) {}
23 template<
typename A,
typename B>
24 void make_coherent(
const A& a,
const B&b)
26 make_coherent_impl<A,B>::run(a.const_cast_derived(), b.const_cast_derived());
29 template<
typename _DerType,
bool Enable>
struct auto_diff_special_op;
59 template<
typename _DerType>
61 :
public internal::auto_diff_special_op
62 <_DerType, !internal::is_same<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar,
63 typename NumTraits<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar>::Real>::value>
66 typedef internal::auto_diff_special_op
67 <_DerType, !internal::is_same<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar,
68 typename NumTraits<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar>::Real>::value> Base;
69 typedef typename internal::remove_all<_DerType>::type DerType;
70 typedef typename internal::traits<DerType>::Scalar Scalar;
71 typedef typename NumTraits<Scalar>::Real Real;
73 using Base::operator+;
74 using Base::operator*;
82 : m_value(value), m_derivatives(DerType::Zero(nbDer))
84 m_derivatives.coeffRef(derNumber) = Scalar(1);
92 if(m_derivatives.size()>0)
93 m_derivatives.setZero();
98 : m_value(value), m_derivatives(der)
101 template<
typename OtherDerType>
103 : m_value(other.value()), m_derivatives(other.derivatives())
106 friend std::ostream & operator << (std::ostream & s,
const AutoDiffScalar& a)
108 return s << a.value();
112 : m_value(other.value()), m_derivatives(other.derivatives())
115 template<
typename OtherDerType>
116 inline AutoDiffScalar& operator=(
const AutoDiffScalar<OtherDerType>& other)
118 m_value = other.value();
119 m_derivatives = other.derivatives();
125 m_value = other.value();
126 m_derivatives = other.derivatives();
133 inline const Scalar& value()
const {
return m_value; }
134 inline Scalar& value() {
return m_value; }
136 inline const DerType& derivatives()
const {
return m_derivatives; }
137 inline DerType& derivatives() {
return m_derivatives; }
139 inline bool operator< (
const Scalar& other)
const {
return m_value < other; }
140 inline bool operator<=(
const Scalar& other)
const {
return m_value <= other; }
141 inline bool operator> (
const Scalar& other)
const {
return m_value > other; }
142 inline bool operator>=(
const Scalar& other)
const {
return m_value >= other; }
143 inline bool operator==(
const Scalar& other)
const {
return m_value == other; }
144 inline bool operator!=(
const Scalar& other)
const {
return m_value != other; }
146 friend inline bool operator< (
const Scalar& a,
const AutoDiffScalar& b) {
return a < b.value(); }
147 friend inline bool operator<=(
const Scalar& a,
const AutoDiffScalar& b) {
return a <= b.value(); }
148 friend inline bool operator> (
const Scalar& a,
const AutoDiffScalar& b) {
return a > b.value(); }
149 friend inline bool operator>=(
const Scalar& a,
const AutoDiffScalar& b) {
return a >= b.value(); }
150 friend inline bool operator==(
const Scalar& a,
const AutoDiffScalar& b) {
return a == b.value(); }
151 friend inline bool operator!=(
const Scalar& a,
const AutoDiffScalar& b) {
return a != b.value(); }
153 template<
typename OtherDerType>
inline bool operator< (const AutoDiffScalar<OtherDerType>& b)
const {
return m_value < b.value(); }
154 template<
typename OtherDerType>
inline bool operator<=(const AutoDiffScalar<OtherDerType>& b)
const {
return m_value <= b.value(); }
155 template<
typename OtherDerType>
inline bool operator> (
const AutoDiffScalar<OtherDerType>& b)
const {
return m_value > b.value(); }
156 template<
typename OtherDerType>
inline bool operator>=(
const AutoDiffScalar<OtherDerType>& b)
const {
return m_value >= b.value(); }
157 template<
typename OtherDerType>
inline bool operator==(
const AutoDiffScalar<OtherDerType>& b)
const {
return m_value == b.value(); }
158 template<
typename OtherDerType>
inline bool operator!=(
const AutoDiffScalar<OtherDerType>& b)
const {
return m_value != b.value(); }
160 inline const AutoDiffScalar<DerType&> operator+(
const Scalar& other)
const
162 return AutoDiffScalar<DerType&>(m_value + other, m_derivatives);
165 friend inline const AutoDiffScalar<DerType&> operator+(
const Scalar& a,
const AutoDiffScalar& b)
167 return AutoDiffScalar<DerType&>(a + b.value(), b.derivatives());
186 template<
typename OtherDerType>
187 inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
const DerType,
const typename internal::remove_all<OtherDerType>::type> >
188 operator+(
const AutoDiffScalar<OtherDerType>& other)
const
190 internal::make_coherent(m_derivatives, other.derivatives());
191 return AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
const DerType,
const typename internal::remove_all<OtherDerType>::type> >(
192 m_value + other.value(),
193 m_derivatives + other.derivatives());
196 template<
typename OtherDerType>
198 operator+=(
const AutoDiffScalar<OtherDerType>& other)
200 (*this) = (*this) + other;
204 inline const AutoDiffScalar<DerType&> operator-(
const Scalar& b)
const
206 return AutoDiffScalar<DerType&>(m_value - b, m_derivatives);
209 friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>,
const DerType> >
212 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>,
const DerType> >
213 (a - b.value(), -b.derivatives());
222 template<
typename OtherDerType>
223 inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
const DerType,
const typename internal::remove_all<OtherDerType>::type> >
224 operator-(
const AutoDiffScalar<OtherDerType>& other)
const
226 internal::make_coherent(m_derivatives, other.derivatives());
227 return AutoDiffScalar<CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
const DerType,
const typename internal::remove_all<OtherDerType>::type> >(
228 m_value - other.value(),
229 m_derivatives - other.derivatives());
232 template<
typename OtherDerType>
234 operator-=(
const AutoDiffScalar<OtherDerType>& other)
236 *
this = *
this - other;
240 inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>,
const DerType> >
243 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>,
const DerType> >(
248 inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >
249 operator*(
const Scalar& other)
const
251 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >(
253 (m_derivatives * other));
256 friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >
259 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >(
261 a.derivatives() * other);
280 inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >
281 operator/(
const Scalar& other)
const
283 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >(
285 (m_derivatives * (Scalar(1)/other)));
288 friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >
291 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType> >(
293 a.derivatives() * (Scalar(-other) / (a.value()*a.value())));
312 template<
typename OtherDerType>
313 inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
314 const CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
315 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType>,
316 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const typename internal::remove_all<OtherDerType>::type > > > >
317 operator/(
const AutoDiffScalar<OtherDerType>& other)
const
319 internal::make_coherent(m_derivatives, other.derivatives());
320 return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
321 const CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
322 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType>,
323 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const typename internal::remove_all<OtherDerType>::type > > > >(
324 m_value / other.value(),
325 ((m_derivatives * other.value()) - (m_value * other.derivatives()))
326 * (Scalar(1)/(other.value()*other.value())));
329 template<
typename OtherDerType>
330 inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
331 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType>,
332 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const typename internal::remove_all<OtherDerType>::type> > >
333 operator*(
const AutoDiffScalar<OtherDerType>& other)
const
335 internal::make_coherent(m_derivatives, other.derivatives());
336 return AutoDiffScalar<const CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
337 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const DerType>,
338 const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
const typename internal::remove_all<OtherDerType>::type > > >(
339 m_value * other.value(),
340 (m_derivatives * other.value()) + (m_value * other.derivatives()));
345 *
this = *
this * other;
349 template<
typename OtherDerType>
350 inline AutoDiffScalar& operator*=(
const AutoDiffScalar<OtherDerType>& other)
352 *
this = *
this * other;
358 *
this = *
this / other;
362 template<
typename OtherDerType>
363 inline AutoDiffScalar& operator/=(
const AutoDiffScalar<OtherDerType>& other)
365 *
this = *
this / other;
371 DerType m_derivatives;
377 template<
typename _DerType>
378 struct auto_diff_special_op<_DerType, true>
382 typedef typename remove_all<_DerType>::type DerType;
383 typedef typename traits<DerType>::Scalar Scalar;
384 typedef typename NumTraits<Scalar>::Real Real;
396 const AutoDiffScalar<_DerType>& derived()
const {
return *
static_cast<const AutoDiffScalar<_DerType>*
>(
this); }
397 AutoDiffScalar<_DerType>& derived() {
return *
static_cast<AutoDiffScalar<_DerType>*
>(
this); }
400 inline const AutoDiffScalar<DerType&> operator+(
const Real& other)
const
402 return AutoDiffScalar<DerType&>(derived().value() + other, derived().derivatives());
405 friend inline const AutoDiffScalar<DerType&> operator+(
const Real& a,
const AutoDiffScalar<_DerType>& b)
407 return AutoDiffScalar<DerType&>(a + b.value(), b.derivatives());
410 inline AutoDiffScalar<_DerType>& operator+=(
const Real& other)
412 derived().value() += other;
417 inline const AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >
418 operator*(
const Real& other)
const
420 return AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >(
421 derived().value() * other,
422 derived().derivatives() * other);
425 friend inline const AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >
426 operator*(
const Real& other,
const AutoDiffScalar<_DerType>& a)
428 return AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >(
430 a.derivatives() * other);
433 inline AutoDiffScalar<_DerType>& operator*=(
const Scalar& other)
435 *
this = *
this * other;
440 template<
typename _DerType>
441 struct auto_diff_special_op<_DerType, false>
443 void operator*()
const;
444 void operator-()
const;
445 void operator+()
const;
448 template<
typename A_Scalar,
int A_Rows,
int A_Cols,
int A_Options,
int A_MaxRows,
int A_MaxCols,
typename B>
449 struct make_coherent_impl<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>, B> {
450 typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> A;
451 static void run(A& a, B& b) {
452 if((A_Rows==Dynamic || A_Cols==Dynamic) && (a.size()==0))
460 template<
typename A,
typename B_Scalar,
int B_Rows,
int B_Cols,
int B_Options,
int B_MaxRows,
int B_MaxCols>
461 struct make_coherent_impl<A, Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> > {
462 typedef Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> B;
463 static void run(A& a, B& b) {
464 if((B_Rows==Dynamic || B_Cols==Dynamic) && (b.size()==0))
472 template<
typename A_Scalar,
int A_Rows,
int A_Cols,
int A_Options,
int A_MaxRows,
int A_MaxCols,
473 typename B_Scalar,
int B_Rows,
int B_Cols,
int B_Options,
int B_MaxRows,
int B_MaxCols>
474 struct make_coherent_impl<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>,
475 Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> > {
476 typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> A;
477 typedef Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> B;
478 static void run(A& a, B& b) {
479 if((A_Rows==Dynamic || A_Cols==Dynamic) && (a.size()==0))
484 else if((B_Rows==Dynamic || B_Cols==Dynamic) && (b.size()==0))
492 template<
typename A_Scalar,
int A_Rows,
int A_Cols,
int A_Options,
int A_MaxRows,
int A_MaxCols>
493 struct scalar_product_traits<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>,A_Scalar>
495 enum { Defined = 1 };
496 typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> ReturnType;
499 template<
typename A_Scalar,
int A_Rows,
int A_Cols,
int A_Options,
int A_MaxRows,
int A_MaxCols>
500 struct scalar_product_traits<A_Scalar, Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> >
502 enum { Defined = 1 };
503 typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> ReturnType;
506 template<
typename DerType>
507 struct scalar_product_traits<AutoDiffScalar<DerType>,typename DerType::Scalar>
509 enum { Defined = 1 };
510 typedef AutoDiffScalar<DerType> ReturnType;
513 template<
typename DerType>
514 struct scalar_product_traits<typename DerType::Scalar,AutoDiffScalar<DerType> >
516 enum { Defined = 1 };
517 typedef AutoDiffScalar<DerType> ReturnType;
522 #define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
523 template<typename DerType> \
524 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> > \
525 FUNC(const Eigen::AutoDiffScalar<DerType>& x) { \
526 using namespace Eigen; \
527 typedef typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar Scalar; \
528 typedef AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > ReturnType; \
532 template<
typename DerType>
533 inline const AutoDiffScalar<DerType>& conj(
const AutoDiffScalar<DerType>& x) {
return x; }
534 template<
typename DerType>
535 inline const AutoDiffScalar<DerType>& real(
const AutoDiffScalar<DerType>& x) {
return x; }
536 template<
typename DerType>
537 inline typename DerType::Scalar imag(
const AutoDiffScalar<DerType>&) {
return 0.; }
538 template<
typename DerType,
typename T>
539 inline AutoDiffScalar<DerType> (min)(
const AutoDiffScalar<DerType>& x,
const T& y) {
return (x <= y ? x : y); }
540 template<
typename DerType,
typename T>
541 inline AutoDiffScalar<DerType> (max)(
const AutoDiffScalar<DerType>& x,
const T& y) {
return (x >= y ? x : y); }
542 template<
typename DerType,
typename T>
543 inline AutoDiffScalar<DerType> (min)(
const T& x,
const AutoDiffScalar<DerType>& y) {
return (x < y ? x : y); }
544 template<
typename DerType,
typename T>
545 inline AutoDiffScalar<DerType> (max)(
const T& x,
const AutoDiffScalar<DerType>& y) {
return (x > y ? x : y); }
547 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs,
549 return ReturnType(abs(x.value()), x.derivatives() * (x.value()<0 ? -1 : 1) );)
551 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs2,
553 return ReturnType(abs2(x.value()), x.derivatives() * (Scalar(2)*x.value()));)
555 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sqrt,
557 Scalar sqrtx = sqrt(x.value());
558 return ReturnType(sqrtx,x.derivatives() * (Scalar(0.5) / sqrtx));)
560 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(cos,
563 return ReturnType(cos(x.value()), x.derivatives() * (-sin(x.value())));)
565 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sin,
568 return ReturnType(sin(x.value()),x.derivatives() * cos(x.value()));)
570 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(exp,
572 Scalar expx = exp(x.value());
573 return ReturnType(expx,x.derivatives() * expx);)
575 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
577 return ReturnType(log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
579 template<typename DerType>
580 inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<
typename Eigen::internal::traits<DerType>::Scalar>,
const DerType> >
583 using namespace Eigen;
584 typedef typename Eigen::internal::traits<DerType>::Scalar Scalar;
586 std::pow(x.value(),y),
587 x.derivatives() * (y * std::pow(x.value(),y-1)));
591 template<
typename DerTypeA,
typename DerTypeB>
597 typedef typename internal::traits<DerTypeA>::Scalar Scalar;
600 ret.value() = atan2(a.value(), b.value());
602 Scalar tmp2 = a.value() * a.value();
603 Scalar tmp3 = b.value() * b.value();
604 Scalar tmp4 = tmp3/(tmp2+tmp3);
607 ret.derivatives() = (a.derivatives() * b.value() - a.value() * b.derivatives()) * (tmp2+tmp3);
612 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(tan,
615 return ReturnType(tan(x.value()),x.derivatives() * (Scalar(1)/numext::abs2(cos(x.value()))));)
617 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(asin,
620 return ReturnType(asin(x.value()),x.derivatives() * (Scalar(1)/sqrt(1-numext::abs2(x.value()))));)
622 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(acos,
625 return ReturnType(acos(x.value()),x.derivatives() * (Scalar(-1)/sqrt(1-numext::abs2(x.value()))));)
627 #undef EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY
630 : NumTraits<
typename NumTraits<typename DerType::Scalar>::Real >
636 RequireInitialization = 1
642 #endif // EIGEN_AUTODIFF_SCALAR_H
A scalar type replacement with automatic differentation capability.
Definition: AutoDiffScalar.h:60
AutoDiffScalar(const Scalar &value, const DerType &der)
Definition: AutoDiffScalar.h:97
AutoDiffScalar()
Definition: AutoDiffScalar.h:77
AutoDiffScalar(const Scalar &value, int nbDer, int derNumber)
Definition: AutoDiffScalar.h:81
AutoDiffScalar(const Real &value)
Definition: AutoDiffScalar.h:89