aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Geometry
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-12-22 22:51:08 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-12-22 22:51:08 +0100
commiteaaba30cacf078335ad41314bf1e4ce993f0b3b4 (patch)
tree61156273b56dc19fc4e661a5aeb42d8e6112a307 /Eigen/src/Geometry
parente182e9c6166840b8db44e0be459a2e26d26fda0f (diff)
parenteeec7d842e05394703c42e7569230835f7842089 (diff)
merge with default branch
Diffstat (limited to 'Eigen/src/Geometry')
-rw-r--r--Eigen/src/Geometry/AlignedBox.h2
-rw-r--r--Eigen/src/Geometry/AngleAxis.h4
-rw-r--r--Eigen/src/Geometry/EulerAngles.h2
-rw-r--r--Eigen/src/Geometry/Homogeneous.h10
-rw-r--r--Eigen/src/Geometry/Hyperplane.h6
-rw-r--r--Eigen/src/Geometry/OrthoMethods.h5
-rw-r--r--Eigen/src/Geometry/ParametrizedLine.h2
-rw-r--r--Eigen/src/Geometry/Quaternion.h58
-rw-r--r--Eigen/src/Geometry/Rotation2D.h2
-rw-r--r--Eigen/src/Geometry/Scaling.h2
-rw-r--r--Eigen/src/Geometry/Transform.h4
-rw-r--r--Eigen/src/Geometry/Translation.h2
-rw-r--r--Eigen/src/Geometry/arch/Geometry_SSE.h6
13 files changed, 54 insertions, 51 deletions
diff --git a/Eigen/src/Geometry/AlignedBox.h b/Eigen/src/Geometry/AlignedBox.h
index 707cfeb4f..d47adb7ab 100644
--- a/Eigen/src/Geometry/AlignedBox.h
+++ b/Eigen/src/Geometry/AlignedBox.h
@@ -171,7 +171,7 @@ EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF_VECTORIZABLE_FIXED_SIZE(_Scalar,_AmbientDim)
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const AlignedBox& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const AlignedBox& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_min.isApprox(other.m_min, prec) && m_max.isApprox(other.m_max, prec); }
protected:
diff --git a/Eigen/src/Geometry/AngleAxis.h b/Eigen/src/Geometry/AngleAxis.h
index d948ea40b..da698bbfe 100644
--- a/Eigen/src/Geometry/AngleAxis.h
+++ b/Eigen/src/Geometry/AngleAxis.h
@@ -146,7 +146,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const AngleAxis& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const AngleAxis& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_axis.isApprox(other.m_axis, prec) && ei_isApprox(m_angle,other.m_angle, prec); }
};
@@ -165,7 +165,7 @@ template<typename QuatDerived>
AngleAxis<Scalar>& AngleAxis<Scalar>::operator=(const QuaternionBase<QuatDerived>& q)
{
Scalar n2 = q.vec().squaredNorm();
- if (n2 < precision<Scalar>()*precision<Scalar>())
+ if (n2 < dummy_precision<Scalar>()*dummy_precision<Scalar>())
{
m_angle = 0;
m_axis << 1, 0, 0;
diff --git a/Eigen/src/Geometry/EulerAngles.h b/Eigen/src/Geometry/EulerAngles.h
index dbcc7ae89..b6dbf8ae9 100644
--- a/Eigen/src/Geometry/EulerAngles.h
+++ b/Eigen/src/Geometry/EulerAngles.h
@@ -50,7 +50,7 @@ MatrixBase<Derived>::eulerAngles(int a0, int a1, int a2) const
Matrix<Scalar,3,1> res;
typedef Matrix<typename Derived::Scalar,2,1> Vector2;
- const Scalar epsilon = precision<Scalar>();
+ const Scalar epsilon = dummy_precision<Scalar>();
const int odd = ((a0+1)%3 == a1) ? 0 : 1;
const int i = a0;
diff --git a/Eigen/src/Geometry/Homogeneous.h b/Eigen/src/Geometry/Homogeneous.h
index 6cf6a381d..a601e29cf 100644
--- a/Eigen/src/Geometry/Homogeneous.h
+++ b/Eigen/src/Geometry/Homogeneous.h
@@ -176,7 +176,7 @@ MatrixBase<Derived>::hnormalized() const
EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived);
return StartMinusOne(derived(),0,0,
ColsAtCompileTime==1?size()-1:1,
- ColsAtCompileTime==1?1:size()-1).nestByValue() / coeff(size()-1);
+ ColsAtCompileTime==1?1:size()-1) / coeff(size()-1);
}
/** \geometry_module
@@ -193,8 +193,7 @@ VectorwiseOp<ExpressionType,Direction>::hnormalized() const
{
return HNormalized_Block(_expression(),0,0,
Direction==Vertical ? _expression().rows()-1 : _expression().rows(),
- Direction==Horizontal ? _expression().cols()-1 : _expression().cols()).nestByValue()
- .cwiseQuotient(
+ Direction==Horizontal ? _expression().cols()-1 : _expression().cols()).cwiseQuotient(
Replicate<NestByValue<HNormalized_Factors>,
Direction==Vertical ? HNormalized_SizeMinusOne : 1,
Direction==Horizontal ? HNormalized_SizeMinusOne : 1>
@@ -202,9 +201,9 @@ VectorwiseOp<ExpressionType,Direction>::hnormalized() const
Direction==Vertical ? _expression().rows()-1:0,
Direction==Horizontal ? _expression().cols()-1:0,
Direction==Vertical ? 1 : _expression().rows(),
- Direction==Horizontal ? 1 : _expression().cols()).nestByValue(),
+ Direction==Horizontal ? 1 : _expression().cols()),
Direction==Vertical ? _expression().rows()-1 : 1,
- Direction==Horizontal ? _expression().cols()-1 : 1).nestByValue());
+ Direction==Horizontal ? _expression().cols()-1 : 1));
}
template<typename MatrixType,typename Lhs>
@@ -281,7 +280,6 @@ struct ei_homogeneous_right_product_impl<Homogeneous<MatrixType,Horizontal>,Rhs>
const typename MatrixType::Nested m_lhs;
const typename Rhs::Nested m_rhs;
-
};
#endif // EIGEN_HOMOGENEOUS_H
diff --git a/Eigen/src/Geometry/Hyperplane.h b/Eigen/src/Geometry/Hyperplane.h
index ab65ca2ae..aab3d5b35 100644
--- a/Eigen/src/Geometry/Hyperplane.h
+++ b/Eigen/src/Geometry/Hyperplane.h
@@ -52,9 +52,9 @@ public:
typedef _Scalar Scalar;
typedef typename NumTraits<Scalar>::Real RealScalar;
typedef Matrix<Scalar,AmbientDimAtCompileTime,1> VectorType;
- typedef Matrix<Scalar,AmbientDimAtCompileTime==Dynamic
+ typedef Matrix<Scalar,int(AmbientDimAtCompileTime)==Dynamic
? Dynamic
- : AmbientDimAtCompileTime+1,1> Coefficients;
+ : int(AmbientDimAtCompileTime)+1,1> Coefficients;
typedef Block<Coefficients,AmbientDimAtCompileTime,1> NormalReturnType;
/** Default constructor without initialization */
@@ -257,7 +257,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const Hyperplane& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const Hyperplane& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_coeffs.isApprox(other.m_coeffs, prec); }
protected:
diff --git a/Eigen/src/Geometry/OrthoMethods.h b/Eigen/src/Geometry/OrthoMethods.h
index 631b82bf5..79baeadd5 100644
--- a/Eigen/src/Geometry/OrthoMethods.h
+++ b/Eigen/src/Geometry/OrthoMethods.h
@@ -90,8 +90,9 @@ MatrixBase<Derived>::cross3(const MatrixBase<OtherDerived>& other) const
const DerivedNested lhs(derived());
const OtherDerivedNested rhs(other.derived());
- return ei_cross3_impl<EiArch,typename ei_cleantype<DerivedNested>::type,
- typename ei_cleantype<OtherDerivedNested>::type>::run(lhs,rhs);
+ return ei_cross3_impl<Architecture::Target,
+ typename ei_cleantype<DerivedNested>::type,
+ typename ei_cleantype<OtherDerivedNested>::type>::run(lhs,rhs);
}
/** \returns a matrix expression of the cross product of each column or row
diff --git a/Eigen/src/Geometry/ParametrizedLine.h b/Eigen/src/Geometry/ParametrizedLine.h
index 6a4fcb1c5..21a5595b9 100644
--- a/Eigen/src/Geometry/ParametrizedLine.h
+++ b/Eigen/src/Geometry/ParametrizedLine.h
@@ -123,7 +123,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const ParametrizedLine& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_origin.isApprox(other.m_origin, prec) && m_direction.isApprox(other.m_direction, prec); }
protected:
diff --git a/Eigen/src/Geometry/Quaternion.h b/Eigen/src/Geometry/Quaternion.h
index 660e10d1e..861eff19c 100644
--- a/Eigen/src/Geometry/Quaternion.h
+++ b/Eigen/src/Geometry/Quaternion.h
@@ -152,9 +152,9 @@ public:
/** \returns the conjugated quaternion */
Quaternion<Scalar> conjugate() const;
- /** \returns an interpolation for a constant motion between \a other and \c *this
+ /** \returns an interpolation for a constant motion between \a other and \c *this
* \a t in [0;1]
- * see http://en.wikipedia.org/wiki/Slerp
+ * see http://en.wikipedia.org/wiki/Slerp
*/
template<class OtherDerived> Quaternion<Scalar> slerp(Scalar t, const QuaternionBase<OtherDerived>& other) const;
@@ -163,7 +163,7 @@ public:
*
* \sa MatrixBase::isApprox() */
template<class OtherDerived>
- bool isApprox(const QuaternionBase<OtherDerived>& other, RealScalar prec = precision<Scalar>()) const
+ bool isApprox(const QuaternionBase<OtherDerived>& other, RealScalar prec = dummy_precision<Scalar>()) const
{ return coeffs().isApprox(other.coeffs(), prec); }
/** return the result vector of \a v through the rotation*/
@@ -221,7 +221,7 @@ struct ei_traits<Quaternion<_Scalar> >
template<typename _Scalar>
class Quaternion : public QuaternionBase<Quaternion<_Scalar> >{
typedef QuaternionBase<Quaternion<_Scalar> > Base;
-public:
+public:
typedef _Scalar Scalar;
EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Quaternion<Scalar>)
@@ -304,27 +304,20 @@ template<typename _Scalar, int PacketAccess>
class Map<Quaternion<_Scalar>, PacketAccess >
: public QuaternionBase<Map<Quaternion<_Scalar>, PacketAccess> >
{
- public:
- typedef _Scalar Scalar;
- typedef Map<Quaternion<Scalar>, PacketAccess > MapQuat;
-
- private:
- Map<Quaternion<Scalar>, PacketAccess >();
- Map<Quaternion<Scalar>, PacketAccess >(const Map<Quaternion<Scalar>, PacketAccess>&);
-
typedef QuaternionBase<Map<Quaternion<_Scalar>, PacketAccess> > Base;
+
public:
- EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(MapQuat)
+ typedef _Scalar Scalar;
+ typedef typename ei_traits<Map>::Coefficients Coefficients;
+ EIGEN_INHERIT_ASSIGNMENT_EQUAL_OPERATOR(Map)
using Base::operator*=;
- typedef typename ei_traits<Map<Quaternion<Scalar>, PacketAccess> >::Coefficients Coefficients;
-
/** Constructs a Mapped Quaternion object from the pointer \a coeffs
*
* The pointer \a coeffs must reference the four coeffecients of Quaternion in the following order:
* \code *coeffs == {x, y, z, w} \endcode
*
- * If the template paramter PacketAccess is set to Aligned, then the pointer coeffs must be aligned. */
+ * If the template parameter PacketAccess is set to Aligned, then the pointer coeffs must be aligned. */
EIGEN_STRONG_INLINE Map(const Scalar* coeffs) : m_coeffs(coeffs) {}
inline Coefficients& coeffs() { return m_coeffs;}
@@ -374,7 +367,7 @@ QuaternionBase<Derived>::operator* (const QuaternionBase<OtherDerived>& other) c
{
EIGEN_STATIC_ASSERT((ei_is_same_type<typename Derived::Scalar, typename OtherDerived::Scalar>::ret),
YOU_MIXED_DIFFERENT_NUMERIC_TYPES__YOU_NEED_TO_USE_THE_CAST_METHOD_OF_MATRIXBASE_TO_CAST_NUMERIC_TYPES_EXPLICITLY)
- return ei_quat_product<EiArch, Derived, OtherDerived,
+ return ei_quat_product<Architecture::Target, Derived, OtherDerived,
typename ei_traits<Derived>::Scalar,
ei_traits<Derived>::PacketAccess && ei_traits<OtherDerived>::PacketAccess>::run(*this, other);
}
@@ -514,7 +507,7 @@ inline Derived& QuaternionBase<Derived>::setFromTwoVectors(const MatrixBase<Deri
// under the constraint:
// ||x|| = 1
// which yields a singular value problem
- if (c < Scalar(-1)+precision<Scalar>())
+ if (c < Scalar(-1)+dummy_precision<Scalar>())
{
c = std::max<Scalar>(c,-1);
Matrix<Scalar,2,3> m; m << v0.transpose(), v1.transpose();
@@ -590,20 +583,29 @@ template <class OtherDerived>
Quaternion<typename ei_traits<Derived>::Scalar>
QuaternionBase<Derived>::slerp(Scalar t, const QuaternionBase<OtherDerived>& other) const
{
- static const Scalar one = Scalar(1) - precision<Scalar>();
+ static const Scalar one = Scalar(1) - epsilon<Scalar>();
Scalar d = this->dot(other);
Scalar absD = ei_abs(d);
- if (absD>=one)
- return Quaternion<Scalar>(derived());
- // theta is the angle between the 2 quaternions
- Scalar theta = std::acos(absD);
- Scalar sinTheta = ei_sin(theta);
+ Scalar scale0;
+ Scalar scale1;
- Scalar scale0 = ei_sin( ( Scalar(1) - t ) * theta) / sinTheta;
- Scalar scale1 = ei_sin( ( t * theta) ) / sinTheta;
- if (d<0)
- scale1 = -scale1;
+ if (absD>=one)
+ {
+ scale0 = Scalar(1) - t;
+ scale1 = t;
+ }
+ else
+ {
+ // theta is the angle between the 2 quaternions
+ Scalar theta = std::acos(absD);
+ Scalar sinTheta = ei_sin(theta);
+
+ scale0 = ei_sin( ( Scalar(1) - t ) * theta) / sinTheta;
+ scale1 = ei_sin( ( t * theta) ) / sinTheta;
+ if (d<0)
+ scale1 = -scale1;
+ }
return Quaternion<Scalar>(scale0 * coeffs() + scale1 * other.coeffs());
}
diff --git a/Eigen/src/Geometry/Rotation2D.h b/Eigen/src/Geometry/Rotation2D.h
index b80fcace0..7f24a3eae 100644
--- a/Eigen/src/Geometry/Rotation2D.h
+++ b/Eigen/src/Geometry/Rotation2D.h
@@ -121,7 +121,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const Rotation2D& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const Rotation2D& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return ei_isApprox(m_angle,other.m_angle, prec); }
};
diff --git a/Eigen/src/Geometry/Scaling.h b/Eigen/src/Geometry/Scaling.h
index ce191b5da..4695914fd 100644
--- a/Eigen/src/Geometry/Scaling.h
+++ b/Eigen/src/Geometry/Scaling.h
@@ -107,7 +107,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const UniformScaling& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const UniformScaling& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return ei_isApprox(m_factor, other.factor(), prec); }
};
diff --git a/Eigen/src/Geometry/Transform.h b/Eigen/src/Geometry/Transform.h
index cf961b98b..89df73505 100644
--- a/Eigen/src/Geometry/Transform.h
+++ b/Eigen/src/Geometry/Transform.h
@@ -187,7 +187,7 @@ public:
/** type of read/write reference to the affine part of the transformation */
typedef typename ei_meta_if<int(Mode)==int(AffineCompact),
MatrixType&,
- NestByValue<Block<MatrixType,Dim,HDim> > >::ret AffinePartNested;
+ Block<MatrixType,Dim,HDim> >::ret AffinePartNested;
/** type of a vector */
typedef Matrix<Scalar,Dim,1> VectorType;
/** type of a read/write reference to the translation part of the rotation */
@@ -424,7 +424,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const Transform& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const Transform& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_matrix.isApprox(other.m_matrix, prec); }
/** Sets the last row to [0 ... 0 1]
diff --git a/Eigen/src/Geometry/Translation.h b/Eigen/src/Geometry/Translation.h
index 1fff03810..b7477df9f 100644
--- a/Eigen/src/Geometry/Translation.h
+++ b/Eigen/src/Geometry/Translation.h
@@ -154,7 +154,7 @@ public:
* determined by \a prec.
*
* \sa MatrixBase::isApprox() */
- bool isApprox(const Translation& other, typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
+ bool isApprox(const Translation& other, typename NumTraits<Scalar>::Real prec = dummy_precision<Scalar>()) const
{ return m_coeffs.isApprox(other.m_coeffs, prec); }
};
diff --git a/Eigen/src/Geometry/arch/Geometry_SSE.h b/Eigen/src/Geometry/arch/Geometry_SSE.h
index a6ed10d82..297932f92 100644
--- a/Eigen/src/Geometry/arch/Geometry_SSE.h
+++ b/Eigen/src/Geometry/arch/Geometry_SSE.h
@@ -26,7 +26,8 @@
#ifndef EIGEN_GEOMETRY_SSE_H
#define EIGEN_GEOMETRY_SSE_H
-template<class Derived, class OtherDerived> struct ei_quat_product<EiArch_SSE, Derived, OtherDerived, float, Aligned>
+template<class Derived, class OtherDerived>
+struct ei_quat_product<Architecture::SSE, Derived, OtherDerived, float, Aligned>
{
inline static Quaternion<float> run(const QuaternionBase<Derived>& _a, const QuaternionBase<OtherDerived>& _b)
{
@@ -48,7 +49,8 @@ template<class Derived, class OtherDerived> struct ei_quat_product<EiArch_SSE, D
};
template<typename VectorLhs,typename VectorRhs>
-struct ei_cross3_impl<EiArch_SSE,VectorLhs,VectorRhs,float,true> {
+struct ei_cross3_impl<Architecture::SSE,VectorLhs,VectorRhs,float,true>
+{
inline static typename ei_plain_matrix_type<VectorLhs>::type
run(const VectorLhs& lhs, const VectorRhs& rhs)
{