13#ifndef EIGEN_HYBRIDNONLINEARSOLVER_H
14#define EIGEN_HYBRIDNONLINEARSOLVER_H
18namespace HybridNonLinearSolverSpace {
42template<
typename FunctorType,
typename Scalar=
double>
100 FunctorType &functor;
109 Scalar pnorm, xnorm, fnorm1;
110 Index nslow1, nslow2;
120template<
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;
141template<
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)
170 if ( functor(
x, fvec) < 0)
172 fnorm = fvec.stableNorm();
184template<
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;
232 if (!useExternalScaling)
233 diag = diag.cwiseMax(wa2);
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;
275 delta =
Scalar(.5) * delta;
280 delta = (std::max)(delta, pnorm /
Scalar(.5));
282 delta = pnorm /
Scalar(.5);
290 wa2 = diag.cwiseProduct(
x);
292 xnorm = wa2.stableNorm();
299 if (actred >=
Scalar(.001))
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;
330 wa2 = (wa2-wa3)/pnorm;
342template<
typename FunctorType,
typename Scalar>
350 status = solveOneStep(
x);
356template<
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);
378template<
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)
410 if ( functor(
x, fvec) < 0)
412 fnorm = fvec.stableNorm();
424template<
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;
475 if (!useExternalScaling)
476 diag = diag.cwiseMax(wa2);
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;
518 delta =
Scalar(.5) * delta;
523 delta = (std::max)(delta, pnorm /
Scalar(.5));
525 delta = pnorm /
Scalar(.5);
533 wa2 = diag.cwiseProduct(
x);
535 xnorm = wa2.stableNorm();
542 if (actred >=
Scalar(.001))
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;
573 wa2 = (wa2-wa3)/pnorm;
585template<
typename FunctorType,
typename Scalar>
593 status = solveNumericalDiffOneStep(
x);
int n
Definition BiCGSTAB_simple.cpp:1
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define eigen_assert(x)
Definition Macros.h:1037
SCALAR Scalar
Definition bench_gemm.cpp:46
Householder QR decomposition of a matrix.
Definition HouseholderQR.h:58
const MatrixType & matrixQR() const
Definition HouseholderQR.h:172
HouseholderSequenceType householderQ() const
Definition HouseholderQR.h:163
TransposeReturnType transpose() const
Transpose of the Householder sequence.
Definition HouseholderSequence.h:236
Finds a zero of a system of n nonlinear functions in n variables by a modification of the Powell hybr...
Definition HybridNonLinearSolver.h:44
HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType &x)
Definition HybridNonLinearSolver.h:587
HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType &x)
Definition HybridNonLinearSolver.h:380
void resetParameters(void)
Definition HybridNonLinearSolver.h:89
FVectorType fvec
Definition HybridNonLinearSolver.h:91
HybridNonLinearSolverSpace::Status solveInit(FVectorType &x)
Definition HybridNonLinearSolver.h:143
Index iter
Definition HybridNonLinearSolver.h:96
bool useExternalScaling
Definition HybridNonLinearSolver.h:98
UpperTriangularType R
Definition HybridNonLinearSolver.h:93
HybridNonLinearSolverSpace::Status solve(FVectorType &x)
Definition HybridNonLinearSolver.h:344
Matrix< Scalar, Dynamic, Dynamic > JacobianType
Definition HybridNonLinearSolver.h:67
Parameters parameters
Definition HybridNonLinearSolver.h:90
HybridNonLinearSolverSpace::Status solveOneStep(FVectorType &x)
Definition HybridNonLinearSolver.h:186
Index nfev
Definition HybridNonLinearSolver.h:94
Matrix< Scalar, Dynamic, 1 > FVectorType
Definition HybridNonLinearSolver.h:66
HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType &x)
Definition HybridNonLinearSolver.h:426
Scalar fnorm
Definition HybridNonLinearSolver.h:97
DenseIndex Index
Definition HybridNonLinearSolver.h:46
HybridNonLinearSolverSpace::Status hybrj1(FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
Definition HybridNonLinearSolver.h:122
HybridNonLinearSolver(FunctorType &_functor)
Definition HybridNonLinearSolver.h:48
JacobianType fjac
Definition HybridNonLinearSolver.h:92
FVectorType qtf
Definition HybridNonLinearSolver.h:91
Index njev
Definition HybridNonLinearSolver.h:95
HybridNonLinearSolverSpace::Status hybrd1(FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
Definition HybridNonLinearSolver.h:358
FVectorType diag
Definition HybridNonLinearSolver.h:91
Matrix< Scalar, Dynamic, Dynamic > UpperTriangularType
Definition HybridNonLinearSolver.h:69
#define abs(x)
Definition datatypes.h:17
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size ratio
Definition gnuplot_common_settings.hh:44
Status
Definition HybridNonLinearSolver.h:19
@ ImproperInputParameters
Definition HybridNonLinearSolver.h:21
@ TolTooSmall
Definition HybridNonLinearSolver.h:24
@ UserAsked
Definition HybridNonLinearSolver.h:27
@ Running
Definition HybridNonLinearSolver.h:20
@ NotMakingProgressIterations
Definition HybridNonLinearSolver.h:26
@ TooManyFunctionEvaluation
Definition HybridNonLinearSolver.h:23
@ NotMakingProgressJacobian
Definition HybridNonLinearSolver.h:25
@ RelativeErrorTooSmall
Definition HybridNonLinearSolver.h:22
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)
Definition fdjac1.h:6
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE float sqrt(const float &x)
Definition MathFunctions.h:177
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition MathFunctions.h:1292
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition Meta.h:66
Definition HybridNonLinearSolver.h:51
Index nb_of_subdiagonals
Definition HybridNonLinearSolver.h:62
Index nb_of_superdiagonals
Definition HybridNonLinearSolver.h:63
Parameters()
Definition HybridNonLinearSolver.h:52
Scalar factor
Definition HybridNonLinearSolver.h:59
Scalar epsfcn
Definition HybridNonLinearSolver.h:64
Scalar xtol
Definition HybridNonLinearSolver.h:61
Index maxfev
Definition HybridNonLinearSolver.h:60
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233
Definition ForwardDeclarations.h:17
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2