AutoDiffScalar.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_AUTODIFF_SCALAR_H
00011 #define EIGEN_AUTODIFF_SCALAR_H
00012 
00013 namespace Eigen {
00014 
00015 namespace internal {
00016 
00017 template<typename A, typename B>
00018 struct make_coherent_impl {
00019   static void run(A&, B&) {}
00020 };
00021 
00022 // resize a to match b is a.size()==0, and conversely.
00023 template<typename A, typename B>
00024 void make_coherent(const A& a, const B&b)
00025 {
00026   make_coherent_impl<A,B>::run(a.const_cast_derived(), b.const_cast_derived());
00027 }
00028 
00029 template<typename _DerType, bool Enable> struct auto_diff_special_op;
00030 
00031 } // end namespace internal
00032 
00059 template<typename _DerType>
00060 class AutoDiffScalar
00061   : public internal::auto_diff_special_op
00062             <_DerType, !internal::is_same<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar,
00063                                         typename NumTraits<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar>::Real>::value>
00064 {
00065   public:
00066     typedef internal::auto_diff_special_op
00067             <_DerType, !internal::is_same<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar,
00068                        typename NumTraits<typename internal::traits<typename internal::remove_all<_DerType>::type>::Scalar>::Real>::value> Base;
00069     typedef typename internal::remove_all<_DerType>::type DerType;
00070     typedef typename internal::traits<DerType>::Scalar Scalar;
00071     typedef typename NumTraits<Scalar>::Real Real;
00072 
00073     using Base::operator+;
00074     using Base::operator*;
00075 
00077     AutoDiffScalar() {}
00078 
00081     AutoDiffScalar(const Scalar& value, int nbDer, int derNumber)
00082       : m_value(value), m_derivatives(DerType::Zero(nbDer))
00083     {
00084       m_derivatives.coeffRef(derNumber) = Scalar(1);
00085     }
00086 
00089     /*explicit*/ AutoDiffScalar(const Real& value)
00090       : m_value(value)
00091     {
00092       if(m_derivatives.size()>0)
00093         m_derivatives.setZero();
00094     }
00095 
00097     AutoDiffScalar(const Scalar& value, const DerType& der)
00098       : m_value(value), m_derivatives(der)
00099     {}
00100 
00101     template<typename OtherDerType>
00102     AutoDiffScalar(const AutoDiffScalar<OtherDerType>& other)
00103       : m_value(other.value()), m_derivatives(other.derivatives())
00104     {}
00105 
00106     friend  std::ostream & operator << (std::ostream & s, const AutoDiffScalar& a)
00107     {
00108       return s << a.value();
00109     }
00110 
00111     AutoDiffScalar(const AutoDiffScalar& other)
00112       : m_value(other.value()), m_derivatives(other.derivatives())
00113     {}
00114 
00115     template<typename OtherDerType>
00116     inline AutoDiffScalar& operator=(const AutoDiffScalar<OtherDerType>& other)
00117     {
00118       m_value = other.value();
00119       m_derivatives = other.derivatives();
00120       return *this;
00121     }
00122 
00123     inline AutoDiffScalar& operator=(const AutoDiffScalar& other)
00124     {
00125       m_value = other.value();
00126       m_derivatives = other.derivatives();
00127       return *this;
00128     }
00129 
00130 //     inline operator const Scalar& () const { return m_value; }
00131 //     inline operator Scalar& () { return m_value; }
00132 
00133     inline const Scalar& value() const { return m_value; }
00134     inline Scalar& value() { return m_value; }
00135 
00136     inline const DerType& derivatives() const { return m_derivatives; }
00137     inline DerType& derivatives() { return m_derivatives; }
00138 
00139     inline bool operator< (const Scalar& other) const  { return m_value <  other; }
00140     inline bool operator<=(const Scalar& other) const  { return m_value <= other; }
00141     inline bool operator> (const Scalar& other) const  { return m_value >  other; }
00142     inline bool operator>=(const Scalar& other) const  { return m_value >= other; }
00143     inline bool operator==(const Scalar& other) const  { return m_value == other; }
00144     inline bool operator!=(const Scalar& other) const  { return m_value != other; }
00145 
00146     friend inline bool operator< (const Scalar& a, const AutoDiffScalar& b) { return a <  b.value(); }
00147     friend inline bool operator<=(const Scalar& a, const AutoDiffScalar& b) { return a <= b.value(); }
00148     friend inline bool operator> (const Scalar& a, const AutoDiffScalar& b) { return a >  b.value(); }
00149     friend inline bool operator>=(const Scalar& a, const AutoDiffScalar& b) { return a >= b.value(); }
00150     friend inline bool operator==(const Scalar& a, const AutoDiffScalar& b) { return a == b.value(); }
00151     friend inline bool operator!=(const Scalar& a, const AutoDiffScalar& b) { return a != b.value(); }
00152 
00153     template<typename OtherDerType> inline bool operator< (const AutoDiffScalar<OtherDerType>& b) const  { return m_value <  b.value(); }
00154     template<typename OtherDerType> inline bool operator<=(const AutoDiffScalar<OtherDerType>& b) const  { return m_value <= b.value(); }
00155     template<typename OtherDerType> inline bool operator> (const AutoDiffScalar<OtherDerType>& b) const  { return m_value >  b.value(); }
00156     template<typename OtherDerType> inline bool operator>=(const AutoDiffScalar<OtherDerType>& b) const  { return m_value >= b.value(); }
00157     template<typename OtherDerType> inline bool operator==(const AutoDiffScalar<OtherDerType>& b) const  { return m_value == b.value(); }
00158     template<typename OtherDerType> inline bool operator!=(const AutoDiffScalar<OtherDerType>& b) const  { return m_value != b.value(); }
00159 
00160     inline const AutoDiffScalar<DerType&> operator+(const Scalar& other) const
00161     {
00162       return AutoDiffScalar<DerType&>(m_value + other, m_derivatives);
00163     }
00164 
00165     friend inline const AutoDiffScalar<DerType&> operator+(const Scalar& a, const AutoDiffScalar& b)
00166     {
00167       return AutoDiffScalar<DerType&>(a + b.value(), b.derivatives());
00168     }
00169 
00170 //     inline const AutoDiffScalar<DerType&> operator+(const Real& other) const
00171 //     {
00172 //       return AutoDiffScalar<DerType&>(m_value + other, m_derivatives);
00173 //     }
00174 
00175 //     friend inline const AutoDiffScalar<DerType&> operator+(const Real& a, const AutoDiffScalar& b)
00176 //     {
00177 //       return AutoDiffScalar<DerType&>(a + b.value(), b.derivatives());
00178 //     }
00179 
00180     inline AutoDiffScalar& operator+=(const Scalar& other)
00181     {
00182       value() += other;
00183       return *this;
00184     }
00185 
00186     template<typename OtherDerType>
00187     inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,const DerType,const typename internal::remove_all<OtherDerType>::type> >
00188     operator+(const AutoDiffScalar<OtherDerType>& other) const
00189     {
00190       internal::make_coherent(m_derivatives, other.derivatives());
00191       return AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,const DerType,const typename internal::remove_all<OtherDerType>::type> >(
00192         m_value + other.value(),
00193         m_derivatives + other.derivatives());
00194     }
00195 
00196     template<typename OtherDerType>
00197     inline AutoDiffScalar&
00198     operator+=(const AutoDiffScalar<OtherDerType>& other)
00199     {
00200       (*this) = (*this) + other;
00201       return *this;
00202     }
00203 
00204     inline const AutoDiffScalar<DerType&> operator-(const Scalar& b) const
00205     {
00206       return AutoDiffScalar<DerType&>(m_value - b, m_derivatives);
00207     }
00208 
00209     friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>, const DerType> >
00210     operator-(const Scalar& a, const AutoDiffScalar& b)
00211     {
00212       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>, const DerType> >
00213             (a - b.value(), -b.derivatives());
00214     }
00215 
00216     inline AutoDiffScalar& operator-=(const Scalar& other)
00217     {
00218       value() -= other;
00219       return *this;
00220     }
00221 
00222     template<typename OtherDerType>
00223     inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_difference_op<Scalar>, const DerType,const typename internal::remove_all<OtherDerType>::type> >
00224     operator-(const AutoDiffScalar<OtherDerType>& other) const
00225     {
00226       internal::make_coherent(m_derivatives, other.derivatives());
00227       return AutoDiffScalar<CwiseBinaryOp<internal::scalar_difference_op<Scalar>, const DerType,const typename internal::remove_all<OtherDerType>::type> >(
00228         m_value - other.value(),
00229         m_derivatives - other.derivatives());
00230     }
00231 
00232     template<typename OtherDerType>
00233     inline AutoDiffScalar&
00234     operator-=(const AutoDiffScalar<OtherDerType>& other)
00235     {
00236       *this = *this - other;
00237       return *this;
00238     }
00239 
00240     inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>, const DerType> >
00241     operator-() const
00242     {
00243       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_opposite_op<Scalar>, const DerType> >(
00244         -m_value,
00245         -m_derivatives);
00246     }
00247 
00248     inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >
00249     operator*(const Scalar& other) const
00250     {
00251       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >(
00252         m_value * other,
00253         (m_derivatives * other));
00254     }
00255 
00256     friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >
00257     operator*(const Scalar& other, const AutoDiffScalar& a)
00258     {
00259       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >(
00260         a.value() * other,
00261         a.derivatives() * other);
00262     }
00263 
00264 //     inline const AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >
00265 //     operator*(const Real& other) const
00266 //     {
00267 //       return AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >(
00268 //         m_value * other,
00269 //         (m_derivatives * other));
00270 //     }
00271 //
00272 //     friend inline const AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >
00273 //     operator*(const Real& other, const AutoDiffScalar& a)
00274 //     {
00275 //       return AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >(
00276 //         a.value() * other,
00277 //         a.derivatives() * other);
00278 //     }
00279 
00280     inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >
00281     operator/(const Scalar& other) const
00282     {
00283       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >(
00284         m_value / other,
00285         (m_derivatives * (Scalar(1)/other)));
00286     }
00287 
00288     friend inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >
00289     operator/(const Scalar& other, const AutoDiffScalar& a)
00290     {
00291       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType> >(
00292         other / a.value(),
00293         a.derivatives() * (Scalar(-other) / (a.value()*a.value())));
00294     }
00295 
00296 //     inline const AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >
00297 //     operator/(const Real& other) const
00298 //     {
00299 //       return AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >(
00300 //         m_value / other,
00301 //         (m_derivatives * (Real(1)/other)));
00302 //     }
00303 //
00304 //     friend inline const AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >
00305 //     operator/(const Real& other, const AutoDiffScalar& a)
00306 //     {
00307 //       return AutoDiffScalar<typename CwiseUnaryOp<internal::scalar_multiple_op<Real>, DerType>::Type >(
00308 //         other / a.value(),
00309 //         a.derivatives() * (-Real(1)/other));
00310 //     }
00311 
00312     template<typename OtherDerType>
00313     inline const AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
00314         const CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
00315           const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType>,
00316           const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const typename internal::remove_all<OtherDerType>::type > > > >
00317     operator/(const AutoDiffScalar<OtherDerType>& other) const
00318     {
00319       internal::make_coherent(m_derivatives, other.derivatives());
00320       return AutoDiffScalar<CwiseUnaryOp<internal::scalar_multiple_op<Scalar>,
00321         const CwiseBinaryOp<internal::scalar_difference_op<Scalar>,
00322           const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType>,
00323           const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const typename internal::remove_all<OtherDerType>::type > > > >(
00324         m_value / other.value(),
00325           ((m_derivatives * other.value()) - (m_value * other.derivatives()))
00326         * (Scalar(1)/(other.value()*other.value())));
00327     }
00328 
00329     template<typename OtherDerType>
00330     inline const AutoDiffScalar<CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
00331         const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType>,
00332         const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const typename internal::remove_all<OtherDerType>::type> > >
00333     operator*(const AutoDiffScalar<OtherDerType>& other) const
00334     {
00335       internal::make_coherent(m_derivatives, other.derivatives());
00336       return AutoDiffScalar<const CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
00337         const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DerType>,
00338         const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const typename internal::remove_all<OtherDerType>::type > > >(
00339         m_value * other.value(),
00340         (m_derivatives * other.value()) + (m_value * other.derivatives()));
00341     }
00342 
00343     inline AutoDiffScalar& operator*=(const Scalar& other)
00344     {
00345       *this = *this * other;
00346       return *this;
00347     }
00348 
00349     template<typename OtherDerType>
00350     inline AutoDiffScalar& operator*=(const AutoDiffScalar<OtherDerType>& other)
00351     {
00352       *this = *this * other;
00353       return *this;
00354     }
00355 
00356     inline AutoDiffScalar& operator/=(const Scalar& other)
00357     {
00358       *this = *this / other;
00359       return *this;
00360     }
00361 
00362     template<typename OtherDerType>
00363     inline AutoDiffScalar& operator/=(const AutoDiffScalar<OtherDerType>& other)
00364     {
00365       *this = *this / other;
00366       return *this;
00367     }
00368 
00369   protected:
00370     Scalar m_value;
00371     DerType m_derivatives;
00372 
00373 };
00374 
00375 namespace internal {
00376 
00377 template<typename _DerType>
00378 struct auto_diff_special_op<_DerType, true>
00379 //   : auto_diff_scalar_op<_DerType, typename NumTraits<Scalar>::Real,
00380 //                            is_same<Scalar,typename NumTraits<Scalar>::Real>::value>
00381 {
00382   typedef typename remove_all<_DerType>::type DerType;
00383   typedef typename traits<DerType>::Scalar Scalar;
00384   typedef typename NumTraits<Scalar>::Real Real;
00385 
00386 //   typedef auto_diff_scalar_op<_DerType, typename NumTraits<Scalar>::Real,
00387 //                            is_same<Scalar,typename NumTraits<Scalar>::Real>::value> Base;
00388 
00389 //   using Base::operator+;
00390 //   using Base::operator+=;
00391 //   using Base::operator-;
00392 //   using Base::operator-=;
00393 //   using Base::operator*;
00394 //   using Base::operator*=;
00395 
00396   const AutoDiffScalar<_DerType>& derived() const { return *static_cast<const AutoDiffScalar<_DerType>*>(this); }
00397   AutoDiffScalar<_DerType>& derived() { return *static_cast<AutoDiffScalar<_DerType>*>(this); }
00398 
00399 
00400   inline const AutoDiffScalar<DerType&> operator+(const Real& other) const
00401   {
00402     return AutoDiffScalar<DerType&>(derived().value() + other, derived().derivatives());
00403   }
00404 
00405   friend inline const AutoDiffScalar<DerType&> operator+(const Real& a, const AutoDiffScalar<_DerType>& b)
00406   {
00407     return AutoDiffScalar<DerType&>(a + b.value(), b.derivatives());
00408   }
00409 
00410   inline AutoDiffScalar<_DerType>& operator+=(const Real& other)
00411   {
00412     derived().value() += other;
00413     return derived();
00414   }
00415 
00416 
00417   inline const AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >
00418   operator*(const Real& other) const
00419   {
00420     return AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >(
00421       derived().value() * other,
00422       derived().derivatives() * other);
00423   }
00424 
00425   friend inline const AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >
00426   operator*(const Real& other, const AutoDiffScalar<_DerType>& a)
00427   {
00428     return AutoDiffScalar<typename CwiseUnaryOp<scalar_multiple2_op<Scalar,Real>, DerType>::Type >(
00429       a.value() * other,
00430       a.derivatives() * other);
00431   }
00432 
00433   inline AutoDiffScalar<_DerType>& operator*=(const Scalar& other)
00434   {
00435     *this = *this * other;
00436     return derived();
00437   }
00438 };
00439 
00440 template<typename _DerType>
00441 struct auto_diff_special_op<_DerType, false>
00442 {
00443   void operator*() const;
00444   void operator-() const;
00445   void operator+() const;
00446 };
00447 
00448 template<typename A_Scalar, int A_Rows, int A_Cols, int A_Options, int A_MaxRows, int A_MaxCols, typename B>
00449 struct make_coherent_impl<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>, B> {
00450   typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> A;
00451   static void run(A& a, B& b) {
00452     if((A_Rows==Dynamic || A_Cols==Dynamic) && (a.size()==0))
00453     {
00454       a.resize(b.size());
00455       a.setZero();
00456     }
00457   }
00458 };
00459 
00460 template<typename A, typename B_Scalar, int B_Rows, int B_Cols, int B_Options, int B_MaxRows, int B_MaxCols>
00461 struct make_coherent_impl<A, Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> > {
00462   typedef Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> B;
00463   static void run(A& a, B& b) {
00464     if((B_Rows==Dynamic || B_Cols==Dynamic) && (b.size()==0))
00465     {
00466       b.resize(a.size());
00467       b.setZero();
00468     }
00469   }
00470 };
00471 
00472 template<typename A_Scalar, int A_Rows, int A_Cols, int A_Options, int A_MaxRows, int A_MaxCols,
00473          typename B_Scalar, int B_Rows, int B_Cols, int B_Options, int B_MaxRows, int B_MaxCols>
00474 struct make_coherent_impl<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>,
00475                              Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> > {
00476   typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> A;
00477   typedef Matrix<B_Scalar, B_Rows, B_Cols, B_Options, B_MaxRows, B_MaxCols> B;
00478   static void run(A& a, B& b) {
00479     if((A_Rows==Dynamic || A_Cols==Dynamic) && (a.size()==0))
00480     {
00481       a.resize(b.size());
00482       a.setZero();
00483     }
00484     else if((B_Rows==Dynamic || B_Cols==Dynamic) && (b.size()==0))
00485     {
00486       b.resize(a.size());
00487       b.setZero();
00488     }
00489   }
00490 };
00491 
00492 template<typename A_Scalar, int A_Rows, int A_Cols, int A_Options, int A_MaxRows, int A_MaxCols>
00493 struct scalar_product_traits<Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols>,A_Scalar>
00494 {
00495   enum { Defined = 1 };
00496   typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> ReturnType;
00497 };
00498 
00499 template<typename A_Scalar, int A_Rows, int A_Cols, int A_Options, int A_MaxRows, int A_MaxCols>
00500 struct scalar_product_traits<A_Scalar, Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> >
00501 {
00502   enum { Defined = 1 };
00503   typedef Matrix<A_Scalar, A_Rows, A_Cols, A_Options, A_MaxRows, A_MaxCols> ReturnType;
00504 };
00505 
00506 template<typename DerType>
00507 struct scalar_product_traits<AutoDiffScalar<DerType>,typename DerType::Scalar>
00508 {
00509   enum { Defined = 1 };
00510   typedef AutoDiffScalar<DerType> ReturnType;
00511 };
00512 
00513 template<typename DerType>
00514 struct scalar_product_traits<typename DerType::Scalar,AutoDiffScalar<DerType> >
00515 {
00516   enum { Defined = 1 };
00517   typedef AutoDiffScalar<DerType> ReturnType;
00518 };
00519 
00520 } // end namespace internal
00521 
00522 #define EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(FUNC,CODE) \
00523   template<typename DerType> \
00524   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> > \
00525   FUNC(const Eigen::AutoDiffScalar<DerType>& x) { \
00526     using namespace Eigen; \
00527     typedef typename Eigen::internal::traits<typename Eigen::internal::remove_all<DerType>::type>::Scalar Scalar; \
00528     typedef AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const typename Eigen::internal::remove_all<DerType>::type> > ReturnType; \
00529     CODE; \
00530   }
00531 
00532 template<typename DerType>
00533 inline const AutoDiffScalar<DerType>& conj(const AutoDiffScalar<DerType>& x)  { return x; }
00534 template<typename DerType>
00535 inline const AutoDiffScalar<DerType>& real(const AutoDiffScalar<DerType>& x)  { return x; }
00536 template<typename DerType>
00537 inline typename DerType::Scalar imag(const AutoDiffScalar<DerType>&)    { return 0.; }
00538 template<typename DerType, typename T>
00539 inline AutoDiffScalar<DerType> (min)(const AutoDiffScalar<DerType>& x, const T& y)    { return (x <= y ? x : y); }
00540 template<typename DerType, typename T>
00541 inline AutoDiffScalar<DerType> (max)(const AutoDiffScalar<DerType>& x, const T& y)    { return (x >= y ? x : y); }
00542 template<typename DerType, typename T>
00543 inline AutoDiffScalar<DerType> (min)(const T& x, const AutoDiffScalar<DerType>& y)    { return (x < y ? x : y); }
00544 template<typename DerType, typename T>
00545 inline AutoDiffScalar<DerType> (max)(const T& x, const AutoDiffScalar<DerType>& y)    { return (x > y ? x : y); }
00546 
00547 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs,
00548   using std::abs;
00549   return ReturnType(abs(x.value()), x.derivatives() * (x.value()<0 ? -1 : 1) );)
00550 
00551 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(abs2,
00552   using numext::abs2;
00553   return ReturnType(abs2(x.value()), x.derivatives() * (Scalar(2)*x.value()));)
00554 
00555 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sqrt,
00556   using std::sqrt;
00557   Scalar sqrtx = sqrt(x.value());
00558   return ReturnType(sqrtx,x.derivatives() * (Scalar(0.5) / sqrtx));)
00559 
00560 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(cos,
00561   using std::cos;
00562   using std::sin;
00563   return ReturnType(cos(x.value()), x.derivatives() * (-sin(x.value())));)
00564 
00565 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(sin,
00566   using std::sin;
00567   using std::cos;
00568   return ReturnType(sin(x.value()),x.derivatives() * cos(x.value()));)
00569 
00570 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(exp,
00571   using std::exp;
00572   Scalar expx = exp(x.value());
00573   return ReturnType(expx,x.derivatives() * expx);)
00574 
00575 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(log,
00576   using std::log;
00577   return ReturnType(log(x.value()),x.derivatives() * (Scalar(1)/x.value()));)
00578 
00579 template<typename DerType>
00580 inline const Eigen::AutoDiffScalar<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<typename Eigen::internal::traits<DerType>::Scalar>, const DerType> >
00581 pow(const Eigen::AutoDiffScalar<DerType>& x, typename Eigen::internal::traits<DerType>::Scalar y)
00582 {
00583   using namespace Eigen;
00584   typedef typename Eigen::internal::traits<DerType>::Scalar Scalar;
00585   return AutoDiffScalar<CwiseUnaryOp<Eigen::internal::scalar_multiple_op<Scalar>, const DerType> >(
00586     std::pow(x.value(),y),
00587     x.derivatives() * (y * std::pow(x.value(),y-1)));
00588 }
00589 
00590 
00591 template<typename DerTypeA,typename DerTypeB>
00592 inline const AutoDiffScalar<Matrix<typename internal::traits<DerTypeA>::Scalar,Dynamic,1> >
00593 atan2(const AutoDiffScalar<DerTypeA>& a, const AutoDiffScalar<DerTypeB>& b)
00594 {
00595   using std::atan2;
00596   using std::max;
00597   typedef typename internal::traits<DerTypeA>::Scalar Scalar;
00598   typedef AutoDiffScalar<Matrix<Scalar,Dynamic,1> > PlainADS;
00599   PlainADS ret;
00600   ret.value() = atan2(a.value(), b.value());
00601   
00602   Scalar tmp2 = a.value() * a.value();
00603   Scalar tmp3 = b.value() * b.value();
00604   Scalar tmp4 = tmp3/(tmp2+tmp3);
00605   
00606   if (tmp4!=0)
00607     ret.derivatives() = (a.derivatives() * b.value() - a.value() * b.derivatives()) * (tmp2+tmp3);
00608 
00609   return ret;
00610 }
00611 
00612 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(tan,
00613   using std::tan;
00614   using std::cos;
00615   return ReturnType(tan(x.value()),x.derivatives() * (Scalar(1)/numext::abs2(cos(x.value()))));)
00616 
00617 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(asin,
00618   using std::sqrt;
00619   using std::asin;
00620   return ReturnType(asin(x.value()),x.derivatives() * (Scalar(1)/sqrt(1-numext::abs2(x.value()))));)
00621   
00622 EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY(acos,
00623   using std::sqrt;
00624   using std::acos;
00625   return ReturnType(acos(x.value()),x.derivatives() * (Scalar(-1)/sqrt(1-numext::abs2(x.value()))));)
00626 
00627 #undef EIGEN_AUTODIFF_DECLARE_GLOBAL_UNARY
00628 
00629 template<typename DerType> struct NumTraits<AutoDiffScalar<DerType> >
00630   : NumTraits< typename NumTraits<typename DerType::Scalar>::Real >
00631 {
00632   typedef AutoDiffScalar<Matrix<typename NumTraits<typename DerType::Scalar>::Real,DerType::RowsAtCompileTime,DerType::ColsAtCompileTime> > Real;
00633   typedef AutoDiffScalar<DerType> NonInteger;
00634   typedef AutoDiffScalar<DerType>& Nested;
00635   enum{
00636     RequireInitialization = 1
00637   };
00638 };
00639 
00640 }
00641 
00642 #endif // EIGEN_AUTODIFF_SCALAR_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 11:57:50