AngleAxis.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) 2008 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_ANGLEAXIS_H
00026 #define EIGEN_ANGLEAXIS_H
00027 
00028 namespace Eigen { 
00029 
00056 namespace internal {
00057 template<typename _Scalar> struct traits<AngleAxis<_Scalar> >
00058 {
00059   typedef _Scalar Scalar;
00060 };
00061 }
00062 
00063 template<typename _Scalar>
00064 class AngleAxis : public RotationBase<AngleAxis<_Scalar>,3>
00065 {
00066   typedef RotationBase<AngleAxis<_Scalar>,3> Base;
00067 
00068 public:
00069 
00070   using Base::operator*;
00071 
00072   enum { Dim = 3 };
00074   typedef _Scalar Scalar;
00075   typedef Matrix<Scalar,3,3> Matrix3;
00076   typedef Matrix<Scalar,3,1> Vector3;
00077   typedef Quaternion<Scalar> QuaternionType;
00078 
00079 protected:
00080 
00081   Vector3 m_axis;
00082   Scalar m_angle;
00083 
00084 public:
00085 
00087   AngleAxis() {}
00093   template<typename Derived>
00094   inline AngleAxis(Scalar angle, const MatrixBase<Derived>& axis) : m_axis(axis), m_angle(angle) {}
00096   template<typename QuatDerived> inline explicit AngleAxis(const QuaternionBase<QuatDerived>& q) { *this = q; }
00098   template<typename Derived>
00099   inline explicit AngleAxis(const MatrixBase<Derived>& m) { *this = m; }
00100 
00101   Scalar angle() const { return m_angle; }
00102   Scalar& angle() { return m_angle; }
00103 
00104   const Vector3& axis() const { return m_axis; }
00105   Vector3& axis() { return m_axis; }
00106 
00108   inline QuaternionType operator* (const AngleAxis& other) const
00109   { return QuaternionType(*this) * QuaternionType(other); }
00110 
00112   inline QuaternionType operator* (const QuaternionType& other) const
00113   { return QuaternionType(*this) * other; }
00114 
00116   friend inline QuaternionType operator* (const QuaternionType& a, const AngleAxis& b)
00117   { return a * QuaternionType(b); }
00118 
00120   AngleAxis inverse() const
00121   { return AngleAxis(-m_angle, m_axis); }
00122 
00123   template<class QuatDerived>
00124   AngleAxis& operator=(const QuaternionBase<QuatDerived>& q);
00125   template<typename Derived>
00126   AngleAxis& operator=(const MatrixBase<Derived>& m);
00127 
00128   template<typename Derived>
00129   AngleAxis& fromRotationMatrix(const MatrixBase<Derived>& m);
00130   Matrix3 toRotationMatrix(void) const;
00131 
00137   template<typename NewScalarType>
00138   inline typename internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type cast() const
00139   { return typename internal::cast_return_type<AngleAxis,AngleAxis<NewScalarType> >::type(*this); }
00140 
00142   template<typename OtherScalarType>
00143   inline explicit AngleAxis(const AngleAxis<OtherScalarType>& other)
00144   {
00145     m_axis = other.axis().template cast<Scalar>();
00146     m_angle = Scalar(other.angle());
00147   }
00148 
00149   static inline const AngleAxis Identity() { return AngleAxis(0, Vector3::UnitX()); }
00150 
00155   bool isApprox(const AngleAxis& other, typename NumTraits<Scalar>::Real prec = NumTraits<Scalar>::dummy_precision()) const
00156   { return m_axis.isApprox(other.m_axis, prec) && internal::isApprox(m_angle,other.m_angle, prec); }
00157 };
00158 
00161 typedef AngleAxis<float> AngleAxisf;
00164 typedef AngleAxis<double> AngleAxisd;
00165 
00172 template<typename Scalar>
00173 template<typename QuatDerived>
00174 AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const QuaternionBase<QuatDerived>& q)
00175 {
00176   using std::acos;
00177   using std::min;
00178   using std::max;
00179   Scalar n2 = q.vec().squaredNorm();
00180   if (n2 < NumTraits<Scalar>::dummy_precision()*NumTraits<Scalar>::dummy_precision())
00181   {
00182     m_angle = 0;
00183     m_axis << 1, 0, 0;
00184   }
00185   else
00186   {
00187     m_angle = Scalar(2)*acos((min)((max)(Scalar(-1),q.w()),Scalar(1)));
00188     m_axis = q.vec() / internal::sqrt(n2);
00189   }
00190   return *this;
00191 }
00192 
00195 template<typename Scalar>
00196 template<typename Derived>
00197 AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const MatrixBase<Derived>& mat)
00198 {
00199   // Since a direct conversion would not be really faster,
00200   // let's use the robust Quaternion implementation:
00201   return *this = QuaternionType(mat);
00202 }
00203 
00207 template<typename Scalar>
00208 template<typename Derived>
00209 AngleAxis<Scalar>& AngleAxis<Scalar>::fromRotationMatrix(const MatrixBase<Derived>& mat)
00210 {
00211   return *this = QuaternionType(mat);
00212 }
00213 
00216 template<typename Scalar>
00217 typename AngleAxis<Scalar>::Matrix3
00218 AngleAxis<Scalar>::toRotationMatrix(void) const
00219 {
00220   Matrix3 res;
00221   Vector3 sin_axis  = internal::sin(m_angle) * m_axis;
00222   Scalar c = internal::cos(m_angle);
00223   Vector3 cos1_axis = (Scalar(1)-c) * m_axis;
00224 
00225   Scalar tmp;
00226   tmp = cos1_axis.x() * m_axis.y();
00227   res.coeffRef(0,1) = tmp - sin_axis.z();
00228   res.coeffRef(1,0) = tmp + sin_axis.z();
00229 
00230   tmp = cos1_axis.x() * m_axis.z();
00231   res.coeffRef(0,2) = tmp + sin_axis.y();
00232   res.coeffRef(2,0) = tmp - sin_axis.y();
00233 
00234   tmp = cos1_axis.y() * m_axis.z();
00235   res.coeffRef(1,2) = tmp - sin_axis.x();
00236   res.coeffRef(2,1) = tmp + sin_axis.x();
00237 
00238   res.diagonal() = (cos1_axis.cwiseProduct(m_axis)).array() + c;
00239 
00240   return res;
00241 }
00242 
00243 } // end namespace Eigen
00244 
00245 #endif // EIGEN_ANGLEAXIS_H