00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #ifndef EIGEN_HYBRIDNONLINEARSOLVER_H
00014 #define EIGEN_HYBRIDNONLINEARSOLVER_H
00015
00016 namespace Eigen {
00017
00018 namespace HybridNonLinearSolverSpace {
00019 enum Status {
00020 Running = -1,
00021 ImproperInputParameters = 0,
00022 RelativeErrorTooSmall = 1,
00023 TooManyFunctionEvaluation = 2,
00024 TolTooSmall = 3,
00025 NotMakingProgressJacobian = 4,
00026 NotMakingProgressIterations = 5,
00027 UserAsked = 6
00028 };
00029 }
00030
00042 template<typename FunctorType, typename Scalar=double>
00043 class HybridNonLinearSolver
00044 {
00045 public:
00046 typedef DenseIndex Index;
00047
00048 HybridNonLinearSolver(FunctorType &_functor)
00049 : functor(_functor) { nfev=njev=iter = 0; fnorm= 0.; useExternalScaling=false;}
00050
00051 struct Parameters {
00052 Parameters()
00053 : factor(Scalar(100.))
00054 , maxfev(1000)
00055 , xtol(internal::sqrt(NumTraits<Scalar>::epsilon()))
00056 , nb_of_subdiagonals(-1)
00057 , nb_of_superdiagonals(-1)
00058 , epsfcn(Scalar(0.)) {}
00059 Scalar factor;
00060 Index maxfev;
00061 Scalar xtol;
00062 Index nb_of_subdiagonals;
00063 Index nb_of_superdiagonals;
00064 Scalar epsfcn;
00065 };
00066 typedef Matrix< Scalar, Dynamic, 1 > FVectorType;
00067 typedef Matrix< Scalar, Dynamic, Dynamic > JacobianType;
00068
00069 typedef Matrix< Scalar, Dynamic, Dynamic > UpperTriangularType;
00070
00071 HybridNonLinearSolverSpace::Status hybrj1(
00072 FVectorType &x,
00073 const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00074 );
00075
00076 HybridNonLinearSolverSpace::Status solveInit(FVectorType &x);
00077 HybridNonLinearSolverSpace::Status solveOneStep(FVectorType &x);
00078 HybridNonLinearSolverSpace::Status solve(FVectorType &x);
00079
00080 HybridNonLinearSolverSpace::Status hybrd1(
00081 FVectorType &x,
00082 const Scalar tol = internal::sqrt(NumTraits<Scalar>::epsilon())
00083 );
00084
00085 HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType &x);
00086 HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType &x);
00087 HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType &x);
00088
00089 void resetParameters(void) { parameters = Parameters(); }
00090 Parameters parameters;
00091 FVectorType fvec, qtf, diag;
00092 JacobianType fjac;
00093 UpperTriangularType R;
00094 Index nfev;
00095 Index njev;
00096 Index iter;
00097 Scalar fnorm;
00098 bool useExternalScaling;
00099 private:
00100 FunctorType &functor;
00101 Index n;
00102 Scalar sum;
00103 bool sing;
00104 Scalar temp;
00105 Scalar delta;
00106 bool jeval;
00107 Index ncsuc;
00108 Scalar ratio;
00109 Scalar pnorm, xnorm, fnorm1;
00110 Index nslow1, nslow2;
00111 Index ncfail;
00112 Scalar actred, prered;
00113 FVectorType wa1, wa2, wa3, wa4;
00114
00115 HybridNonLinearSolver& operator=(const HybridNonLinearSolver&);
00116 };
00117
00118
00119
00120 template<typename FunctorType, typename Scalar>
00121 HybridNonLinearSolverSpace::Status
00122 HybridNonLinearSolver<FunctorType,Scalar>::hybrj1(
00123 FVectorType &x,
00124 const Scalar tol
00125 )
00126 {
00127 n = x.size();
00128
00129
00130 if (n <= 0 || tol < 0.)
00131 return HybridNonLinearSolverSpace::ImproperInputParameters;
00132
00133 resetParameters();
00134 parameters.maxfev = 100*(n+1);
00135 parameters.xtol = tol;
00136 diag.setConstant(n, 1.);
00137 useExternalScaling = true;
00138 return solve(x);
00139 }
00140
00141 template<typename FunctorType, typename Scalar>
00142 HybridNonLinearSolverSpace::Status
00143 HybridNonLinearSolver<FunctorType,Scalar>::solveInit(FVectorType &x)
00144 {
00145 n = x.size();
00146
00147 wa1.resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
00148 fvec.resize(n);
00149 qtf.resize(n);
00150 fjac.resize(n, n);
00151 if (!useExternalScaling)
00152 diag.resize(n);
00153 assert( (!useExternalScaling || diag.size()==n) || "When useExternalScaling is set, the caller must provide a valid 'diag'");
00154
00155
00156 nfev = 0;
00157 njev = 0;
00158
00159
00160 if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0. )
00161 return HybridNonLinearSolverSpace::ImproperInputParameters;
00162 if (useExternalScaling)
00163 for (Index j = 0; j < n; ++j)
00164 if (diag[j] <= 0.)
00165 return HybridNonLinearSolverSpace::ImproperInputParameters;
00166
00167
00168
00169 nfev = 1;
00170 if ( functor(x, fvec) < 0)
00171 return HybridNonLinearSolverSpace::UserAsked;
00172 fnorm = fvec.stableNorm();
00173
00174
00175 iter = 1;
00176 ncsuc = 0;
00177 ncfail = 0;
00178 nslow1 = 0;
00179 nslow2 = 0;
00180
00181 return HybridNonLinearSolverSpace::Running;
00182 }
00183
00184 template<typename FunctorType, typename Scalar>
00185 HybridNonLinearSolverSpace::Status
00186 HybridNonLinearSolver<FunctorType,Scalar>::solveOneStep(FVectorType &x)
00187 {
00188 assert(x.size()==n);
00189
00190 Index j;
00191 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00192
00193 jeval = true;
00194
00195
00196 if ( functor.df(x, fjac) < 0)
00197 return HybridNonLinearSolverSpace::UserAsked;
00198 ++njev;
00199
00200 wa2 = fjac.colwise().blueNorm();
00201
00202
00203
00204 if (iter == 1) {
00205 if (!useExternalScaling)
00206 for (j = 0; j < n; ++j)
00207 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
00208
00209
00210
00211 xnorm = diag.cwiseProduct(x).stableNorm();
00212 delta = parameters.factor * xnorm;
00213 if (delta == 0.)
00214 delta = parameters.factor;
00215 }
00216
00217
00218 HouseholderQR<JacobianType> qrfac(fjac);
00219
00220
00221 R = qrfac.matrixQR();
00222
00223
00224 fjac = qrfac.householderQ();
00225
00226
00227 qtf = fjac.transpose() * fvec;
00228
00229
00230 if (!useExternalScaling)
00231 diag = diag.cwiseMax(wa2);
00232
00233 while (true) {
00234
00235 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00236
00237
00238 wa1 = -wa1;
00239 wa2 = x + wa1;
00240 pnorm = diag.cwiseProduct(wa1).stableNorm();
00241
00242
00243 if (iter == 1)
00244 delta = (std::min)(delta,pnorm);
00245
00246
00247 if ( functor(wa2, wa4) < 0)
00248 return HybridNonLinearSolverSpace::UserAsked;
00249 ++nfev;
00250 fnorm1 = wa4.stableNorm();
00251
00252
00253 actred = -1.;
00254 if (fnorm1 < fnorm)
00255 actred = 1. - internal::abs2(fnorm1 / fnorm);
00256
00257
00258 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00259 temp = wa3.stableNorm();
00260 prered = 0.;
00261 if (temp < fnorm)
00262 prered = 1. - internal::abs2(temp / fnorm);
00263
00264
00265 ratio = 0.;
00266 if (prered > 0.)
00267 ratio = actred / prered;
00268
00269
00270 if (ratio < Scalar(.1)) {
00271 ncsuc = 0;
00272 ++ncfail;
00273 delta = Scalar(.5) * delta;
00274 } else {
00275 ncfail = 0;
00276 ++ncsuc;
00277 if (ratio >= Scalar(.5) || ncsuc > 1)
00278 delta = (std::max)(delta, pnorm / Scalar(.5));
00279 if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00280 delta = pnorm / Scalar(.5);
00281 }
00282 }
00283
00284
00285 if (ratio >= Scalar(1e-4)) {
00286
00287 x = wa2;
00288 wa2 = diag.cwiseProduct(x);
00289 fvec = wa4;
00290 xnorm = wa2.stableNorm();
00291 fnorm = fnorm1;
00292 ++iter;
00293 }
00294
00295
00296 ++nslow1;
00297 if (actred >= Scalar(.001))
00298 nslow1 = 0;
00299 if (jeval)
00300 ++nslow2;
00301 if (actred >= Scalar(.1))
00302 nslow2 = 0;
00303
00304
00305 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00306 return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00307
00308
00309 if (nfev >= parameters.maxfev)
00310 return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00311 if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00312 return HybridNonLinearSolverSpace::TolTooSmall;
00313 if (nslow2 == 5)
00314 return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00315 if (nslow1 == 10)
00316 return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00317
00318
00319 if (ncfail == 2)
00320 break;
00321
00322
00323
00324 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00325 wa2 = fjac.transpose() * wa4;
00326 if (ratio >= Scalar(1e-4))
00327 qtf = wa2;
00328 wa2 = (wa2-wa3)/pnorm;
00329
00330
00331 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00332 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00333 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00334
00335 jeval = false;
00336 }
00337 return HybridNonLinearSolverSpace::Running;
00338 }
00339
00340 template<typename FunctorType, typename Scalar>
00341 HybridNonLinearSolverSpace::Status
00342 HybridNonLinearSolver<FunctorType,Scalar>::solve(FVectorType &x)
00343 {
00344 HybridNonLinearSolverSpace::Status status = solveInit(x);
00345 if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00346 return status;
00347 while (status==HybridNonLinearSolverSpace::Running)
00348 status = solveOneStep(x);
00349 return status;
00350 }
00351
00352
00353
00354 template<typename FunctorType, typename Scalar>
00355 HybridNonLinearSolverSpace::Status
00356 HybridNonLinearSolver<FunctorType,Scalar>::hybrd1(
00357 FVectorType &x,
00358 const Scalar tol
00359 )
00360 {
00361 n = x.size();
00362
00363
00364 if (n <= 0 || tol < 0.)
00365 return HybridNonLinearSolverSpace::ImproperInputParameters;
00366
00367 resetParameters();
00368 parameters.maxfev = 200*(n+1);
00369 parameters.xtol = tol;
00370
00371 diag.setConstant(n, 1.);
00372 useExternalScaling = true;
00373 return solveNumericalDiff(x);
00374 }
00375
00376 template<typename FunctorType, typename Scalar>
00377 HybridNonLinearSolverSpace::Status
00378 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffInit(FVectorType &x)
00379 {
00380 n = x.size();
00381
00382 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00383 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00384
00385 wa1.resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
00386 qtf.resize(n);
00387 fjac.resize(n, n);
00388 fvec.resize(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
00393
00394 nfev = 0;
00395 njev = 0;
00396
00397
00398 if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.nb_of_subdiagonals< 0 || parameters.nb_of_superdiagonals< 0 || parameters.factor <= 0. )
00399 return HybridNonLinearSolverSpace::ImproperInputParameters;
00400 if (useExternalScaling)
00401 for (Index j = 0; j < n; ++j)
00402 if (diag[j] <= 0.)
00403 return HybridNonLinearSolverSpace::ImproperInputParameters;
00404
00405
00406
00407 nfev = 1;
00408 if ( functor(x, fvec) < 0)
00409 return HybridNonLinearSolverSpace::UserAsked;
00410 fnorm = fvec.stableNorm();
00411
00412
00413 iter = 1;
00414 ncsuc = 0;
00415 ncfail = 0;
00416 nslow1 = 0;
00417 nslow2 = 0;
00418
00419 return HybridNonLinearSolverSpace::Running;
00420 }
00421
00422 template<typename FunctorType, typename Scalar>
00423 HybridNonLinearSolverSpace::Status
00424 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(FVectorType &x)
00425 {
00426 assert(x.size()==n);
00427
00428 Index j;
00429 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
00430
00431 jeval = true;
00432 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
00433 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
00434
00435
00436 if (internal::fdjac1(functor, x, fvec, fjac, parameters.nb_of_subdiagonals, parameters.nb_of_superdiagonals, parameters.epsfcn) <0)
00437 return HybridNonLinearSolverSpace::UserAsked;
00438 nfev += (std::min)(parameters.nb_of_subdiagonals+parameters.nb_of_superdiagonals+ 1, n);
00439
00440 wa2 = fjac.colwise().blueNorm();
00441
00442
00443
00444 if (iter == 1) {
00445 if (!useExternalScaling)
00446 for (j = 0; j < n; ++j)
00447 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
00448
00449
00450
00451 xnorm = diag.cwiseProduct(x).stableNorm();
00452 delta = parameters.factor * xnorm;
00453 if (delta == 0.)
00454 delta = parameters.factor;
00455 }
00456
00457
00458 HouseholderQR<JacobianType> qrfac(fjac);
00459
00460
00461 R = qrfac.matrixQR();
00462
00463
00464 fjac = qrfac.householderQ();
00465
00466
00467 qtf = fjac.transpose() * fvec;
00468
00469
00470 if (!useExternalScaling)
00471 diag = diag.cwiseMax(wa2);
00472
00473 while (true) {
00474
00475 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
00476
00477
00478 wa1 = -wa1;
00479 wa2 = x + wa1;
00480 pnorm = diag.cwiseProduct(wa1).stableNorm();
00481
00482
00483 if (iter == 1)
00484 delta = (std::min)(delta,pnorm);
00485
00486
00487 if ( functor(wa2, wa4) < 0)
00488 return HybridNonLinearSolverSpace::UserAsked;
00489 ++nfev;
00490 fnorm1 = wa4.stableNorm();
00491
00492
00493 actred = -1.;
00494 if (fnorm1 < fnorm)
00495 actred = 1. - internal::abs2(fnorm1 / fnorm);
00496
00497
00498 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
00499 temp = wa3.stableNorm();
00500 prered = 0.;
00501 if (temp < fnorm)
00502 prered = 1. - internal::abs2(temp / fnorm);
00503
00504
00505 ratio = 0.;
00506 if (prered > 0.)
00507 ratio = actred / prered;
00508
00509
00510 if (ratio < Scalar(.1)) {
00511 ncsuc = 0;
00512 ++ncfail;
00513 delta = Scalar(.5) * delta;
00514 } else {
00515 ncfail = 0;
00516 ++ncsuc;
00517 if (ratio >= Scalar(.5) || ncsuc > 1)
00518 delta = (std::max)(delta, pnorm / Scalar(.5));
00519 if (internal::abs(ratio - 1.) <= Scalar(.1)) {
00520 delta = pnorm / Scalar(.5);
00521 }
00522 }
00523
00524
00525 if (ratio >= Scalar(1e-4)) {
00526
00527 x = wa2;
00528 wa2 = diag.cwiseProduct(x);
00529 fvec = wa4;
00530 xnorm = wa2.stableNorm();
00531 fnorm = fnorm1;
00532 ++iter;
00533 }
00534
00535
00536 ++nslow1;
00537 if (actred >= Scalar(.001))
00538 nslow1 = 0;
00539 if (jeval)
00540 ++nslow2;
00541 if (actred >= Scalar(.1))
00542 nslow2 = 0;
00543
00544
00545 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
00546 return HybridNonLinearSolverSpace::RelativeErrorTooSmall;
00547
00548
00549 if (nfev >= parameters.maxfev)
00550 return HybridNonLinearSolverSpace::TooManyFunctionEvaluation;
00551 if (Scalar(.1) * (std::max)(Scalar(.1) * delta, pnorm) <= NumTraits<Scalar>::epsilon() * xnorm)
00552 return HybridNonLinearSolverSpace::TolTooSmall;
00553 if (nslow2 == 5)
00554 return HybridNonLinearSolverSpace::NotMakingProgressJacobian;
00555 if (nslow1 == 10)
00556 return HybridNonLinearSolverSpace::NotMakingProgressIterations;
00557
00558
00559 if (ncfail == 2)
00560 break;
00561
00562
00563
00564 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
00565 wa2 = fjac.transpose() * wa4;
00566 if (ratio >= Scalar(1e-4))
00567 qtf = wa2;
00568 wa2 = (wa2-wa3)/pnorm;
00569
00570
00571 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
00572 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
00573 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
00574
00575 jeval = false;
00576 }
00577 return HybridNonLinearSolverSpace::Running;
00578 }
00579
00580 template<typename FunctorType, typename Scalar>
00581 HybridNonLinearSolverSpace::Status
00582 HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiff(FVectorType &x)
00583 {
00584 HybridNonLinearSolverSpace::Status status = solveNumericalDiffInit(x);
00585 if (status==HybridNonLinearSolverSpace::ImproperInputParameters)
00586 return status;
00587 while (status==HybridNonLinearSolverSpace::Running)
00588 status = solveNumericalDiffOneStep(x);
00589 return status;
00590 }
00591
00592 }
00593
00594 #endif // EIGEN_HYBRIDNONLINEARSOLVER_H
00595
00596