aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Eigen2Support
diff options
context:
space:
mode:
authorGravatar Igor Krivenko <anonymous@invalid.net>2011-12-09 23:38:41 +0100
committerGravatar Igor Krivenko <anonymous@invalid.net>2011-12-09 23:38:41 +0100
commit36457178f9d12481e699d07c0ac8d25a7e6a854f (patch)
treee9b12b3ed7e4102fcf1d0c4384fec3c1d4c94f37 /Eigen/src/Eigen2Support
parentd400a6245ec0d86b0e3dd2af30e49349ae28b0f0 (diff)
bug #352:properly cast constants
Diffstat (limited to 'Eigen/src/Eigen2Support')
-rw-r--r--Eigen/src/Eigen2Support/Geometry/AlignedBox.h6
-rw-r--r--Eigen/src/Eigen2Support/Geometry/Quaternion.h12
-rw-r--r--Eigen/src/Eigen2Support/SVD.h2
3 files changed, 10 insertions, 10 deletions
diff --git a/Eigen/src/Eigen2Support/Geometry/AlignedBox.h b/Eigen/src/Eigen2Support/Geometry/AlignedBox.h
index 2116f9c83..e4e555715 100644
--- a/Eigen/src/Eigen2Support/Geometry/AlignedBox.h
+++ b/Eigen/src/Eigen2Support/Geometry/AlignedBox.h
@@ -157,13 +157,13 @@ protected:
template<typename Scalar,int AmbiantDim>
inline Scalar AlignedBox<Scalar,AmbiantDim>::squaredExteriorDistance(const VectorType& p) const
{
- Scalar dist2 = 0.;
+ Scalar dist2(0);
Scalar aux;
for (int k=0; k<dim(); ++k)
{
- if ((aux = (p[k]-m_min[k]))<0.)
+ if ((aux = (p[k]-m_min[k]))<Scalar(0))
dist2 += aux*aux;
- else if ( (aux = (m_max[k]-p[k]))<0. )
+ else if ( (aux = (m_max[k]-p[k]))<Scalar(0))
dist2 += aux*aux;
}
return dist2;
diff --git a/Eigen/src/Eigen2Support/Geometry/Quaternion.h b/Eigen/src/Eigen2Support/Geometry/Quaternion.h
index a75fa42ae..bb6326416 100644
--- a/Eigen/src/Eigen2Support/Geometry/Quaternion.h
+++ b/Eigen/src/Eigen2Support/Geometry/Quaternion.h
@@ -314,9 +314,9 @@ Quaternion<Scalar>::toRotationMatrix(void) const
// it has to be inlined, and so the return by value is not an issue
Matrix3 res;
- const Scalar tx = 2*this->x();
- const Scalar ty = 2*this->y();
- const Scalar tz = 2*this->z();
+ const Scalar tx = Scalar(2)*this->x();
+ const Scalar ty = Scalar(2)*this->y();
+ const Scalar tz = Scalar(2)*this->z();
const Scalar twx = tx*this->w();
const Scalar twy = ty*this->w();
const Scalar twz = tz*this->w();
@@ -327,15 +327,15 @@ Quaternion<Scalar>::toRotationMatrix(void) const
const Scalar tyz = tz*this->y();
const Scalar tzz = tz*this->z();
- res.coeffRef(0,0) = 1-(tyy+tzz);
+ res.coeffRef(0,0) = Scalar(1)-(tyy+tzz);
res.coeffRef(0,1) = txy-twz;
res.coeffRef(0,2) = txz+twy;
res.coeffRef(1,0) = txy+twz;
- res.coeffRef(1,1) = 1-(txx+tzz);
+ res.coeffRef(1,1) = Scalar(1)-(txx+tzz);
res.coeffRef(1,2) = tyz-twx;
res.coeffRef(2,0) = txz-twy;
res.coeffRef(2,1) = tyz+twx;
- res.coeffRef(2,2) = 1-(txx+tyy);
+ res.coeffRef(2,2) = Scalar(1)-(txx+tyy);
return res;
}
diff --git a/Eigen/src/Eigen2Support/SVD.h b/Eigen/src/Eigen2Support/SVD.h
index 16b4b488f..44b72850f 100644
--- a/Eigen/src/Eigen2Support/SVD.h
+++ b/Eigen/src/Eigen2Support/SVD.h
@@ -390,7 +390,7 @@ void SVD<MatrixType>::compute(const MatrixType& matrix)
Scalar ek = e[k]/scale;
Scalar b = ((spm1 + sp)*(spm1 - sp) + epm1*epm1)/Scalar(2);
Scalar c = (sp*epm1)*(sp*epm1);
- Scalar shift = 0.0;
+ Scalar shift(0);
if ((b != 0.0) || (c != 0.0))
{
shift = ei_sqrt(b*b + c);