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


libicr
Author(s): Robert Krug
autogenerated on Mon Jan 6 2014 11:32:28