LevenbergMarquardt.h
Go to the documentation of this file.
00001 // -*- coding: utf-8
00002 // vim: set fileencoding=utf-8
00003 
00004 // This file is part of Eigen, a lightweight C++ template library
00005 // for linear algebra.
00006 //
00007 // Copyright (C) 2009 Thomas Capricelli <orzel@freehackers.org>
00008 //
00009 // This Source Code Form is subject to the terms of the Mozilla
00010 // Public License v. 2.0. If a copy of the MPL was not distributed
00011 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00012 
00013 #ifndef EIGEN_LEVENBERGMARQUARDT__H
00014 #define EIGEN_LEVENBERGMARQUARDT__H
00015 
00016 namespace Eigen { 
00017 
00018 namespace LevenbergMarquardtSpace {
00019     enum Status {
00020         NotStarted = -2,
00021         Running = -1,
00022         ImproperInputParameters = 0,
00023         RelativeReductionTooSmall = 1,
00024         RelativeErrorTooSmall = 2,
00025         RelativeErrorAndReductionTooSmall = 3,
00026         CosinusTooSmall = 4,
00027         TooManyFunctionEvaluation = 5,
00028         FtolTooSmall = 6,
00029         XtolTooSmall = 7,
00030         GtolTooSmall = 8,
00031         UserAsked = 9
00032     };
00033 }
00034 
00035 
00036 
00045 template<typename FunctorType, typename Scalar=double>
00046 class LevenbergMarquardt
00047 {
00048 public:
00049     LevenbergMarquardt(FunctorType &_functor)
00050         : functor(_functor) { nfev = njev = iter = 0;  fnorm = gnorm = 0.; useExternalScaling=false; }
00051 
00052     typedef DenseIndex Index;
00053 
00054     struct Parameters {
00055         Parameters()
00056             : factor(Scalar(100.))
00057             , maxfev(400)
00058             , ftol(internal::sqrt(NumTraits<Scalar>::epsilon()))
00059             , xtol(internal::sqrt(NumTraits<Scalar>::epsilon()))
00060             , gtol(Scalar(0.))
00061             , epsfcn(Scalar(0.)) {}
00062         Scalar factor;
00063         Index maxfev;   // maximum number of function evaluation
00064         Scalar ftol;
00065         Scalar xtol;
00066         Scalar gtol;
00067         Scalar epsfcn;
00068     };
00069 
00070     typedef Matrix< Scalar, Dynamic, 1 > FVectorType;
00071     typedef Matrix< Scalar, Dynamic, Dynamic > JacobianType;
00072 
00073     LevenbergMarquardtSpace::Status lmder1(
00074             FVectorType &x,
00075             const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00076             );
00077 
00078     LevenbergMarquardtSpace::Status minimize(FVectorType &x);
00079     LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x);
00080     LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x);
00081 
00082     static LevenbergMarquardtSpace::Status lmdif1(
00083             FunctorType &functor,
00084             FVectorType &x,
00085             Index *nfev,
00086             const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00087             );
00088 
00089     LevenbergMarquardtSpace::Status lmstr1(
00090             FVectorType  &x,
00091             const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00092             );
00093 
00094     LevenbergMarquardtSpace::Status minimizeOptimumStorage(FVectorType  &x);
00095     LevenbergMarquardtSpace::Status minimizeOptimumStorageInit(FVectorType  &x);
00096     LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep(FVectorType  &x);
00097 
00098     void resetParameters(void) { parameters = Parameters(); }
00099 
00100     Parameters parameters;
00101     FVectorType  fvec, qtf, diag;
00102     JacobianType fjac;
00103     PermutationMatrix<Dynamic,Dynamic> permutation;
00104     Index nfev;
00105     Index njev;
00106     Index iter;
00107     Scalar fnorm, gnorm;
00108     bool useExternalScaling; 
00109 
00110     Scalar lm_param(void) { return par; }
00111 private:
00112     FunctorType &functor;
00113     Index n;
00114     Index m;
00115     FVectorType wa1, wa2, wa3, wa4;
00116 
00117     Scalar par, sum;
00118     Scalar temp, temp1, temp2;
00119     Scalar delta;
00120     Scalar ratio;
00121     Scalar pnorm, xnorm, fnorm1, actred, dirder, prered;
00122 
00123     LevenbergMarquardt& operator=(const LevenbergMarquardt&);
00124 };
00125 
00126 template<typename FunctorType, typename Scalar>
00127 LevenbergMarquardtSpace::Status
00128 LevenbergMarquardt<FunctorType,Scalar>::lmder1(
00129         FVectorType  &x,
00130         const Scalar tol
00131         )
00132 {
00133     n = x.size();
00134     m = functor.values();
00135 
00136     /* check the input parameters for errors. */
00137     if (n <= 0 || m < n || tol < 0.)
00138         return LevenbergMarquardtSpace::ImproperInputParameters;
00139 
00140     resetParameters();
00141     parameters.ftol = tol;
00142     parameters.xtol = tol;
00143     parameters.maxfev = 100*(n+1);
00144 
00145     return minimize(x);
00146 }
00147 
00148 
00149 template<typename FunctorType, typename Scalar>
00150 LevenbergMarquardtSpace::Status
00151 LevenbergMarquardt<FunctorType,Scalar>::minimize(FVectorType  &x)
00152 {
00153     LevenbergMarquardtSpace::Status status = minimizeInit(x);
00154     if (status==LevenbergMarquardtSpace::ImproperInputParameters)
00155         return status;
00156     do {
00157         status = minimizeOneStep(x);
00158     } while (status==LevenbergMarquardtSpace::Running);
00159     return status;
00160 }
00161 
00162 template<typename FunctorType, typename Scalar>
00163 LevenbergMarquardtSpace::Status
00164 LevenbergMarquardt<FunctorType,Scalar>::minimizeInit(FVectorType  &x)
00165 {
00166     n = x.size();
00167     m = functor.values();
00168 
00169     wa1.resize(n); wa2.resize(n); wa3.resize(n);
00170     wa4.resize(m);
00171     fvec.resize(m);
00172     fjac.resize(m, n);
00173     if (!useExternalScaling)
00174         diag.resize(n);
00175     assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00176     qtf.resize(n);
00177 
00178     /* Function Body */
00179     nfev = 0;
00180     njev = 0;
00181 
00182     /*     check the input parameters for errors. */
00183     if (n <= 0 || m < n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0.)
00184         return LevenbergMarquardtSpace::ImproperInputParameters;
00185 
00186     if (useExternalScaling)
00187         for (Index j = 0; j < n; ++j)
00188             if (diag[j] <= 0.)
00189                 return LevenbergMarquardtSpace::ImproperInputParameters;
00190 
00191     /*     evaluate the function at the starting point */
00192     /*     and calculate its norm. */
00193     nfev = 1;
00194     if ( functor(x, fvec) < 0)
00195         return LevenbergMarquardtSpace::UserAsked;
00196     fnorm = fvec.stableNorm();
00197 
00198     /*     initialize levenberg-marquardt parameter and iteration counter. */
00199     par = 0.;
00200     iter = 1;
00201 
00202     return LevenbergMarquardtSpace::NotStarted;
00203 }
00204 
00205 template<typename FunctorType, typename Scalar>
00206 LevenbergMarquardtSpace::Status
00207 LevenbergMarquardt<FunctorType,Scalar>::minimizeOneStep(FVectorType  &x)
00208 {
00209     assert(x.size()==n); // check the caller is not cheating us
00210 
00211     /* calculate the jacobian matrix. */
00212     Index df_ret = functor.df(x, fjac);
00213     if (df_ret<0)
00214         return LevenbergMarquardtSpace::UserAsked;
00215     if (df_ret>0)
00216         // numerical diff, we evaluated the function df_ret times
00217         nfev += df_ret;
00218     else njev++;
00219 
00220     /* compute the qr factorization of the jacobian. */
00221     wa2 = fjac.colwise().blueNorm();
00222     ColPivHouseholderQR<JacobianType> qrfac(fjac);
00223     fjac = qrfac.matrixQR();
00224     permutation = qrfac.colsPermutation();
00225 
00226     /* on the first iteration and if external scaling is not used, scale according */
00227     /* to the norms of the columns of the initial jacobian. */
00228     if (iter == 1) {
00229         if (!useExternalScaling)
00230             for (Index j = 0; j < n; ++j)
00231                 diag[j] = (wa2[j]==0.)? 1. : wa2[j];
00232 
00233         /* on the first iteration, calculate the norm of the scaled x */
00234         /* and initialize the step bound delta. */
00235         xnorm = diag.cwiseProduct(x).stableNorm();
00236         delta = parameters.factor * xnorm;
00237         if (delta == 0.)
00238             delta = parameters.factor;
00239     }
00240 
00241     /* form (q transpose)*fvec and store the first n components in */
00242     /* qtf. */
00243     wa4 = fvec;
00244     wa4.applyOnTheLeft(qrfac.householderQ().adjoint());
00245     qtf = wa4.head(n);
00246 
00247     /* compute the norm of the scaled gradient. */
00248     gnorm = 0.;
00249     if (fnorm != 0.)
00250         for (Index j = 0; j < n; ++j)
00251             if (wa2[permutation.indices()[j]] != 0.)
00252                 gnorm = (std::max)(gnorm, internal::abs( fjac.col(j).head(j+1).dot(qtf.head(j+1)/fnorm) / wa2[permutation.indices()[j]]));
00253 
00254     /* test for convergence of the gradient norm. */
00255     if (gnorm <= parameters.gtol)
00256         return LevenbergMarquardtSpace::CosinusTooSmall;
00257 
00258     /* rescale if necessary. */
00259     if (!useExternalScaling)
00260         diag = diag.cwiseMax(wa2);
00261 
00262     do {
00263 
00264         /* determine the levenberg-marquardt parameter. */
00265         internal::lmpar2<Scalar>(qrfac, diag, qtf, delta, par, wa1);
00266 
00267         /* store the direction p and x + p. calculate the norm of p. */
00268         wa1 = -wa1;
00269         wa2 = x + wa1;
00270         pnorm = diag.cwiseProduct(wa1).stableNorm();
00271 
00272         /* on the first iteration, adjust the initial step bound. */
00273         if (iter == 1)
00274             delta = (std::min)(delta,pnorm);
00275 
00276         /* evaluate the function at x + p and calculate its norm. */
00277         if ( functor(wa2, wa4) < 0)
00278             return LevenbergMarquardtSpace::UserAsked;
00279         ++nfev;
00280         fnorm1 = wa4.stableNorm();
00281 
00282         /* compute the scaled actual reduction. */
00283         actred = -1.;
00284         if (Scalar(.1) * fnorm1 < fnorm)
00285             actred = 1. - internal::abs2(fnorm1 / fnorm);
00286 
00287         /* compute the scaled predicted reduction and */
00288         /* the scaled directional derivative. */
00289         wa3 = fjac.template triangularView<Upper>() * (qrfac.colsPermutation().inverse() *wa1);
00290         temp1 = internal::abs2(wa3.stableNorm() / fnorm);
00291         temp2 = internal::abs2(internal::sqrt(par) * pnorm / fnorm);
00292         prered = temp1 + temp2 / Scalar(.5);
00293         dirder = -(temp1 + temp2);
00294 
00295         /* compute the ratio of the actual to the predicted */
00296         /* reduction. */
00297         ratio = 0.;
00298         if (prered != 0.)
00299             ratio = actred / prered;
00300 
00301         /* update the step bound. */
00302         if (ratio <= Scalar(.25)) {
00303             if (actred >= 0.)
00304                 temp = Scalar(.5);
00305             if (actred < 0.)
00306                 temp = Scalar(.5) * dirder / (dirder + Scalar(.5) * actred);
00307             if (Scalar(.1) * fnorm1 >= fnorm || temp < Scalar(.1))
00308                 temp = Scalar(.1);
00309             /* Computing MIN */
00310             delta = temp * (std::min)(delta, pnorm / Scalar(.1));
00311             par /= temp;
00312         } else if (!(par != 0. && ratio < Scalar(.75))) {
00313             delta = pnorm / Scalar(.5);
00314             par = Scalar(.5) * par;
00315         }
00316 
00317         /* test for successful iteration. */
00318         if (ratio >= Scalar(1e-4)) {
00319             /* successful iteration. update x, fvec, and their norms. */
00320             x = wa2;
00321             wa2 = diag.cwiseProduct(x);
00322             fvec = wa4;
00323             xnorm = wa2.stableNorm();
00324             fnorm = fnorm1;
00325             ++iter;
00326         }
00327 
00328         /* tests for convergence. */
00329         if (internal::abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1. && delta <= parameters.xtol * xnorm)
00330             return LevenbergMarquardtSpace::RelativeErrorAndReductionTooSmall;
00331         if (internal::abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1.)
00332             return LevenbergMarquardtSpace::RelativeReductionTooSmall;
00333         if (delta <= parameters.xtol * xnorm)
00334             return LevenbergMarquardtSpace::RelativeErrorTooSmall;
00335 
00336         /* tests for termination and stringent tolerances. */
00337         if (nfev >= parameters.maxfev)
00338             return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
00339         if (internal::abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
00340             return LevenbergMarquardtSpace::FtolTooSmall;
00341         if (delta <= NumTraits<Scalar>::epsilon() * xnorm)
00342             return LevenbergMarquardtSpace::XtolTooSmall;
00343         if (gnorm <= NumTraits<Scalar>::epsilon())
00344             return LevenbergMarquardtSpace::GtolTooSmall;
00345 
00346     } while (ratio < Scalar(1e-4));
00347 
00348     return LevenbergMarquardtSpace::Running;
00349 }
00350 
00351 template<typename FunctorType, typename Scalar>
00352 LevenbergMarquardtSpace::Status
00353 LevenbergMarquardt<FunctorType,Scalar>::lmstr1(
00354         FVectorType  &x,
00355         const Scalar tol
00356         )
00357 {
00358     n = x.size();
00359     m = functor.values();
00360 
00361     /* check the input parameters for errors. */
00362     if (n <= 0 || m < n || tol < 0.)
00363         return LevenbergMarquardtSpace::ImproperInputParameters;
00364 
00365     resetParameters();
00366     parameters.ftol = tol;
00367     parameters.xtol = tol;
00368     parameters.maxfev = 100*(n+1);
00369 
00370     return minimizeOptimumStorage(x);
00371 }
00372 
00373 template<typename FunctorType, typename Scalar>
00374 LevenbergMarquardtSpace::Status
00375 LevenbergMarquardt<FunctorType,Scalar>::minimizeOptimumStorageInit(FVectorType  &x)
00376 {
00377     n = x.size();
00378     m = functor.values();
00379 
00380     wa1.resize(n); wa2.resize(n); wa3.resize(n);
00381     wa4.resize(m);
00382     fvec.resize(m);
00383     // Only R is stored in fjac. Q is only used to compute 'qtf', which is
00384     // Q.transpose()*rhs. qtf will be updated using givens rotation,
00385     // instead of storing them in Q.
00386     // The purpose it to only use a nxn matrix, instead of mxn here, so
00387     // that we can handle cases where m>>n :
00388     fjac.resize(n, n);
00389     if (!useExternalScaling)
00390         diag.resize(n);
00391     assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00392     qtf.resize(n);
00393 
00394     /* Function Body */
00395     nfev = 0;
00396     njev = 0;
00397 
00398     /*     check the input parameters for errors. */
00399     if (n <= 0 || m < n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0.)
00400         return LevenbergMarquardtSpace::ImproperInputParameters;
00401 
00402     if (useExternalScaling)
00403         for (Index j = 0; j < n; ++j)
00404             if (diag[j] <= 0.)
00405                 return LevenbergMarquardtSpace::ImproperInputParameters;
00406 
00407     /*     evaluate the function at the starting point */
00408     /*     and calculate its norm. */
00409     nfev = 1;
00410     if ( functor(x, fvec) < 0)
00411         return LevenbergMarquardtSpace::UserAsked;
00412     fnorm = fvec.stableNorm();
00413 
00414     /*     initialize levenberg-marquardt parameter and iteration counter. */
00415     par = 0.;
00416     iter = 1;
00417 
00418     return LevenbergMarquardtSpace::NotStarted;
00419 }
00420 
00421 
00422 template<typename FunctorType, typename Scalar>
00423 LevenbergMarquardtSpace::Status
00424 LevenbergMarquardt<FunctorType,Scalar>::minimizeOptimumStorageOneStep(FVectorType  &x)
00425 {
00426     assert(x.size()==n); // check the caller is not cheating us
00427 
00428     Index i, j;
00429     bool sing;
00430 
00431     /* compute the qr factorization of the jacobian matrix */
00432     /* calculated one row at a time, while simultaneously */
00433     /* forming (q transpose)*fvec and storing the first */
00434     /* n components in qtf. */
00435     qtf.fill(0.);
00436     fjac.fill(0.);
00437     Index rownb = 2;
00438     for (i = 0; i < m; ++i) {
00439         if (functor.df(x, wa3, rownb) < 0) return LevenbergMarquardtSpace::UserAsked;
00440         internal::rwupdt<Scalar>(fjac, wa3, qtf, fvec[i]);
00441         ++rownb;
00442     }
00443     ++njev;
00444 
00445     /* if the jacobian is rank deficient, call qrfac to */
00446     /* reorder its columns and update the components of qtf. */
00447     sing = false;
00448     for (j = 0; j < n; ++j) {
00449         if (fjac(j,j) == 0.)
00450             sing = true;
00451         wa2[j] = fjac.col(j).head(j).stableNorm();
00452     }
00453     permutation.setIdentity(n);
00454     if (sing) {
00455         wa2 = fjac.colwise().blueNorm();
00456         // TODO We have no unit test covering this code path, do not modify
00457         // until it is carefully tested
00458         ColPivHouseholderQR<JacobianType> qrfac(fjac);
00459         fjac = qrfac.matrixQR();
00460         wa1 = fjac.diagonal();
00461         fjac.diagonal() = qrfac.hCoeffs();
00462         permutation = qrfac.colsPermutation();
00463         // TODO : avoid this:
00464         for(Index ii=0; ii< fjac.cols(); ii++) fjac.col(ii).segment(ii+1, fjac.rows()-ii-1) *= fjac(ii,ii); // rescale vectors
00465 
00466         for (j = 0; j < n; ++j) {
00467             if (fjac(j,j) != 0.) {
00468                 sum = 0.;
00469                 for (i = j; i < n; ++i)
00470                     sum += fjac(i,j) * qtf[i];
00471                 temp = -sum / fjac(j,j);
00472                 for (i = j; i < n; ++i)
00473                     qtf[i] += fjac(i,j) * temp;
00474             }
00475             fjac(j,j) = wa1[j];
00476         }
00477     }
00478 
00479     /* on the first iteration and if external scaling is not used, scale according */
00480     /* to the norms of the columns of the initial jacobian. */
00481     if (iter == 1) {
00482         if (!useExternalScaling)
00483             for (j = 0; j < n; ++j)
00484                 diag[j] = (wa2[j]==0.)? 1. : wa2[j];
00485 
00486         /* on the first iteration, calculate the norm of the scaled x */
00487         /* and initialize the step bound delta. */
00488         xnorm = diag.cwiseProduct(x).stableNorm();
00489         delta = parameters.factor * xnorm;
00490         if (delta == 0.)
00491             delta = parameters.factor;
00492     }
00493 
00494     /* compute the norm of the scaled gradient. */
00495     gnorm = 0.;
00496     if (fnorm != 0.)
00497         for (j = 0; j < n; ++j)
00498             if (wa2[permutation.indices()[j]] != 0.)
00499                 gnorm = (std::max)(gnorm, internal::abs( fjac.col(j).head(j+1).dot(qtf.head(j+1)/fnorm) / wa2[permutation.indices()[j]]));
00500 
00501     /* test for convergence of the gradient norm. */
00502     if (gnorm <= parameters.gtol)
00503         return LevenbergMarquardtSpace::CosinusTooSmall;
00504 
00505     /* rescale if necessary. */
00506     if (!useExternalScaling)
00507         diag = diag.cwiseMax(wa2);
00508 
00509     do {
00510 
00511         /* determine the levenberg-marquardt parameter. */
00512         internal::lmpar<Scalar>(fjac, permutation.indices(), diag, qtf, delta, par, wa1);
00513 
00514         /* store the direction p and x + p. calculate the norm of p. */
00515         wa1 = -wa1;
00516         wa2 = x + wa1;
00517         pnorm = diag.cwiseProduct(wa1).stableNorm();
00518 
00519         /* on the first iteration, adjust the initial step bound. */
00520         if (iter == 1)
00521             delta = (std::min)(delta,pnorm);
00522 
00523         /* evaluate the function at x + p and calculate its norm. */
00524         if ( functor(wa2, wa4) < 0)
00525             return LevenbergMarquardtSpace::UserAsked;
00526         ++nfev;
00527         fnorm1 = wa4.stableNorm();
00528 
00529         /* compute the scaled actual reduction. */
00530         actred = -1.;
00531         if (Scalar(.1) * fnorm1 < fnorm)
00532             actred = 1. - internal::abs2(fnorm1 / fnorm);
00533 
00534         /* compute the scaled predicted reduction and */
00535         /* the scaled directional derivative. */
00536         wa3 = fjac.topLeftCorner(n,n).template triangularView<Upper>() * (permutation.inverse() * wa1);
00537         temp1 = internal::abs2(wa3.stableNorm() / fnorm);
00538         temp2 = internal::abs2(internal::sqrt(par) * pnorm / fnorm);
00539         prered = temp1 + temp2 / Scalar(.5);
00540         dirder = -(temp1 + temp2);
00541 
00542         /* compute the ratio of the actual to the predicted */
00543         /* reduction. */
00544         ratio = 0.;
00545         if (prered != 0.)
00546             ratio = actred / prered;
00547 
00548         /* update the step bound. */
00549         if (ratio <= Scalar(.25)) {
00550             if (actred >= 0.)
00551                 temp = Scalar(.5);
00552             if (actred < 0.)
00553                 temp = Scalar(.5) * dirder / (dirder + Scalar(.5) * actred);
00554             if (Scalar(.1) * fnorm1 >= fnorm || temp < Scalar(.1))
00555                 temp = Scalar(.1);
00556             /* Computing MIN */
00557             delta = temp * (std::min)(delta, pnorm / Scalar(.1));
00558             par /= temp;
00559         } else if (!(par != 0. && ratio < Scalar(.75))) {
00560             delta = pnorm / Scalar(.5);
00561             par = Scalar(.5) * par;
00562         }
00563 
00564         /* test for successful iteration. */
00565         if (ratio >= Scalar(1e-4)) {
00566             /* successful iteration. update x, fvec, and their norms. */
00567             x = wa2;
00568             wa2 = diag.cwiseProduct(x);
00569             fvec = wa4;
00570             xnorm = wa2.stableNorm();
00571             fnorm = fnorm1;
00572             ++iter;
00573         }
00574 
00575         /* tests for convergence. */
00576         if (internal::abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1. && delta <= parameters.xtol * xnorm)
00577             return LevenbergMarquardtSpace::RelativeErrorAndReductionTooSmall;
00578         if (internal::abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1.)
00579             return LevenbergMarquardtSpace::RelativeReductionTooSmall;
00580         if (delta <= parameters.xtol * xnorm)
00581             return LevenbergMarquardtSpace::RelativeErrorTooSmall;
00582 
00583         /* tests for termination and stringent tolerances. */
00584         if (nfev >= parameters.maxfev)
00585             return LevenbergMarquardtSpace::TooManyFunctionEvaluation;
00586         if (internal::abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
00587             return LevenbergMarquardtSpace::FtolTooSmall;
00588         if (delta <= NumTraits<Scalar>::epsilon() * xnorm)
00589             return LevenbergMarquardtSpace::XtolTooSmall;
00590         if (gnorm <= NumTraits<Scalar>::epsilon())
00591             return LevenbergMarquardtSpace::GtolTooSmall;
00592 
00593     } while (ratio < Scalar(1e-4));
00594 
00595     return LevenbergMarquardtSpace::Running;
00596 }
00597 
00598 template<typename FunctorType, typename Scalar>
00599 LevenbergMarquardtSpace::Status
00600 LevenbergMarquardt<FunctorType,Scalar>::minimizeOptimumStorage(FVectorType  &x)
00601 {
00602     LevenbergMarquardtSpace::Status status = minimizeOptimumStorageInit(x);
00603     if (status==LevenbergMarquardtSpace::ImproperInputParameters)
00604         return status;
00605     do {
00606         status = minimizeOptimumStorageOneStep(x);
00607     } while (status==LevenbergMarquardtSpace::Running);
00608     return status;
00609 }
00610 
00611 template<typename FunctorType, typename Scalar>
00612 LevenbergMarquardtSpace::Status
00613 LevenbergMarquardt<FunctorType,Scalar>::lmdif1(
00614         FunctorType &functor,
00615         FVectorType  &x,
00616         Index *nfev,
00617         const Scalar tol
00618         )
00619 {
00620     Index n = x.size();
00621     Index m = functor.values();
00622 
00623     /* check the input parameters for errors. */
00624     if (n <= 0 || m < n || tol < 0.)
00625         return LevenbergMarquardtSpace::ImproperInputParameters;
00626 
00627     NumericalDiff<FunctorType> numDiff(functor);
00628     // embedded LevenbergMarquardt
00629     LevenbergMarquardt<NumericalDiff<FunctorType>, Scalar > lm(numDiff);
00630     lm.parameters.ftol = tol;
00631     lm.parameters.xtol = tol;
00632     lm.parameters.maxfev = 200*(n+1);
00633 
00634     LevenbergMarquardtSpace::Status info = LevenbergMarquardtSpace::Status(lm.minimize(x));
00635     if (nfev)
00636         * nfev = lm.nfev;
00637     return info;
00638 }
00639 
00640 } // end namespace Eigen
00641 
00642 #endif // EIGEN_LEVENBERGMARQUARDT__H
00643 
00644 //vim: ai ts=4 sts=4 et sw=4


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:11:07