13 #ifndef EIGEN_HYBRIDNONLINEARSOLVER_H 14 #define EIGEN_HYBRIDNONLINEARSOLVER_H 18 namespace HybridNonLinearSolverSpace {
42 template<
typename FunctorType,
typename Scalar=
double>
49 : functor(_functor) { nfev=njev=iter = 0; fnorm= 0.; useExternalScaling=
false;}
53 : factor(Scalar(100.))
56 , nb_of_subdiagonals(-1)
57 , nb_of_superdiagonals(-1)
58 , epsfcn(Scalar(0.)) {}
91 FVectorType fvec,
qtf, diag;
93 UpperTriangularType
R;
113 FVectorType wa1, wa2, wa3,
wa4;
120 template<
typename FunctorType,
typename Scalar>
130 if (n <= 0 || tol < 0.)
134 parameters.maxfev = 100*(n+1);
135 parameters.xtol = tol;
136 diag.setConstant(n, 1.);
137 useExternalScaling =
true;
141 template<
typename FunctorType,
typename Scalar>
147 wa1.
resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
151 if (!useExternalScaling)
153 eigen_assert( (!useExternalScaling || diag.size()==n) &&
"When useExternalScaling is set, the caller must provide a valid 'diag'");
160 if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0. )
162 if (useExternalScaling)
163 for (Index j = 0; j < n; ++j)
170 if ( functor(x, fvec) < 0)
172 fnorm = fvec.stableNorm();
184 template<
typename FunctorType,
typename Scalar>
193 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
198 if ( functor.df(x, fjac) < 0)
202 wa2 = fjac.colwise().blueNorm();
207 if (!useExternalScaling)
208 for (j = 0; j < n; ++j)
209 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
213 xnorm = diag.cwiseProduct(x).stableNorm();
214 delta = parameters.factor * xnorm;
216 delta = parameters.factor;
229 qtf = fjac.transpose() * fvec;
232 if (!useExternalScaling)
233 diag = diag.cwiseMax(wa2);
237 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
242 pnorm = diag.cwiseProduct(wa1).stableNorm();
246 delta = (std::min)(delta,pnorm);
249 if ( functor(wa2, wa4) < 0)
252 fnorm1 = wa4.stableNorm();
260 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
261 temp = wa3.stableNorm();
269 ratio = actred / prered;
272 if (ratio < Scalar(.1)) {
275 delta = Scalar(.5) * delta;
279 if (ratio >= Scalar(.5) || ncsuc > 1)
280 delta = (
std::max)(delta, pnorm / Scalar(.5));
281 if (
abs(ratio - 1.) <= Scalar(.1)) {
282 delta = pnorm / Scalar(.5);
287 if (ratio >= Scalar(1e-4)) {
290 wa2 = diag.cwiseProduct(x);
292 xnorm = wa2.stableNorm();
299 if (actred >= Scalar(.001))
303 if (actred >= Scalar(.1))
307 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
311 if (nfev >= parameters.maxfev)
326 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
327 wa2 = fjac.transpose() * wa4;
328 if (ratio >= Scalar(1e-4))
330 wa2 = (wa2-wa3)/pnorm;
333 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
334 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
335 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
342 template<
typename FunctorType,
typename Scalar>
350 status = solveOneStep(x);
356 template<
typename FunctorType,
typename Scalar>
366 if (n <= 0 || tol < 0.)
370 parameters.maxfev = 200*(n+1);
371 parameters.xtol = tol;
373 diag.setConstant(n, 1.);
374 useExternalScaling =
true;
375 return solveNumericalDiff(x);
378 template<
typename FunctorType,
typename Scalar>
384 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
385 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
387 wa1.
resize(n); wa2.resize(n); wa3.resize(n); wa4.resize(n);
391 if (!useExternalScaling)
393 eigen_assert( (!useExternalScaling || diag.size()==n) &&
"When useExternalScaling is set, the caller must provide a valid 'diag'");
400 if (n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.nb_of_subdiagonals< 0 || parameters.nb_of_superdiagonals< 0 || parameters.factor <= 0. )
402 if (useExternalScaling)
403 for (Index j = 0; j < n; ++j)
410 if ( functor(x, fvec) < 0)
412 fnorm = fvec.stableNorm();
424 template<
typename FunctorType,
typename Scalar>
434 std::vector<JacobiRotation<Scalar> > v_givens(n), w_givens(n);
437 if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
438 if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
441 if (
internal::fdjac1(functor, x, fvec, fjac, parameters.nb_of_subdiagonals, parameters.nb_of_superdiagonals, parameters.epsfcn) <0)
443 nfev += (std::min)(parameters.nb_of_subdiagonals+parameters.nb_of_superdiagonals+ 1, n);
445 wa2 = fjac.colwise().blueNorm();
450 if (!useExternalScaling)
451 for (j = 0; j < n; ++j)
452 diag[j] = (wa2[j]==0.) ? 1. : wa2[j];
456 xnorm = diag.cwiseProduct(x).stableNorm();
457 delta = parameters.factor * xnorm;
459 delta = parameters.factor;
472 qtf = fjac.transpose() * fvec;
475 if (!useExternalScaling)
476 diag = diag.cwiseMax(wa2);
480 internal::dogleg<Scalar>(R, diag, qtf, delta, wa1);
485 pnorm = diag.cwiseProduct(wa1).stableNorm();
489 delta = (std::min)(delta,pnorm);
492 if ( functor(wa2, wa4) < 0)
495 fnorm1 = wa4.stableNorm();
503 wa3 = R.template triangularView<Upper>()*wa1 + qtf;
504 temp = wa3.stableNorm();
512 ratio = actred / prered;
515 if (ratio < Scalar(.1)) {
518 delta = Scalar(.5) * delta;
522 if (ratio >= Scalar(.5) || ncsuc > 1)
523 delta = (
std::max)(delta, pnorm / Scalar(.5));
524 if (
abs(ratio - 1.) <= Scalar(.1)) {
525 delta = pnorm / Scalar(.5);
530 if (ratio >= Scalar(1e-4)) {
533 wa2 = diag.cwiseProduct(x);
535 xnorm = wa2.stableNorm();
542 if (actred >= Scalar(.001))
546 if (actred >= Scalar(.1))
550 if (delta <= parameters.xtol * xnorm || fnorm == 0.)
554 if (nfev >= parameters.maxfev)
569 wa1 = diag.cwiseProduct( diag.cwiseProduct(wa1)/pnorm );
570 wa2 = fjac.transpose() * wa4;
571 if (ratio >= Scalar(1e-4))
573 wa2 = (wa2-wa3)/pnorm;
576 internal::r1updt<Scalar>(R, wa1, v_givens, w_givens, wa2, wa3, &sing);
577 internal::r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
578 internal::r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
585 template<
typename FunctorType,
typename Scalar>
593 status = solveNumericalDiffOneStep(x);
599 #endif // EIGEN_HYBRIDNONLINEARSOLVER_H Matrix< Scalar, Dynamic, Dynamic > UpperTriangularType
HouseholderSequenceType householderQ() const
HybridNonLinearSolver(FunctorType &_functor)
HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType &x)
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
void resetParameters(void)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const AbsReturnType abs() const
HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType &x)
HybridNonLinearSolverSpace::Status solve(FVectorType &x)
HybridNonLinearSolverSpace::Status solveOneStep(FVectorType &x)
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half() max(const half &a, const half &b)
Matrix< Scalar, Dynamic, Dynamic > JacobianType
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Index nb_of_superdiagonals
HybridNonLinearSolverSpace::Status hybrd1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType &x)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
HybridNonLinearSolverSpace::Status hybrj1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
HybridNonLinearSolverSpace::Status solveInit(FVectorType &x)
Householder QR decomposition of a matrix.
Matrix< Scalar, Dynamic, 1 > FVectorType
Finds a zero of a system of n nonlinear functions in n variables by a modification of the Powell hybr...
const MatrixType & matrixQR() const
DenseIndex fdjac1(const FunctorType &Functor, Matrix< Scalar, Dynamic, 1 > &x, Matrix< Scalar, Dynamic, 1 > &fvec, Matrix< Scalar, Dynamic, Dynamic > &fjac, DenseIndex ml, DenseIndex mu, Scalar epsfcn)