aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/src/NonLinearOptimization
diff options
context:
space:
mode:
authorGravatar Thomas Capricelli <orzel@freehackers.org>2010-02-11 01:43:03 +0100
committerGravatar Thomas Capricelli <orzel@freehackers.org>2010-02-11 01:43:03 +0100
commitccc58e935ea2123590013f0cdc8e6ae25bca4f28 (patch)
treebf855673783a71acb65d07b189d26f5a1b4e7eba /unsupported/Eigen/src/NonLinearOptimization
parentfc5fa7774340da06cb4de7337c0928e76dfe4946 (diff)
fix usage of epsilon wrt to latest API change
Diffstat (limited to 'unsupported/Eigen/src/NonLinearOptimization')
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/HybridNonLinearSolver.h10
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/LevenbergMarquardt.h22
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/chkder.h6
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/covar.h2
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/dogleg.h2
-rw-r--r--unsupported/Eigen/src/NonLinearOptimization/fdjac1.h2
6 files changed, 22 insertions, 22 deletions
diff --git a/unsupported/Eigen/src/NonLinearOptimization/HybridNonLinearSolver.h b/unsupported/Eigen/src/NonLinearOptimization/HybridNonLinearSolver.h
index 5c832b73d..0754a426b 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/HybridNonLinearSolver.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/HybridNonLinearSolver.h
@@ -63,7 +63,7 @@ public:
Parameters()
: factor(Scalar(100.))
, maxfev(1000)
- , xtol(ei_sqrt(epsilon<Scalar>()))
+ , xtol(ei_sqrt(NumTraits<Scalar>::epsilon()))
, nb_of_subdiagonals(-1)
, nb_of_superdiagonals(-1)
, epsfcn(Scalar(0.)) {}
@@ -81,7 +81,7 @@ public:
HybridNonLinearSolverSpace::Status hybrj1(
FVectorType &x,
- const Scalar tol = ei_sqrt(epsilon<Scalar>())
+ const Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon())
);
HybridNonLinearSolverSpace::Status solveInit(
@@ -99,7 +99,7 @@ public:
HybridNonLinearSolverSpace::Status hybrd1(
FVectorType &x,
- const Scalar tol = ei_sqrt(epsilon<Scalar>())
+ const Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon())
);
HybridNonLinearSolverSpace::Status solveNumericalDiffInit(
@@ -341,7 +341,7 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveOneStep(
/* tests for termination and stringent tolerances. */
if (nfev >= parameters.maxfev)
return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
- if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= epsilon<Scalar>() * xnorm)
+ if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
return HybridNonLinearSolverSpace::TolTooSmall;
if (nslow2 == 5)
return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
@@ -590,7 +590,7 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(
/* tests for termination and stringent tolerances. */
if (nfev >= parameters.maxfev)
return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
- if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= epsilon<Scalar>() * xnorm)
+ if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
return HybridNonLinearSolverSpace::TolTooSmall;
if (nslow2 == 5)
return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
diff --git a/unsupported/Eigen/src/NonLinearOptimization/LevenbergMarquardt.h b/unsupported/Eigen/src/NonLinearOptimization/LevenbergMarquardt.h
index 49b63bed8..f21331e3e 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/LevenbergMarquardt.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/LevenbergMarquardt.h
@@ -67,8 +67,8 @@ public:
Parameters()
: factor(Scalar(100.))
, maxfev(400)
- , ftol(ei_sqrt(epsilon<Scalar>()))
- , xtol(ei_sqrt(epsilon<Scalar>()))
+ , ftol(ei_sqrt(NumTraits<Scalar>::epsilon()))
+ , xtol(ei_sqrt(NumTraits<Scalar>::epsilon()))
, gtol(Scalar(0.))
, epsfcn(Scalar(0.)) {}
Scalar factor;
@@ -84,7 +84,7 @@ public:
LevenbergMarquardtSpace::Status lmder1(
FVectorType &x,
- const Scalar tol = ei_sqrt(epsilon<Scalar>())
+ const Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon())
);
LevenbergMarquardtSpace::Status minimize(
@@ -104,12 +104,12 @@ public:
FunctorType &functor,
FVectorType &x,
int *nfev,
- const Scalar tol = ei_sqrt(epsilon<Scalar>())
+ const Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon())
);
LevenbergMarquardtSpace::Status lmstr1(
FVectorType &x,
- const Scalar tol = ei_sqrt(epsilon<Scalar>())
+ const Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon())
);
LevenbergMarquardtSpace::Status minimizeOptimumStorage(
@@ -370,11 +370,11 @@ LevenbergMarquardt<FunctorType,Scalar>::minimizeOneStep(
/* tests for termination and stringent tolerances. */
if (nfev >= parameters.maxfev)
return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
- if (ei_abs(actred) <= epsilon<Scalar>() && prered <= epsilon<Scalar>() && Scalar(.5) * ratio <= 1.)
+ if (ei_abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
return LevenbergMarquardtSpace::FtolTooSmall;
- if (delta <= epsilon<Scalar>() * xnorm)
+ if (delta <= NumTraits<Scalar>::epsilon() * xnorm)
return LevenbergMarquardtSpace::XtolTooSmall;
- if (gnorm <= epsilon<Scalar>())
+ if (gnorm <= NumTraits<Scalar>::epsilon())
return LevenbergMarquardtSpace::GtolTooSmall;
} while (ratio < Scalar(1e-4));
@@ -621,11 +621,11 @@ LevenbergMarquardt<FunctorType,Scalar>::minimizeOptimumStorageOneStep(
/* tests for termination and stringent tolerances. */
if (nfev >= parameters.maxfev)
return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
- if (ei_abs(actred) <= epsilon<Scalar>() && prered <= epsilon<Scalar>() && Scalar(.5) * ratio <= 1.)
+ if (ei_abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
return LevenbergMarquardtSpace::FtolTooSmall;
- if (delta <= epsilon<Scalar>() * xnorm)
+ if (delta <= NumTraits<Scalar>::epsilon() * xnorm)
return LevenbergMarquardtSpace::XtolTooSmall;
- if (gnorm <= epsilon<Scalar>())
+ if (gnorm <= NumTraits<Scalar>::epsilon())
return LevenbergMarquardtSpace::GtolTooSmall;
} while (ratio < Scalar(1e-4));
diff --git a/unsupported/Eigen/src/NonLinearOptimization/chkder.h b/unsupported/Eigen/src/NonLinearOptimization/chkder.h
index be8115230..591e8bef7 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/chkder.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/chkder.h
@@ -13,8 +13,8 @@ void ei_chkder(
Matrix< Scalar, Dynamic, 1 > &err
)
{
- const Scalar eps = ei_sqrt(epsilon<Scalar>());
- const Scalar epsf = chkder_factor * epsilon<Scalar>();
+ const Scalar eps = ei_sqrt(NumTraits<Scalar>::epsilon());
+ const Scalar epsf = chkder_factor * NumTraits<Scalar>::epsilon();
const Scalar epslog = chkder_log10e * ei_log(eps);
Scalar temp;
@@ -44,7 +44,7 @@ void ei_chkder(
if (fvec[i] != 0. && fvecp[i] != 0. && ei_abs(fvecp[i] - fvec[i]) >= epsf * ei_abs(fvec[i]))
temp = eps * ei_abs((fvecp[i] - fvec[i]) / eps - err[i]) / (ei_abs(fvec[i]) + ei_abs(fvecp[i]));
err[i] = 1.;
- if (temp > epsilon<Scalar>() && temp < eps)
+ if (temp > NumTraits<Scalar>::epsilon() && temp < eps)
err[i] = (chkder_log10e * ei_log(temp) - epslog) / epslog;
if (temp >= eps)
err[i] = 0.;
diff --git a/unsupported/Eigen/src/NonLinearOptimization/covar.h b/unsupported/Eigen/src/NonLinearOptimization/covar.h
index 15c72d6ed..2f983a958 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/covar.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/covar.h
@@ -3,7 +3,7 @@ template <typename Scalar>
void ei_covar(
Matrix< Scalar, Dynamic, Dynamic > &r,
const VectorXi &ipvt,
- Scalar tol = ei_sqrt(epsilon<Scalar>()) )
+ Scalar tol = ei_sqrt(NumTraits<Scalar>::epsilon()) )
{
/* Local variables */
int i, j, k, l, ii, jj;
diff --git a/unsupported/Eigen/src/NonLinearOptimization/dogleg.h b/unsupported/Eigen/src/NonLinearOptimization/dogleg.h
index b3c4fbb96..9c1d38dea 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/dogleg.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/dogleg.h
@@ -14,7 +14,7 @@ void ei_dogleg(
Scalar sgnorm;
/* Function Body */
- const Scalar epsmch = epsilon<Scalar>();
+ const Scalar epsmch = NumTraits<Scalar>::epsilon();
const int n = qrfac.cols();
assert(n==qtb.size());
assert(n==x.size());
diff --git a/unsupported/Eigen/src/NonLinearOptimization/fdjac1.h b/unsupported/Eigen/src/NonLinearOptimization/fdjac1.h
index 9564a784b..3dc1e8070 100644
--- a/unsupported/Eigen/src/NonLinearOptimization/fdjac1.h
+++ b/unsupported/Eigen/src/NonLinearOptimization/fdjac1.h
@@ -17,7 +17,7 @@ int ei_fdjac1(
int start, length;
/* Function Body */
- const Scalar epsmch = epsilon<Scalar>();
+ const Scalar epsmch = NumTraits<Scalar>::epsilon();
const int n = x.size();
assert(fvec.size()==n);
Matrix< Scalar, Dynamic, 1 > wa1(n);