13 #ifndef EIGEN_HYBRIDNONLINEARSOLVER_H
14 #define EIGEN_HYBRIDNONLINEARSOLVER_H
21 namespace HybridNonLinearSolverSpace {
45 template <
typename FunctorType,
typename Scalar =
double>
120 template <
typename FunctorType,
typename Scalar>
129 parameters.maxfev = 100 * (
n + 1);
130 parameters.xtol = tol;
131 diag.setConstant(
n, 1.);
132 useExternalScaling =
true;
136 template <
typename FunctorType,
typename Scalar>
147 if (!useExternalScaling)
diag.resize(
n);
149 "When useExternalScaling is set, the caller must provide a valid 'diag'");
156 if (
n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0.)
158 if (useExternalScaling)
166 fnorm = fvec.stableNorm();
178 template <
typename FunctorType,
typename Scalar>
185 std::vector<JacobiRotation<Scalar> > v_givens(
n), w_givens(
n);
193 wa2 = fjac.colwise().blueNorm();
198 if (!useExternalScaling)
199 for (
j = 0;
j <
n; ++
j)
diag[
j] = (wa2[
j] == 0.) ? 1. : wa2[
j];
203 xnorm =
diag.cwiseProduct(
x).stableNorm();
204 delta = parameters.factor * xnorm;
221 if (!useExternalScaling)
diag =
diag.cwiseMax(wa2);
225 internal::dogleg<Scalar>(
R,
diag, qtf,
delta, wa1);
230 pnorm =
diag.cwiseProduct(wa1).stableNorm();
238 fnorm1 = wa4.stableNorm();
246 wa3 =
R.template triangularView<Upper>() * wa1 + qtf;
247 temp = wa3.stableNorm();
254 if (prered > 0.) ratio = actred / prered;
274 wa2 =
diag.cwiseProduct(
x);
276 xnorm = wa2.stableNorm();
283 if (actred >=
Scalar(.001)) nslow1 = 0;
285 if (actred >=
Scalar(.1)) nslow2 = 0;
298 if (ncfail == 2)
break;
302 wa1 =
diag.cwiseProduct(
diag.cwiseProduct(wa1) / pnorm);
303 wa2 = fjac.transpose() * wa4;
304 if (ratio >=
Scalar(1
e-4)) qtf = wa2;
305 wa2 = (wa2 - wa3) / pnorm;
308 internal::r1updt<Scalar>(
R, wa1, v_givens, w_givens, wa2, wa3, &sing);
309 internal::r1mpyq<Scalar>(
n,
n, fjac.data(), v_givens, w_givens);
310 internal::r1mpyq<Scalar>(1,
n, qtf.data(), v_givens, w_givens);
317 template <
typename FunctorType,
typename Scalar>
325 template <
typename FunctorType,
typename Scalar>
334 parameters.maxfev = 200 * (
n + 1);
335 parameters.xtol = tol;
337 diag.setConstant(
n, 1.);
338 useExternalScaling =
true;
339 return solveNumericalDiff(
x);
342 template <
typename FunctorType,
typename Scalar>
346 if (parameters.nb_of_subdiagonals < 0) parameters.nb_of_subdiagonals =
n - 1;
347 if (parameters.nb_of_superdiagonals < 0) parameters.nb_of_superdiagonals =
n - 1;
356 if (!useExternalScaling)
diag.resize(
n);
358 "When useExternalScaling is set, the caller must provide a valid 'diag'");
365 if (
n <= 0 || parameters.xtol < 0. || parameters.maxfev <= 0 || parameters.nb_of_subdiagonals < 0 ||
366 parameters.nb_of_superdiagonals < 0 || parameters.factor <= 0.)
368 if (useExternalScaling)
376 fnorm = fvec.stableNorm();
388 template <
typename FunctorType,
typename Scalar>
397 std::vector<JacobiRotation<Scalar> > v_givens(
n), w_givens(
n);
400 if (parameters.nb_of_subdiagonals < 0) parameters.nb_of_subdiagonals =
n - 1;
401 if (parameters.nb_of_superdiagonals < 0) parameters.nb_of_superdiagonals =
n - 1;
404 if (
internal::fdjac1(functor,
x, fvec, fjac, parameters.nb_of_subdiagonals, parameters.nb_of_superdiagonals,
405 parameters.epsfcn) < 0)
407 nfev += (
std::min)(parameters.nb_of_subdiagonals + parameters.nb_of_superdiagonals + 1,
n);
409 wa2 = fjac.colwise().blueNorm();
414 if (!useExternalScaling)
415 for (
j = 0;
j <
n; ++
j)
diag[
j] = (wa2[
j] == 0.) ? 1. : wa2[
j];
419 xnorm =
diag.cwiseProduct(
x).stableNorm();
420 delta = parameters.factor * xnorm;
437 if (!useExternalScaling)
diag =
diag.cwiseMax(wa2);
441 internal::dogleg<Scalar>(
R,
diag, qtf,
delta, wa1);
446 pnorm =
diag.cwiseProduct(wa1).stableNorm();
454 fnorm1 = wa4.stableNorm();
462 wa3 =
R.template triangularView<Upper>() * wa1 + qtf;
463 temp = wa3.stableNorm();
470 if (prered > 0.) ratio = actred / prered;
490 wa2 =
diag.cwiseProduct(
x);
492 xnorm = wa2.stableNorm();
499 if (actred >=
Scalar(.001)) nslow1 = 0;
501 if (actred >=
Scalar(.1)) nslow2 = 0;
514 if (ncfail == 2)
break;
518 wa1 =
diag.cwiseProduct(
diag.cwiseProduct(wa1) / pnorm);
519 wa2 = fjac.transpose() * wa4;
520 if (ratio >=
Scalar(1
e-4)) qtf = wa2;
521 wa2 = (wa2 - wa3) / pnorm;
524 internal::r1updt<Scalar>(
R, wa1, v_givens, w_givens, wa2, wa3, &sing);
525 internal::r1mpyq<Scalar>(
n,
n, fjac.data(), v_givens, w_givens);
526 internal::r1mpyq<Scalar>(1,
n, qtf.data(), v_givens, w_givens);
533 template <
typename FunctorType,
typename Scalar>
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define eigen_assert(x)
Definition: Macros.h:910
@ R
Definition: StatisticsVector.h:21
SCALAR Scalar
Definition: bench_gemm.cpp:45
Householder QR decomposition of a matrix.
Definition: HouseholderQR.h:59
HouseholderSequenceType householderQ() const
Definition: HouseholderQR.h:160
const MatrixType & matrixQR() const
Definition: HouseholderQR.h:168
TransposeReturnType transpose() const
Transpose of the Householder sequence.
Definition: HouseholderSequence.h:216
Finds a zero of a system of n nonlinear functions in n variables by a modification of the Powell hybr...
Definition: HybridNonLinearSolver.h:46
HybridNonLinearSolverSpace::Status solveNumericalDiff(FVectorType &x)
Definition: HybridNonLinearSolver.h:534
HybridNonLinearSolverSpace::Status solveNumericalDiffInit(FVectorType &x)
Definition: HybridNonLinearSolver.h:343
void resetParameters(void)
Definition: HybridNonLinearSolver.h:90
FVectorType fvec
Definition: HybridNonLinearSolver.h:92
FVectorType wa2
Definition: HybridNonLinearSolver.h:115
Scalar delta
Definition: HybridNonLinearSolver.h:107
HybridNonLinearSolverSpace::Status solveInit(FVectorType &x)
Definition: HybridNonLinearSolver.h:137
Scalar xnorm
Definition: HybridNonLinearSolver.h:111
Index iter
Definition: HybridNonLinearSolver.h:97
Scalar fnorm1
Definition: HybridNonLinearSolver.h:111
HybridNonLinearSolver & operator=(const HybridNonLinearSolver &)
bool useExternalScaling
Definition: HybridNonLinearSolver.h:99
Matrix< Scalar, Dynamic, Dynamic > JacobianType
Definition: HybridNonLinearSolver.h:72
Index ncsuc
Definition: HybridNonLinearSolver.h:109
UpperTriangularType R
Definition: HybridNonLinearSolver.h:94
Index nslow1
Definition: HybridNonLinearSolver.h:112
Scalar sum
Definition: HybridNonLinearSolver.h:104
HybridNonLinearSolverSpace::Status solve(FVectorType &x)
Definition: HybridNonLinearSolver.h:318
Parameters parameters
Definition: HybridNonLinearSolver.h:91
FVectorType wa4
Definition: HybridNonLinearSolver.h:115
Scalar prered
Definition: HybridNonLinearSolver.h:114
HybridNonLinearSolverSpace::Status solveOneStep(FVectorType &x)
Definition: HybridNonLinearSolver.h:179
Index nfev
Definition: HybridNonLinearSolver.h:95
HybridNonLinearSolverSpace::Status solveNumericalDiffOneStep(FVectorType &x)
Definition: HybridNonLinearSolver.h:389
Index nslow2
Definition: HybridNonLinearSolver.h:112
Scalar fnorm
Definition: HybridNonLinearSolver.h:98
FVectorType wa1
Definition: HybridNonLinearSolver.h:115
Matrix< Scalar, Dynamic, Dynamic > UpperTriangularType
Definition: HybridNonLinearSolver.h:74
Scalar temp
Definition: HybridNonLinearSolver.h:106
DenseIndex Index
Definition: HybridNonLinearSolver.h:48
Index n
Definition: HybridNonLinearSolver.h:103
Scalar ratio
Definition: HybridNonLinearSolver.h:110
FVectorType wa3
Definition: HybridNonLinearSolver.h:115
bool sing
Definition: HybridNonLinearSolver.h:105
Index ncfail
Definition: HybridNonLinearSolver.h:113
bool jeval
Definition: HybridNonLinearSolver.h:108
HybridNonLinearSolverSpace::Status hybrj1(FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
Definition: HybridNonLinearSolver.h:121
HybridNonLinearSolver(FunctorType &_functor)
Definition: HybridNonLinearSolver.h:50
JacobianType fjac
Definition: HybridNonLinearSolver.h:93
Scalar actred
Definition: HybridNonLinearSolver.h:114
FVectorType qtf
Definition: HybridNonLinearSolver.h:92
Index njev
Definition: HybridNonLinearSolver.h:96
Matrix< Scalar, Dynamic, 1 > FVectorType
Definition: HybridNonLinearSolver.h:71
HybridNonLinearSolverSpace::Status hybrd1(FVectorType &x, const Scalar tol=numext::sqrt(NumTraits< Scalar >::epsilon()))
Definition: HybridNonLinearSolver.h:326
FunctorType & functor
Definition: HybridNonLinearSolver.h:102
FVectorType diag
Definition: HybridNonLinearSolver.h:92
Scalar pnorm
Definition: HybridNonLinearSolver.h:111
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
const char const char const char * diag
Definition: level2_impl.h:86
Status
Definition: HybridNonLinearSolver.h:22
@ ImproperInputParameters
Definition: HybridNonLinearSolver.h:24
@ TolTooSmall
Definition: HybridNonLinearSolver.h:27
@ UserAsked
Definition: HybridNonLinearSolver.h:30
@ Running
Definition: HybridNonLinearSolver.h:23
@ NotMakingProgressIterations
Definition: HybridNonLinearSolver.h:29
@ TooManyFunctionEvaluation
Definition: HybridNonLinearSolver.h:26
@ NotMakingProgressJacobian
Definition: HybridNonLinearSolver.h:28
@ RelativeErrorTooSmall
Definition: HybridNonLinearSolver.h:25
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:9
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE float sqrt(const float &x)
Definition: arch/SSE/MathFunctions.h:69
EIGEN_DEVICE_FUNC bool abs2(bool x)
Definition: MathFunctions.h:1102
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:75
int delta
Definition: MultiOpt.py:96
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
list x
Definition: plotDoE.py:28
Update the problem specs before solve
Definition: steady_axisym_advection_diffusion.cc:353
Definition: HybridNonLinearSolver.h:56
Index nb_of_subdiagonals
Definition: HybridNonLinearSolver.h:67
Index nb_of_superdiagonals
Definition: HybridNonLinearSolver.h:68
Parameters()
Definition: HybridNonLinearSolver.h:57
Scalar factor
Definition: HybridNonLinearSolver.h:64
Scalar epsfcn
Definition: HybridNonLinearSolver.h:69
Scalar xtol
Definition: HybridNonLinearSolver.h:66
Index maxfev
Definition: HybridNonLinearSolver.h:65
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2