HybridNonLinearSolver.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 // Eigen is free software; you can redistribute it and/or
00010 // modify it under the terms of the GNU Lesser General Public
00011 // License as published by the Free Software Foundation; either
00012 // version 3 of the License, or (at your option) any later version.
00013 //
00014 // Alternatively, you can redistribute it and/or
00015 // modify it under the terms of the GNU General Public License as
00016 // published by the Free Software Foundation; either version 2 of
00017 // the License, or (at your option) any later version.
00018 //
00019 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00020 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00021 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00022 // GNU General Public License for more details.
00023 //
00024 // You should have received a copy of the GNU Lesser General Public
00025 // License and a copy of the GNU General Public License along with
00026 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00027 
00028 #ifndef EIGEN_HYBRIDNONLINEARSOLVER_H
00029 #define EIGEN_HYBRIDNONLINEARSOLVER_H
00030 
00031 namespace HybridNonLinearSolverSpace { 
00032     enum Status {
00033         Running = -1,
00034         ImproperInputParameters = 0,
00035         RelativeErrorTooSmall = 1,
00036         TooManyFunctionEvaluation = 2,
00037         TolTooSmall = 3,
00038         NotMakingProgressJacobian = 4,
00039         NotMakingProgressIterations = 5,
00040         UserAksed = 6
00041     };
00042 }
00043 
00055 template<typename FunctorType, typename Scalar=double>
00056 class HybridNonLinearSolver
00057 {
00058 public:
00059     typedef DenseIndex Index;
00060 
00061     HybridNonLinearSolver(FunctorType &_functor)
00062         : functor(_functor) { nfev=njev=iter = 0;  fnorm= 0.; useExternalScaling=false;}
00063 
00064     struct Parameters {
00065         Parameters()
00066             : factor(Scalar(100.))
00067             , maxfev(1000)
00068             , xtol(internal::sqrt(NumTraits<Scalar>::epsilon()))
00069             , nb_of_subdiagonals(-1)
00070             , nb_of_superdiagonals(-1)
00071             , epsfcn(Scalar(0.)) {}
00072         Scalar factor;
00073         Index maxfev;   // maximum number of function evaluation
00074         Scalar xtol;
00075         Index nb_of_subdiagonals;
00076         Index nb_of_superdiagonals;
00077         Scalar epsfcn;
00078     };
00079     typedef Matrix< Scalar, Dynamic, 1 > FVectorType;
00080     typedef Matrix< Scalar, Dynamic, Dynamic > JacobianType;
00081     /* TODO: if eigen provides a triangular storage, use it here */
00082     typedef Matrix< Scalar, Dynamic, Dynamic > UpperTriangularType;
00083 
00084     HybridNonLinearSolverSpace::Status hybrj1(
00085             FVectorType  &x,
00086             const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00087             );
00088 
00089     HybridNonLinearSolverSpace::Status solveInit(FVectorType  &x);
00090     HybridNonLinearSolverSpace::Status solveOneStep(FVectorType  &x);
00091     HybridNonLinearSolverSpace::Status solve(FVectorType  &x);
00092 
00093     HybridNonLinearSolverSpace::Status hybrd1(
00094             FVectorType  &x,
00095             const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00096             );
00097 
00098     HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType  &x);
00099     HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType  &x);
00100     HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType  &x);
00101 
00102     void resetParameters(void) { parameters = Parameters(); }
00103     Parameters parameters;
00104     FVectorType  fvec, qtf, diag;
00105     JacobianType fjac;
00106     UpperTriangularType R;
00107     Index nfev;
00108     Index njev;
00109     Index iter;
00110     Scalar fnorm;
00111     bool useExternalScaling; 
00112 private:
00113     FunctorType &functor;
00114     Index n;
00115     Scalar sum;
00116     bool sing;
00117     Scalar temp;
00118     Scalar delta;
00119     bool jeval;
00120     Index ncsuc;
00121     Scalar ratio;
00122     Scalar pnorm, xnorm, fnorm1;
00123     Index nslow1, nslow2;
00124     Index ncfail;
00125     Scalar actred, prered;
00126     FVectorType wa1, wa2, wa3, wa4;
00127 
00128     HybridNonLinearSolver& operator=(const HybridNonLinearSolver&);
00129 };
00130 
00131 
00132 
00133 template<typename FunctorType, typename Scalar>
00134 HybridNonLinearSolverSpace::Status
00135 HybridNonLinearSolver<FunctorType,Scalar>::hybrj1(
00136         FVectorType  &x,
00137         const Scalar tol
00138         )
00139 {
00140     n = x.size();
00141 
00142     /* check the input parameters for errors. */
00143     if (n <= 0 || tol < 0.)
00144         return HybridNonLinearSolverSpace::ImproperInputParameters;
00145 
00146     resetParameters();
00147     parameters.maxfev = 100*(n+1);
00148     parameters.xtol = tol;
00149     diag.setConstant(n, 1.);
00150     useExternalScaling = true;
00151     return solve(x);
00152 }
00153 
00154 template<typename FunctorType, typename Scalar>
00155 HybridNonLinearSolverSpace::Status
00156 HybridNonLinearSolver<FunctorType,Scalar>::solveInit(FVectorType  &x)
00157 {
00158     n = x.size();
00159 
00160     wa1.resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
00161     fvec.resize(n);
00162     qtf.resize(n);
00163     fjac.resize(n, n);
00164     if (!useExternalScaling)
00165         diag.resize(n);
00166     assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00167 
00168     /* Function Body */
00169     nfev = 0;
00170     njev = 0;
00171 
00172     /*     check the input parameters for errors. */
00173     if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0. )
00174         return HybridNonLinearSolverSpace::ImproperInputParameters;
00175     if (useExternalScaling)
00176         for (Index j = 0; j < n; ++j)
00177             if (diag[j] <= 0.)
00178                 return HybridNonLinearSolverSpace::ImproperInputParameters;
00179 
00180     /*     evaluate the function at the starting point */
00181     /*     and calculate its norm. */
00182     nfev = 1;
00183     if ( functor(x, fvec) < 0)
00184         return HybridNonLinearSolverSpace::UserAksed;
00185     fnorm = fvec.stableNorm();
00186 
00187     /*     initialize iteration counter and monitors. */
00188     iter = 1;
00189     ncsuc = 0;
00190     ncfail = 0;
00191     nslow1 = 0;
00192     nslow2 = 0;
00193 
00194     return HybridNonLinearSolverSpace::Running;
00195 }
00196 
00197 template<typename FunctorType, typename Scalar>
00198 HybridNonLinearSolverSpace::Status
00199 HybridNonLinearSolver<FunctorType,Scalar>::solveOneStep(FVectorType  &x)
00200 {
00201     assert(x.size()==n); // check the caller is not cheating us
00202 
00203     Index j;
00204     std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00205 
00206     jeval = true;
00207 
00208     /* calculate the jacobian matrix. */
00209     if ( functor.df(x, fjac) < 0)
00210         return HybridNonLinearSolverSpace::UserAksed;
00211     ++njev;
00212 
00213     wa2 = fjac.colwise().blueNorm();
00214 
00215     /* on the first iteration and if external scaling is not used, scale according */
00216     /* to the norms of the columns of the initial jacobian. */
00217     if (iter == 1) {
00218         if (!useExternalScaling)
00219             for (j = 0; j < n; ++j)
00220                 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
00221 
00222         /* on the first iteration, calculate the norm of the scaled x */
00223         /* and initialize the step bound delta. */
00224         xnorm = diag.cwiseProduct(x).stableNorm();
00225         delta = parameters.factor * xnorm;
00226         if (delta == 0.)
00227             delta = parameters.factor;
00228     }
00229 
00230     /* compute the qr factorization of the jacobian. */
00231     wa2 = fjac.colwise().blueNorm();
00232     HouseholderQR<JacobianType> qrfac(fjac); // no pivoting:
00233 
00234     /* copy the triangular factor of the qr factorization into r. */
00235     R = qrfac.matrixQR();
00236 
00237     /* accumulate the orthogonal factor in fjac. */
00238     fjac = qrfac.householderQ();
00239 
00240     /* form (q transpose)*fvec and store in qtf. */
00241     qtf = fjac.transpose() * fvec;
00242 
00243     /* rescale if necessary. */
00244     if (!useExternalScaling)
00245         diag = diag.cwiseMax(wa2);
00246 
00247     while (true) {
00248         /* determine the direction p. */
00249         internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00250 
00251         /* store the direction p and x + p. calculate the norm of p. */
00252         wa1 = -wa1;
00253         wa2 = x + wa1;
00254         pnorm = diag.cwiseProduct(wa1).stableNorm();
00255 
00256         /* on the first iteration, adjust the initial step bound. */
00257         if (iter == 1)
00258             delta = std::min(delta,pnorm);
00259 
00260         /* evaluate the function at x + p and calculate its norm. */
00261         if ( functor(wa2, wa4) < 0)
00262             return HybridNonLinearSolverSpace::UserAksed;
00263         ++nfev;
00264         fnorm1 = wa4.stableNorm();
00265 
00266         /* compute the scaled actual reduction. */
00267         actred = -1.;
00268         if (fnorm1 < fnorm) /* Computing 2nd power */
00269             actred = 1. - internal::abs2(fnorm1 / fnorm);
00270 
00271         /* compute the scaled predicted reduction. */
00272         wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00273         temp = wa3.stableNorm();
00274         prered = 0.;
00275         if (temp < fnorm) /* Computing 2nd power */
00276             prered = 1. - internal::abs2(temp / fnorm);
00277 
00278         /* compute the ratio of the actual to the predicted reduction. */
00279         ratio = 0.;
00280         if (prered > 0.)
00281             ratio = actred / prered;
00282 
00283         /* update the step bound. */
00284         if (ratio < Scalar(.1)) {
00285             ncsuc = 0;
00286             ++ncfail;
00287             delta = Scalar(.5) * delta;
00288         } else {
00289             ncfail = 0;
00290             ++ncsuc;
00291             if (ratio >= Scalar(.5) || ncsuc > 1)
00292                 delta = std::max(delta, pnorm / Scalar(.5));
00293             if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00294                 delta = pnorm / Scalar(.5);
00295             }
00296         }
00297 
00298         /* test for successful iteration. */
00299         if (ratio >= Scalar(1e-4)) {
00300             /* successful iteration. update x, fvec, and their norms. */
00301             x = wa2;
00302             wa2 = diag.cwiseProduct(x);
00303             fvec = wa4;
00304             xnorm = wa2.stableNorm();
00305             fnorm = fnorm1;
00306             ++iter;
00307         }
00308 
00309         /* determine the progress of the iteration. */
00310         ++nslow1;
00311         if (actred >= Scalar(.001))
00312             nslow1 = 0;
00313         if (jeval)
00314             ++nslow2;
00315         if (actred >= Scalar(.1))
00316             nslow2 = 0;
00317 
00318         /* test for convergence. */
00319         if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00320             return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00321 
00322         /* tests for termination and stringent tolerances. */
00323         if (nfev >= parameters.maxfev)
00324             return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00325         if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00326             return HybridNonLinearSolverSpace::TolTooSmall;
00327         if (nslow2 == 5)
00328             return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00329         if (nslow1 == 10)
00330             return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00331 
00332         /* criterion for recalculating jacobian. */
00333         if (ncfail == 2)
00334             break; // leave inner loop and go for the next outer loop iteration
00335 
00336         /* calculate the rank one modification to the jacobian */
00337         /* and update qtf if necessary. */
00338         wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00339         wa2 = fjac.transpose() * wa4;
00340         if (ratio >= Scalar(1e-4))
00341             qtf = wa2;
00342         wa2 = (wa2-wa3)/pnorm;
00343 
00344         /* compute the qr factorization of the updated jacobian. */
00345         internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00346         internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00347         internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00348 
00349         jeval = false;
00350     }
00351     return HybridNonLinearSolverSpace::Running;
00352 }
00353 
00354 template<typename FunctorType, typename Scalar>
00355 HybridNonLinearSolverSpace::Status
00356 HybridNonLinearSolver<FunctorType,Scalar>::solve(FVectorType  &x)
00357 {
00358     HybridNonLinearSolverSpace::Status status = solveInit(x);
00359     if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00360         return status;
00361     while (status==HybridNonLinearSolverSpace::Running)
00362         status = solveOneStep(x);
00363     return status;
00364 }
00365 
00366 
00367 
00368 template<typename FunctorType, typename Scalar>
00369 HybridNonLinearSolverSpace::Status
00370 HybridNonLinearSolver<FunctorType,Scalar>::hybrd1(
00371         FVectorType  &x,
00372         const Scalar tol
00373         )
00374 {
00375     n = x.size();
00376 
00377     /* check the input parameters for errors. */
00378     if (n <= 0 || tol < 0.)
00379         return HybridNonLinearSolverSpace::ImproperInputParameters;
00380 
00381     resetParameters();
00382     parameters.maxfev = 200*(n+1);
00383     parameters.xtol = tol;
00384 
00385     diag.setConstant(n, 1.);
00386     useExternalScaling = true;
00387     return solveNumericalDiff(x);
00388 }
00389 
00390 template<typename FunctorType, typename Scalar>
00391 HybridNonLinearSolverSpace::Status
00392 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffInit(FVectorType  &x)
00393 {
00394     n = x.size();
00395 
00396     if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00397     if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00398 
00399     wa1.resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
00400     qtf.resize(n);
00401     fjac.resize(n, n);
00402     fvec.resize(n);
00403     if (!useExternalScaling)
00404         diag.resize(n);
00405     assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00406 
00407     /* Function Body */
00408     nfev = 0;
00409     njev = 0;
00410 
00411     /*     check the input parameters for errors. */
00412     if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.nb_of_subdiagonals< 0 || parameters.nb_of_superdiagonals< 0 || parameters.factor <= 0. )
00413         return HybridNonLinearSolverSpace::ImproperInputParameters;
00414     if (useExternalScaling)
00415         for (Index j = 0; j < n; ++j)
00416             if (diag[j] <= 0.)
00417                 return HybridNonLinearSolverSpace::ImproperInputParameters;
00418 
00419     /*     evaluate the function at the starting point */
00420     /*     and calculate its norm. */
00421     nfev = 1;
00422     if ( functor(x, fvec) < 0)
00423         return HybridNonLinearSolverSpace::UserAksed;
00424     fnorm = fvec.stableNorm();
00425 
00426     /*     initialize iteration counter and monitors. */
00427     iter = 1;
00428     ncsuc = 0;
00429     ncfail = 0;
00430     nslow1 = 0;
00431     nslow2 = 0;
00432 
00433     return HybridNonLinearSolverSpace::Running;
00434 }
00435 
00436 template<typename FunctorType, typename Scalar>
00437 HybridNonLinearSolverSpace::Status
00438 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(FVectorType  &x)
00439 {
00440     assert(x.size()==n); // check the caller is not cheating us
00441 
00442     Index j;
00443     std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00444 
00445     jeval = true;
00446     if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00447     if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00448 
00449     /* calculate the jacobian matrix. */
00450     if (internal::fdjac1(functor, x, fvec, fjac, parameters.nb_of_subdiagonals, parameters.nb_of_superdiagonals, parameters.epsfcn) <0)
00451         return HybridNonLinearSolverSpace::UserAksed;
00452     nfev += std::min(parameters.nb_of_subdiagonals+parameters.nb_of_superdiagonals+ 1, n);
00453 
00454     wa2 = fjac.colwise().blueNorm();
00455 
00456     /* on the first iteration and if external scaling is not used, scale according */
00457     /* to the norms of the columns of the initial jacobian. */
00458     if (iter == 1) {
00459         if (!useExternalScaling)
00460             for (j = 0; j < n; ++j)
00461                 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
00462 
00463         /* on the first iteration, calculate the norm of the scaled x */
00464         /* and initialize the step bound delta. */
00465         xnorm = diag.cwiseProduct(x).stableNorm();
00466         delta = parameters.factor * xnorm;
00467         if (delta == 0.)
00468             delta = parameters.factor;
00469     }
00470 
00471     /* compute the qr factorization of the jacobian. */
00472     wa2 = fjac.colwise().blueNorm();
00473     HouseholderQR<JacobianType> qrfac(fjac); // no pivoting:
00474 
00475     /* copy the triangular factor of the qr factorization into r. */
00476     R = qrfac.matrixQR();
00477 
00478     /* accumulate the orthogonal factor in fjac. */
00479     fjac = qrfac.householderQ();
00480 
00481     /* form (q transpose)*fvec and store in qtf. */
00482     qtf = fjac.transpose() * fvec;
00483 
00484     /* rescale if necessary. */
00485     if (!useExternalScaling)
00486         diag = diag.cwiseMax(wa2);
00487 
00488     while (true) {
00489         /* determine the direction p. */
00490         internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00491 
00492         /* store the direction p and x + p. calculate the norm of p. */
00493         wa1 = -wa1;
00494         wa2 = x + wa1;
00495         pnorm = diag.cwiseProduct(wa1).stableNorm();
00496 
00497         /* on the first iteration, adjust the initial step bound. */
00498         if (iter == 1)
00499             delta = std::min(delta,pnorm);
00500 
00501         /* evaluate the function at x + p and calculate its norm. */
00502         if ( functor(wa2, wa4) < 0)
00503             return HybridNonLinearSolverSpace::UserAksed;
00504         ++nfev;
00505         fnorm1 = wa4.stableNorm();
00506 
00507         /* compute the scaled actual reduction. */
00508         actred = -1.;
00509         if (fnorm1 < fnorm) /* Computing 2nd power */
00510             actred = 1. - internal::abs2(fnorm1 / fnorm);
00511 
00512         /* compute the scaled predicted reduction. */
00513         wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00514         temp = wa3.stableNorm();
00515         prered = 0.;
00516         if (temp < fnorm) /* Computing 2nd power */
00517             prered = 1. - internal::abs2(temp / fnorm);
00518 
00519         /* compute the ratio of the actual to the predicted reduction. */
00520         ratio = 0.;
00521         if (prered > 0.)
00522             ratio = actred / prered;
00523 
00524         /* update the step bound. */
00525         if (ratio < Scalar(.1)) {
00526             ncsuc = 0;
00527             ++ncfail;
00528             delta = Scalar(.5) * delta;
00529         } else {
00530             ncfail = 0;
00531             ++ncsuc;
00532             if (ratio >= Scalar(.5) || ncsuc > 1)
00533                 delta = std::max(delta, pnorm / Scalar(.5));
00534             if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00535                 delta = pnorm / Scalar(.5);
00536             }
00537         }
00538 
00539         /* test for successful iteration. */
00540         if (ratio >= Scalar(1e-4)) {
00541             /* successful iteration. update x, fvec, and their norms. */
00542             x = wa2;
00543             wa2 = diag.cwiseProduct(x);
00544             fvec = wa4;
00545             xnorm = wa2.stableNorm();
00546             fnorm = fnorm1;
00547             ++iter;
00548         }
00549 
00550         /* determine the progress of the iteration. */
00551         ++nslow1;
00552         if (actred >= Scalar(.001))
00553             nslow1 = 0;
00554         if (jeval)
00555             ++nslow2;
00556         if (actred >= Scalar(.1))
00557             nslow2 = 0;
00558 
00559         /* test for convergence. */
00560         if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00561             return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00562 
00563         /* tests for termination and stringent tolerances. */
00564         if (nfev >= parameters.maxfev)
00565             return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00566         if (Scalar(.1) * std::max(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00567             return HybridNonLinearSolverSpace::TolTooSmall;
00568         if (nslow2 == 5)
00569             return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00570         if (nslow1 == 10)
00571             return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00572 
00573         /* criterion for recalculating jacobian. */
00574         if (ncfail == 2)
00575             break; // leave inner loop and go for the next outer loop iteration
00576 
00577         /* calculate the rank one modification to the jacobian */
00578         /* and update qtf if necessary. */
00579         wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00580         wa2 = fjac.transpose() * wa4;
00581         if (ratio >= Scalar(1e-4))
00582             qtf = wa2;
00583         wa2 = (wa2-wa3)/pnorm;
00584 
00585         /* compute the qr factorization of the updated jacobian. */
00586         internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00587         internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00588         internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00589 
00590         jeval = false;
00591     }
00592     return HybridNonLinearSolverSpace::Running;
00593 }
00594 
00595 template<typename FunctorType, typename Scalar>
00596 HybridNonLinearSolverSpace::Status
00597 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiff(FVectorType  &x)
00598 {
00599     HybridNonLinearSolverSpace::Status status = solveNumericalDiffInit(x);
00600     if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00601         return status;
00602     while (status==HybridNonLinearSolverSpace::Running)
00603         status = solveNumericalDiffOneStep(x);
00604     return status;
00605 }
00606 
00607 //vim: ai ts=4 sts=4 et sw=4
00608 #endif // EIGEN_HYBRIDNONLINEARSOLVER_H
00609 


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:31:24