AutoDiffVector.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 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_AUTODIFF_VECTOR_H
00026 #define EIGEN_AUTODIFF_VECTOR_H
00027 
00028 namespace Eigen {
00029 
00030 /* \class AutoDiffScalar
00031   * \brief A scalar type replacement with automatic differentation capability
00032   *
00033   * \param DerType the vector type used to store/represent the derivatives (e.g. Vector3f)
00034   *
00035   * This class represents a scalar value while tracking its respective derivatives.
00036   *
00037   * It supports the following list of global math function:
00038   *  - std::abs, std::sqrt, std::pow, std::exp, std::log, std::sin, std::cos,
00039   *  - internal::abs, internal::sqrt, internal::pow, internal::exp, internal::log, internal::sin, internal::cos,
00040   *  - internal::conj, internal::real, internal::imag, internal::abs2.
00041   *
00042   * AutoDiffScalar can be used as the scalar type of an Eigen::Matrix object. However,
00043   * in that case, the expression template mechanism only occurs at the top Matrix level,
00044   * while derivatives are computed right away.
00045   *
00046   */
00047 template<typename ValueType, typename JacobianType>
00048 class AutoDiffVector
00049 {
00050   public:
00051     //typedef typename internal::traits<ValueType>::Scalar Scalar;
00052     typedef typename internal::traits<ValueType>::Scalar BaseScalar;
00053     typedef AutoDiffScalar<Matrix<BaseScalar,JacobianType::RowsAtCompileTime,1> > ActiveScalar;
00054     typedef ActiveScalar Scalar;
00055     typedef AutoDiffScalar<typename JacobianType::ColXpr> CoeffType;
00056     typedef typename JacobianType::Index Index;
00057 
00058     inline AutoDiffVector() {}
00059 
00060     inline AutoDiffVector(const ValueType& values)
00061       : m_values(values)
00062     {
00063       m_jacobian.setZero();
00064     }
00065 
00066 
00067     CoeffType operator[] (Index i) { return CoeffType(m_values[i], m_jacobian.col(i)); }
00068     const CoeffType operator[] (Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
00069 
00070     CoeffType operator() (Index i) { return CoeffType(m_values[i], m_jacobian.col(i)); }
00071     const CoeffType operator() (Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
00072 
00073     CoeffType coeffRef(Index i) { return CoeffType(m_values[i], m_jacobian.col(i)); }
00074     const CoeffType coeffRef(Index i) const { return CoeffType(m_values[i], m_jacobian.col(i)); }
00075 
00076     Index size() const { return m_values.size(); }
00077 
00078     // FIXME here we could return an expression of the sum
00079     Scalar sum() const { /*std::cerr << "sum \n\n";*/ /*std::cerr << m_jacobian.rowwise().sum() << "\n\n";*/ return Scalar(m_values.sum(), m_jacobian.rowwise().sum()); }
00080 
00081 
00082     inline AutoDiffVector(const ValueType& values, const JacobianType& jac)
00083       : m_values(values), m_jacobian(jac)
00084     {}
00085 
00086     template<typename OtherValueType, typename OtherJacobianType>
00087     inline AutoDiffVector(const AutoDiffVector<OtherValueType, OtherJacobianType>& other)
00088       : m_values(other.values()), m_jacobian(other.jacobian())
00089     {}
00090 
00091     inline AutoDiffVector(const AutoDiffVector& other)
00092       : m_values(other.values()), m_jacobian(other.jacobian())
00093     {}
00094 
00095     template<typename OtherValueType, typename OtherJacobianType>
00096     inline AutoDiffVector& operator=(const AutoDiffVector<OtherValueType, OtherJacobianType>& other)
00097     {
00098       m_values = other.values();
00099       m_jacobian = other.jacobian();
00100       return *this;
00101     }
00102 
00103     inline AutoDiffVector& operator=(const AutoDiffVector& other)
00104     {
00105       m_values = other.values();
00106       m_jacobian = other.jacobian();
00107       return *this;
00108     }
00109 
00110     inline const ValueType& values() const { return m_values; }
00111     inline ValueType& values() { return m_values; }
00112 
00113     inline const JacobianType& jacobian() const { return m_jacobian; }
00114     inline JacobianType& jacobian() { return m_jacobian; }
00115 
00116     template<typename OtherValueType,typename OtherJacobianType>
00117     inline const AutoDiffVector<
00118       typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,ValueType,OtherValueType>::Type,
00119       typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,JacobianType,OtherJacobianType>::Type >
00120     operator+(const AutoDiffVector<OtherValueType,OtherJacobianType>& other) const
00121     {
00122       return AutoDiffVector<
00123       typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,ValueType,OtherValueType>::Type,
00124       typename MakeCwiseBinaryOp<internal::scalar_sum_op<BaseScalar>,JacobianType,OtherJacobianType>::Type >(
00125         m_values + other.values(),
00126         m_jacobian + other.jacobian());
00127     }
00128 
00129     template<typename OtherValueType, typename OtherJacobianType>
00130     inline AutoDiffVector&
00131     operator+=(const AutoDiffVector<OtherValueType,OtherJacobianType>& other)
00132     {
00133       m_values += other.values();
00134       m_jacobian += other.jacobian();
00135       return *this;
00136     }
00137 
00138     template<typename OtherValueType,typename OtherJacobianType>
00139     inline const AutoDiffVector<
00140       typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,ValueType,OtherValueType>::Type,
00141       typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,JacobianType,OtherJacobianType>::Type >
00142     operator-(const AutoDiffVector<OtherValueType,OtherJacobianType>& other) const
00143     {
00144       return AutoDiffVector<
00145         typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,ValueType,OtherValueType>::Type,
00146         typename MakeCwiseBinaryOp<internal::scalar_difference_op<Scalar>,JacobianType,OtherJacobianType>::Type >(
00147           m_values - other.values(),
00148           m_jacobian - other.jacobian());
00149     }
00150 
00151     template<typename OtherValueType, typename OtherJacobianType>
00152     inline AutoDiffVector&
00153     operator-=(const AutoDiffVector<OtherValueType,OtherJacobianType>& other)
00154     {
00155       m_values -= other.values();
00156       m_jacobian -= other.jacobian();
00157       return *this;
00158     }
00159 
00160     inline const AutoDiffVector<
00161       typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, ValueType>::Type,
00162       typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type >
00163     operator-() const
00164     {
00165       return AutoDiffVector<
00166         typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, ValueType>::Type,
00167         typename MakeCwiseUnaryOp<internal::scalar_opposite_op<Scalar>, JacobianType>::Type >(
00168           -m_values,
00169           -m_jacobian);
00170     }
00171 
00172     inline const AutoDiffVector<
00173       typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
00174       typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type>
00175     operator*(const BaseScalar& other) const
00176     {
00177       return AutoDiffVector<
00178         typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
00179         typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >(
00180           m_values * other,
00181           m_jacobian * other);
00182     }
00183 
00184     friend inline const AutoDiffVector<
00185       typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
00186       typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >
00187     operator*(const Scalar& other, const AutoDiffVector& v)
00188     {
00189       return AutoDiffVector<
00190         typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, ValueType>::Type,
00191         typename MakeCwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>::Type >(
00192           v.values() * other,
00193           v.jacobian() * other);
00194     }
00195 
00196 //     template<typename OtherValueType,typename OtherJacobianType>
00197 //     inline const AutoDiffVector<
00198 //       CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
00199 //       CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
00200 //         CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
00201 //         CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >
00202 //     operator*(const AutoDiffVector<OtherValueType,OtherJacobianType>& other) const
00203 //     {
00204 //       return AutoDiffVector<
00205 //         CwiseBinaryOp<internal::scalar_multiple_op<Scalar>, ValueType, OtherValueType>
00206 //         CwiseBinaryOp<internal::scalar_sum_op<Scalar>,
00207 //           CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, JacobianType>,
00208 //           CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, OtherJacobianType> > >(
00209 //             m_values.cwise() * other.values(),
00210 //             (m_jacobian * other.values()) + (m_values * other.jacobian()));
00211 //     }
00212 
00213     inline AutoDiffVector& operator*=(const Scalar& other)
00214     {
00215       m_values *= other;
00216       m_jacobian *= other;
00217       return *this;
00218     }
00219 
00220     template<typename OtherValueType,typename OtherJacobianType>
00221     inline AutoDiffVector& operator*=(const AutoDiffVector<OtherValueType,OtherJacobianType>& other)
00222     {
00223       *this = *this * other;
00224       return *this;
00225     }
00226 
00227   protected:
00228     ValueType m_values;
00229     JacobianType m_jacobian;
00230 
00231 };
00232 
00233 }
00234 
00235 #endif // EIGEN_AUTODIFF_VECTOR_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:30:45