00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
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 UserAsked = 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;
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
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
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
00169 nfev = 0;
00170 njev = 0;
00171
00172
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
00181
00182 nfev = 1;
00183 if ( functor(x, fvec) < 0)
00184 return HybridNonLinearSolverSpace::UserAsked;
00185 fnorm = fvec.stableNorm();
00186
00187
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);
00202
00203 Index j;
00204 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00205
00206 jeval = true;
00207
00208
00209 if ( functor.df(x, fjac) < 0)
00210 return HybridNonLinearSolverSpace::UserAsked;
00211 ++njev;
00212
00213 wa2 = fjac.colwise().blueNorm();
00214
00215
00216
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
00223
00224 xnorm = diag.cwiseProduct(x).stableNorm();
00225 delta = parameters.factor * xnorm;
00226 if (delta == 0.)
00227 delta = parameters.factor;
00228 }
00229
00230
00231 HouseholderQR<JacobianType> qrfac(fjac);
00232
00233
00234 R = qrfac.matrixQR();
00235
00236
00237 fjac = qrfac.householderQ();
00238
00239
00240 qtf = fjac.transpose() * fvec;
00241
00242
00243 if (!useExternalScaling)
00244 diag = diag.cwiseMax(wa2);
00245
00246 while (true) {
00247
00248 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00249
00250
00251 wa1 = -wa1;
00252 wa2 = x + wa1;
00253 pnorm = diag.cwiseProduct(wa1).stableNorm();
00254
00255
00256 if (iter == 1)
00257 delta = (std::min)(delta,pnorm);
00258
00259
00260 if ( functor(wa2, wa4) < 0)
00261 return HybridNonLinearSolverSpace::UserAsked;
00262 ++nfev;
00263 fnorm1 = wa4.stableNorm();
00264
00265
00266 actred = -1.;
00267 if (fnorm1 < fnorm)
00268 actred = 1. - internal::abs2(fnorm1 / fnorm);
00269
00270
00271 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00272 temp = wa3.stableNorm();
00273 prered = 0.;
00274 if (temp < fnorm)
00275 prered = 1. - internal::abs2(temp / fnorm);
00276
00277
00278 ratio = 0.;
00279 if (prered > 0.)
00280 ratio = actred / prered;
00281
00282
00283 if (ratio < Scalar(.1)) {
00284 ncsuc = 0;
00285 ++ncfail;
00286 delta = Scalar(.5) * delta;
00287 } else {
00288 ncfail = 0;
00289 ++ncsuc;
00290 if (ratio >= Scalar(.5) || ncsuc > 1)
00291 delta = (std::max)(delta, pnorm / Scalar(.5));
00292 if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00293 delta = pnorm / Scalar(.5);
00294 }
00295 }
00296
00297
00298 if (ratio >= Scalar(1e-4)) {
00299
00300 x = wa2;
00301 wa2 = diag.cwiseProduct(x);
00302 fvec = wa4;
00303 xnorm = wa2.stableNorm();
00304 fnorm = fnorm1;
00305 ++iter;
00306 }
00307
00308
00309 ++nslow1;
00310 if (actred >= Scalar(.001))
00311 nslow1 = 0;
00312 if (jeval)
00313 ++nslow2;
00314 if (actred >= Scalar(.1))
00315 nslow2 = 0;
00316
00317
00318 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00319 return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00320
00321
00322 if (nfev >= parameters.maxfev)
00323 return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00324 if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00325 return HybridNonLinearSolverSpace::TolTooSmall;
00326 if (nslow2 == 5)
00327 return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00328 if (nslow1 == 10)
00329 return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00330
00331
00332 if (ncfail == 2)
00333 break;
00334
00335
00336
00337 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00338 wa2 = fjac.transpose() * wa4;
00339 if (ratio >= Scalar(1e-4))
00340 qtf = wa2;
00341 wa2 = (wa2-wa3)/pnorm;
00342
00343
00344 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00345 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00346 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00347
00348 jeval = false;
00349 }
00350 return HybridNonLinearSolverSpace::Running;
00351 }
00352
00353 template<typename FunctorType, typename Scalar>
00354 HybridNonLinearSolverSpace::Status
00355 HybridNonLinearSolver<FunctorType,Scalar>::solve(FVectorType &x)
00356 {
00357 HybridNonLinearSolverSpace::Status status = solveInit(x);
00358 if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00359 return status;
00360 while (status==HybridNonLinearSolverSpace::Running)
00361 status = solveOneStep(x);
00362 return status;
00363 }
00364
00365
00366
00367 template<typename FunctorType, typename Scalar>
00368 HybridNonLinearSolverSpace::Status
00369 HybridNonLinearSolver<FunctorType,Scalar>::hybrd1(
00370 FVectorType &x,
00371 const Scalar tol
00372 )
00373 {
00374 n = x.size();
00375
00376
00377 if (n <= 0 || tol < 0.)
00378 return HybridNonLinearSolverSpace::ImproperInputParameters;
00379
00380 resetParameters();
00381 parameters.maxfev = 200*(n+1);
00382 parameters.xtol = tol;
00383
00384 diag.setConstant(n, 1.);
00385 useExternalScaling = true;
00386 return solveNumericalDiff(x);
00387 }
00388
00389 template<typename FunctorType, typename Scalar>
00390 HybridNonLinearSolverSpace::Status
00391 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffInit(FVectorType &x)
00392 {
00393 n = x.size();
00394
00395 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00396 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00397
00398 wa1.resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
00399 qtf.resize(n);
00400 fjac.resize(n, n);
00401 fvec.resize(n);
00402 if (!useExternalScaling)
00403 diag.resize(n);
00404 assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00405
00406
00407 nfev = 0;
00408 njev = 0;
00409
00410
00411 if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.nb_of_subdiagonals< 0 || parameters.nb_of_superdiagonals< 0 || parameters.factor <= 0. )
00412 return HybridNonLinearSolverSpace::ImproperInputParameters;
00413 if (useExternalScaling)
00414 for (Index j = 0; j < n; ++j)
00415 if (diag[j] <= 0.)
00416 return HybridNonLinearSolverSpace::ImproperInputParameters;
00417
00418
00419
00420 nfev = 1;
00421 if ( functor(x, fvec) < 0)
00422 return HybridNonLinearSolverSpace::UserAsked;
00423 fnorm = fvec.stableNorm();
00424
00425
00426 iter = 1;
00427 ncsuc = 0;
00428 ncfail = 0;
00429 nslow1 = 0;
00430 nslow2 = 0;
00431
00432 return HybridNonLinearSolverSpace::Running;
00433 }
00434
00435 template<typename FunctorType, typename Scalar>
00436 HybridNonLinearSolverSpace::Status
00437 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(FVectorType &x)
00438 {
00439 assert(x.size()==n);
00440
00441 Index j;
00442 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00443
00444 jeval = true;
00445 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00446 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00447
00448
00449 if (internal::fdjac1(functor, x, fvec, fjac, parameters.nb_of_subdiagonals, parameters.nb_of_superdiagonals, parameters.epsfcn) <0)
00450 return HybridNonLinearSolverSpace::UserAsked;
00451 nfev += (std::min)(parameters.nb_of_subdiagonals+parameters.nb_of_superdiagonals+ 1, n);
00452
00453 wa2 = fjac.colwise().blueNorm();
00454
00455
00456
00457 if (iter == 1) {
00458 if (!useExternalScaling)
00459 for (j = 0; j < n; ++j)
00460 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
00461
00462
00463
00464 xnorm = diag.cwiseProduct(x).stableNorm();
00465 delta = parameters.factor * xnorm;
00466 if (delta == 0.)
00467 delta = parameters.factor;
00468 }
00469
00470
00471 HouseholderQR<JacobianType> qrfac(fjac);
00472
00473
00474 R = qrfac.matrixQR();
00475
00476
00477 fjac = qrfac.householderQ();
00478
00479
00480 qtf = fjac.transpose() * fvec;
00481
00482
00483 if (!useExternalScaling)
00484 diag = diag.cwiseMax(wa2);
00485
00486 while (true) {
00487
00488 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00489
00490
00491 wa1 = -wa1;
00492 wa2 = x + wa1;
00493 pnorm = diag.cwiseProduct(wa1).stableNorm();
00494
00495
00496 if (iter == 1)
00497 delta = (std::min)(delta,pnorm);
00498
00499
00500 if ( functor(wa2, wa4) < 0)
00501 return HybridNonLinearSolverSpace::UserAsked;
00502 ++nfev;
00503 fnorm1 = wa4.stableNorm();
00504
00505
00506 actred = -1.;
00507 if (fnorm1 < fnorm)
00508 actred = 1. - internal::abs2(fnorm1 / fnorm);
00509
00510
00511 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00512 temp = wa3.stableNorm();
00513 prered = 0.;
00514 if (temp < fnorm)
00515 prered = 1. - internal::abs2(temp / fnorm);
00516
00517
00518 ratio = 0.;
00519 if (prered > 0.)
00520 ratio = actred / prered;
00521
00522
00523 if (ratio < Scalar(.1)) {
00524 ncsuc = 0;
00525 ++ncfail;
00526 delta = Scalar(.5) * delta;
00527 } else {
00528 ncfail = 0;
00529 ++ncsuc;
00530 if (ratio >= Scalar(.5) || ncsuc > 1)
00531 delta = (std::max)(delta, pnorm / Scalar(.5));
00532 if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00533 delta = pnorm / Scalar(.5);
00534 }
00535 }
00536
00537
00538 if (ratio >= Scalar(1e-4)) {
00539
00540 x = wa2;
00541 wa2 = diag.cwiseProduct(x);
00542 fvec = wa4;
00543 xnorm = wa2.stableNorm();
00544 fnorm = fnorm1;
00545 ++iter;
00546 }
00547
00548
00549 ++nslow1;
00550 if (actred >= Scalar(.001))
00551 nslow1 = 0;
00552 if (jeval)
00553 ++nslow2;
00554 if (actred >= Scalar(.1))
00555 nslow2 = 0;
00556
00557
00558 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00559 return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00560
00561
00562 if (nfev >= parameters.maxfev)
00563 return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00564 if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00565 return HybridNonLinearSolverSpace::TolTooSmall;
00566 if (nslow2 == 5)
00567 return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00568 if (nslow1 == 10)
00569 return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00570
00571
00572 if (ncfail == 2)
00573 break;
00574
00575
00576
00577 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00578 wa2 = fjac.transpose() * wa4;
00579 if (ratio >= Scalar(1e-4))
00580 qtf = wa2;
00581 wa2 = (wa2-wa3)/pnorm;
00582
00583
00584 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00585 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00586 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00587
00588 jeval = false;
00589 }
00590 return HybridNonLinearSolverSpace::Running;
00591 }
00592
00593 template<typename FunctorType, typename Scalar>
00594 HybridNonLinearSolverSpace::Status
00595 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiff(FVectorType &x)
00596 {
00597 HybridNonLinearSolverSpace::Status status = solveNumericalDiffInit(x);
00598 if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00599 return status;
00600 while (status==HybridNonLinearSolverSpace::Running)
00601 status = solveNumericalDiffOneStep(x);
00602 return status;
00603 }
00604
00605
00606 #endif // EIGEN_HYBRIDNONLINEARSOLVER_H
00607