13 #ifndef EIGEN_LEVENBERGMARQUARDT__H
14 #define EIGEN_LEVENBERGMARQUARDT__H
21 namespace LevenbergMarquardtSpace {
46 template <
typename FunctorType,
typename Scalar =
double>
47 class LevenbergMarquardt {
125 template <
typename FunctorType,
typename Scalar>
128 m = functor.values();
134 parameters.ftol = tol;
135 parameters.xtol = tol;
136 parameters.maxfev = 100 * (
n + 1);
141 template <
typename FunctorType,
typename Scalar>
146 status = minimizeOneStep(
x);
151 template <
typename FunctorType,
typename Scalar>
154 m = functor.values();
162 if (!useExternalScaling)
diag.resize(
n);
164 "When useExternalScaling is set, the caller must provide a valid 'diag'");
172 if (
n <= 0 ||
m <
n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. ||
173 parameters.maxfev <= 0 || parameters.factor <= 0.)
176 if (useExternalScaling)
184 fnorm = fvec.stableNorm();
193 template <
typename FunctorType,
typename Scalar>
201 Index df_ret = functor.df(
x, fjac);
210 wa2 = fjac.colwise().blueNorm();
211 ColPivHouseholderQR<JacobianType> qrfac(fjac);
212 fjac = qrfac.matrixQR();
213 permutation = qrfac.colsPermutation();
218 if (!useExternalScaling)
223 xnorm =
diag.cwiseProduct(
x).stableNorm();
224 delta = parameters.factor * xnorm;
231 wa4.applyOnTheLeft(qrfac.householderQ().adjoint());
238 if (wa2[permutation.indices()[
j]] != 0.)
240 abs(fjac.col(
j).head(
j + 1).dot(qtf.head(
j + 1) / fnorm) / wa2[permutation.indices()[
j]]));
246 if (!useExternalScaling)
diag =
diag.cwiseMax(wa2);
250 internal::lmpar2<Scalar>(qrfac,
diag, qtf,
delta,
par, wa1);
255 pnorm =
diag.cwiseProduct(wa1).stableNorm();
263 fnorm1 = wa4.stableNorm();
271 wa3 = fjac.template triangularView<Upper>() * (qrfac.colsPermutation().inverse() * wa1);
274 prered = temp1 + temp2 /
Scalar(.5);
275 dirder = -(temp1 + temp2);
280 if (prered != 0.) ratio = actred / prered;
283 if (ratio <=
Scalar(.25)) {
284 if (actred >= 0.) temp =
Scalar(.5);
285 if (actred < 0.) temp =
Scalar(.5) * dirder / (dirder +
Scalar(.5) * actred);
290 }
else if (!(
par != 0. && ratio <
Scalar(.75))) {
299 wa2 =
diag.cwiseProduct(
x);
301 xnorm = wa2.stableNorm();
307 if (
abs(actred) <= parameters.ftol && prered <= parameters.ftol &&
Scalar(.5) * ratio <= 1. &&
308 delta <= parameters.xtol * xnorm)
310 if (
abs(actred) <= parameters.ftol && prered <= parameters.ftol &&
Scalar(.5) * ratio <= 1.)
322 }
while (ratio <
Scalar(1
e-4));
327 template <
typename FunctorType,
typename Scalar>
330 m = functor.values();
336 parameters.ftol = tol;
337 parameters.xtol = tol;
338 parameters.maxfev = 100 * (
n + 1);
340 return minimizeOptimumStorage(
x);
343 template <
typename FunctorType,
typename Scalar>
346 m = functor.values();
359 if (!useExternalScaling)
diag.resize(
n);
361 "When useExternalScaling is set, the caller must provide a valid 'diag'");
369 if (
n <= 0 ||
m <
n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. ||
370 parameters.maxfev <= 0 || parameters.factor <= 0.)
373 if (useExternalScaling)
381 fnorm = fvec.stableNorm();
390 template <
typename FunctorType,
typename Scalar>
407 for (
i = 0;
i <
m; ++
i) {
409 internal::rwupdt<Scalar>(fjac, wa3, qtf, fvec[
i]);
417 for (
j = 0;
j <
n; ++
j) {
418 if (fjac(
j,
j) == 0.) sing =
true;
419 wa2[
j] = fjac.col(
j).head(
j).stableNorm();
421 permutation.setIdentity(
n);
423 wa2 = fjac.colwise().blueNorm();
428 wa1 = fjac.diagonal();
429 fjac.diagonal() = qrfac.
hCoeffs();
432 for (
Index ii = 0; ii < fjac.cols(); ii++)
433 fjac.col(ii).segment(ii + 1, fjac.rows() - ii - 1) *= fjac(ii, ii);
435 for (
j = 0;
j <
n; ++
j) {
436 if (fjac(
j,
j) != 0.) {
438 for (
i =
j;
i <
n; ++
i) sum += fjac(
i,
j) * qtf[
i];
439 temp = -sum / fjac(
j,
j);
440 for (
i =
j;
i <
n; ++
i) qtf[
i] += fjac(
i,
j) * temp;
449 if (!useExternalScaling)
450 for (
j = 0;
j <
n; ++
j)
diag[
j] = (wa2[
j] == 0.) ? 1. : wa2[
j];
454 xnorm =
diag.cwiseProduct(
x).stableNorm();
455 delta = parameters.factor * xnorm;
462 for (
j = 0;
j <
n; ++
j)
463 if (wa2[permutation.indices()[
j]] != 0.)
465 abs(fjac.col(
j).head(
j + 1).dot(qtf.head(
j + 1) / fnorm) / wa2[permutation.indices()[
j]]));
471 if (!useExternalScaling)
diag =
diag.cwiseMax(wa2);
475 internal::lmpar<Scalar>(fjac, permutation.indices(),
diag, qtf,
delta,
par, wa1);
480 pnorm =
diag.cwiseProduct(wa1).stableNorm();
488 fnorm1 = wa4.stableNorm();
496 wa3 = fjac.topLeftCorner(
n,
n).template triangularView<Upper>() * (permutation.inverse() * wa1);
499 prered = temp1 + temp2 /
Scalar(.5);
500 dirder = -(temp1 + temp2);
505 if (prered != 0.) ratio = actred / prered;
508 if (ratio <=
Scalar(.25)) {
509 if (actred >= 0.) temp =
Scalar(.5);
510 if (actred < 0.) temp =
Scalar(.5) * dirder / (dirder +
Scalar(.5) * actred);
515 }
else if (!(
par != 0. && ratio <
Scalar(.75))) {
524 wa2 =
diag.cwiseProduct(
x);
526 xnorm = wa2.stableNorm();
532 if (
abs(actred) <= parameters.ftol && prered <= parameters.ftol &&
Scalar(.5) * ratio <= 1. &&
533 delta <= parameters.xtol * xnorm)
535 if (
abs(actred) <= parameters.ftol && prered <= parameters.ftol &&
Scalar(.5) * ratio <= 1.)
547 }
while (ratio <
Scalar(1
e-4));
552 template <
typename FunctorType,
typename Scalar>
557 status = minimizeOptimumStorageOneStep(
x);
562 template <
typename FunctorType,
typename Scalar>
566 Index m = functor.values();
574 lm.parameters.ftol = tol;
575 lm.parameters.xtol = tol;
576 lm.parameters.maxfev = 200 * (
n + 1);
579 if (nfev) *nfev = lm.nfev;
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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
SCALAR Scalar
Definition: bench_gemm.cpp:45
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
Definition: ColPivHouseholderQR.h:54
const HCoeffsType & hCoeffs() const
Definition: ColPivHouseholderQR.h:333
const PermutationType & colsPermutation() const
Definition: ColPivHouseholderQR.h:192
const MatrixType & matrixQR() const
Definition: ColPivHouseholderQR.h:169
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
Definition: LevenbergMarquardt/LevenbergMarquardt.h:102
static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
Definition: LevenbergMarquardt/LevenbergMarquardt.h:340
JacobianType::Scalar Scalar
Definition: LevenbergMarquardt/LevenbergMarquardt.h:107
FVectorType diag
Definition: NonLinearOptimization/LevenbergMarquardt.h:99
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol=sqrt_epsilon())
FunctorType & functor
Definition: NonLinearOptimization/LevenbergMarquardt.h:111
bool useExternalScaling
Definition: NonLinearOptimization/LevenbergMarquardt.h:106
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)
static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=sqrt_epsilon())
Index m
Definition: LevenbergMarquardt/LevenbergMarquardt.h:239
LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep(FVectorType &x)
Definition: NonLinearOptimization/LevenbergMarquardt.h:391
Scalar temp2
Definition: NonLinearOptimization/LevenbergMarquardt.h:117
FVectorType wa4
Definition: NonLinearOptimization/LevenbergMarquardt.h:114
Scalar pnorm
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
JacobianType fjac
Definition: NonLinearOptimization/LevenbergMarquardt.h:100
FVectorType fvec
Definition: NonLinearOptimization/LevenbergMarquardt.h:99
FunctorType_ FunctorType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:104
LevenbergMarquardtSpace::Status minimize(FVectorType &x)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:261
Index iter
Definition: NonLinearOptimization/LevenbergMarquardt.h:104
Scalar delta
Definition: NonLinearOptimization/LevenbergMarquardt.h:118
LevenbergMarquardtSpace::Status minimize(FVectorType &x)
LevenbergMarquardtSpace::Status lmstr1(FVectorType &x, const Scalar tol=sqrt_epsilon())
Definition: NonLinearOptimization/LevenbergMarquardt.h:328
Scalar lm_param(void)
Definition: NonLinearOptimization/LevenbergMarquardt.h:108
Scalar dirder
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
Index njev
Definition: NonLinearOptimization/LevenbergMarquardt.h:103
Scalar temp
Definition: NonLinearOptimization/LevenbergMarquardt.h:117
Matrix< Scalar, Dynamic, 1 > FVectorType
Definition: NonLinearOptimization/LevenbergMarquardt.h:78
Matrix< Scalar, Dynamic, Dynamic > JacobianType
Definition: NonLinearOptimization/LevenbergMarquardt.h:79
FVectorType wa1
Definition: NonLinearOptimization/LevenbergMarquardt.h:114
LevenbergMarquardtSpace::Status lmder1(FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
Definition: LevenbergMarquardt/LevenbergMarquardt.h:324
Index nfev
Definition: NonLinearOptimization/LevenbergMarquardt.h:102
LevenbergMarquardtSpace::Status minimizeOptimumStorage(FVectorType &x)
Definition: NonLinearOptimization/LevenbergMarquardt.h:553
Index n
Definition: LevenbergMarquardt/LevenbergMarquardt.h:238
PermutationMatrix< Dynamic, Dynamic > permutation
Definition: NonLinearOptimization/LevenbergMarquardt.h:101
Scalar xnorm
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
Scalar prered
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
Scalar fnorm1
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
LevenbergMarquardt & operator=(const LevenbergMarquardt &)
void resetParameters(void)
Definition: NonLinearOptimization/LevenbergMarquardt.h:96
DenseIndex Index
Definition: NonLinearOptimization/LevenbergMarquardt.h:60
Scalar temp1
Definition: NonLinearOptimization/LevenbergMarquardt.h:117
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:276
Scalar ratio
Definition: NonLinearOptimization/LevenbergMarquardt.h:119
Parameters parameters
Definition: NonLinearOptimization/LevenbergMarquardt.h:98
FVectorType wa2
Definition: NonLinearOptimization/LevenbergMarquardt.h:114
Scalar gnorm
Definition: NonLinearOptimization/LevenbergMarquardt.h:105
LevenbergMarquardtSpace::Status minimizeOptimumStorageInit(FVectorType &x)
Definition: NonLinearOptimization/LevenbergMarquardt.h:344
LevenbergMarquardt(FunctorType &_functor)
Definition: NonLinearOptimization/LevenbergMarquardt.h:54
FunctorType::JacobianType JacobianType
Definition: LevenbergMarquardt/LevenbergMarquardt.h:106
FVectorType qtf
Definition: NonLinearOptimization/LevenbergMarquardt.h:99
Scalar sum
Definition: NonLinearOptimization/LevenbergMarquardt.h:116
Scalar fnorm
Definition: NonLinearOptimization/LevenbergMarquardt.h:105
static Scalar sqrt_epsilon()
Definition: NonLinearOptimization/LevenbergMarquardt.h:48
Scalar actred
Definition: NonLinearOptimization/LevenbergMarquardt.h:120
Scalar par
Definition: NonLinearOptimization/LevenbergMarquardt.h:116
FVectorType wa3
Definition: NonLinearOptimization/LevenbergMarquardt.h:114
Definition: NumericalDiff.h:35
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
int * m
Definition: level2_cplx_impl.h:294
int info
Definition: level2_cplx_impl.h:39
const char const char const char * diag
Definition: level2_impl.h:86
Status
Definition: LevenbergMarquardt/LevenbergMarquardt.h:27
@ RelativeReductionTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:31
@ GtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:38
@ NotStarted
Definition: LevenbergMarquardt/LevenbergMarquardt.h:28
@ UserAsked
Definition: LevenbergMarquardt/LevenbergMarquardt.h:39
@ Running
Definition: LevenbergMarquardt/LevenbergMarquardt.h:29
@ FtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:36
@ TooManyFunctionEvaluation
Definition: LevenbergMarquardt/LevenbergMarquardt.h:35
@ XtolTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:37
@ CosinusTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:34
@ RelativeErrorTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:32
@ ImproperInputParameters
Definition: LevenbergMarquardt/LevenbergMarquardt.h:30
@ RelativeErrorAndReductionTooSmall
Definition: LevenbergMarquardt/LevenbergMarquardt.h:33
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 Index
The Index type as used for the API.
Definition: Meta.h:83
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:75
int delta
Definition: MultiOpt.py:96
string par
Definition: calibrate.py:135
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
list x
Definition: plotDoE.py:28
Definition: NonLinearOptimization/LevenbergMarquardt.h:62
Scalar factor
Definition: NonLinearOptimization/LevenbergMarquardt.h:70
Parameters()
Definition: NonLinearOptimization/LevenbergMarquardt.h:63
Index maxfev
Definition: NonLinearOptimization/LevenbergMarquardt.h:71
Scalar ftol
Definition: NonLinearOptimization/LevenbergMarquardt.h:72
Scalar xtol
Definition: NonLinearOptimization/LevenbergMarquardt.h:73
Scalar gtol
Definition: NonLinearOptimization/LevenbergMarquardt.h:74
Scalar epsfcn
Definition: NonLinearOptimization/LevenbergMarquardt.h:75
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