mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-12 11:49:02 +08:00
merge both c methods lmdif/lmdif1 into one class
LevenbergMarquardtNumericalDiff with two methods.
This commit is contained in:
parent
a736378331
commit
3f1b81e129
@ -54,7 +54,6 @@ namespace Eigen {
|
|||||||
#include "src/NonLinear/lmdif.h"
|
#include "src/NonLinear/lmdif.h"
|
||||||
#include "src/NonLinear/hybrj.h"
|
#include "src/NonLinear/hybrj.h"
|
||||||
#include "src/NonLinear/lmstr1.h"
|
#include "src/NonLinear/lmstr1.h"
|
||||||
#include "src/NonLinear/lmdif1.h"
|
|
||||||
#include "src/NonLinear/chkder.h"
|
#include "src/NonLinear/chkder.h"
|
||||||
|
|
||||||
//@}
|
//@}
|
||||||
|
@ -1,7 +1,19 @@
|
|||||||
|
|
||||||
|
|
||||||
template<typename FunctorType, typename Scalar>
|
template<typename FunctorType, typename Scalar>
|
||||||
int ei_lmdif(
|
class LevenbergMarquardtNumericalDiff
|
||||||
const FunctorType &Functor,
|
{
|
||||||
|
public:
|
||||||
|
LevenbergMarquardtNumericalDiff(const FunctorType &_functor)
|
||||||
|
: functor(_functor) {}
|
||||||
|
|
||||||
|
int minimize(
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &fvec,
|
||||||
|
const Scalar tol = ei_sqrt(epsilon<Scalar>())
|
||||||
|
);
|
||||||
|
|
||||||
|
int minimize(
|
||||||
Matrix< Scalar, Dynamic, 1 > &x,
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
Matrix< Scalar, Dynamic, 1 > &fvec,
|
Matrix< Scalar, Dynamic, 1 > &fvec,
|
||||||
int &nfev,
|
int &nfev,
|
||||||
@ -17,6 +29,61 @@ int ei_lmdif(
|
|||||||
Scalar gtol = Scalar(0.),
|
Scalar gtol = Scalar(0.),
|
||||||
Scalar epsfcn = Scalar(0.),
|
Scalar epsfcn = Scalar(0.),
|
||||||
int nprint=0
|
int nprint=0
|
||||||
|
);
|
||||||
|
|
||||||
|
private:
|
||||||
|
const FunctorType &functor;
|
||||||
|
};
|
||||||
|
|
||||||
|
|
||||||
|
template<typename FunctorType, typename Scalar>
|
||||||
|
int LevenbergMarquardtNumericalDiff<FunctorType,Scalar>::minimize(
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &fvec,
|
||||||
|
Scalar tol
|
||||||
|
)
|
||||||
|
{
|
||||||
|
const int n = x.size(), m=fvec.size();
|
||||||
|
int info, nfev=0;
|
||||||
|
Matrix< Scalar, Dynamic, Dynamic > fjac(m, n);
|
||||||
|
Matrix< Scalar, Dynamic, 1> diag, qtf;
|
||||||
|
VectorXi ipvt;
|
||||||
|
|
||||||
|
/* check the input parameters for errors. */
|
||||||
|
if (n <= 0 || m < n || tol < 0.) {
|
||||||
|
printf("ei_lmder1 bad args : m,n,tol,...");
|
||||||
|
return 0;
|
||||||
|
}
|
||||||
|
|
||||||
|
info = minimize(
|
||||||
|
x, fvec,
|
||||||
|
nfev,
|
||||||
|
fjac, ipvt, qtf, diag,
|
||||||
|
1,
|
||||||
|
100.,
|
||||||
|
(n+1)*200,
|
||||||
|
tol, tol, Scalar(0.), Scalar(0.)
|
||||||
|
);
|
||||||
|
return (info==8)?4:info;
|
||||||
|
}
|
||||||
|
|
||||||
|
template<typename FunctorType, typename Scalar>
|
||||||
|
int LevenbergMarquardtNumericalDiff<FunctorType,Scalar>::minimize(
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &fvec,
|
||||||
|
int &nfev,
|
||||||
|
Matrix< Scalar, Dynamic, Dynamic > &fjac,
|
||||||
|
VectorXi &ipvt,
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &qtf,
|
||||||
|
Matrix< Scalar, Dynamic, 1 > &diag,
|
||||||
|
int mode,
|
||||||
|
Scalar factor,
|
||||||
|
int maxfev,
|
||||||
|
Scalar ftol,
|
||||||
|
Scalar xtol,
|
||||||
|
Scalar gtol,
|
||||||
|
Scalar epsfcn,
|
||||||
|
int nprint
|
||||||
)
|
)
|
||||||
{
|
{
|
||||||
const int m = fvec.size(), n = x.size();
|
const int m = fvec.size(), n = x.size();
|
||||||
@ -55,7 +122,7 @@ int ei_lmdif(
|
|||||||
/* evaluate the function at the starting point */
|
/* evaluate the function at the starting point */
|
||||||
/* and calculate its norm. */
|
/* and calculate its norm. */
|
||||||
|
|
||||||
iflag = Functor.f(x, fvec);
|
iflag = functor.f(x, fvec);
|
||||||
nfev = 1;
|
nfev = 1;
|
||||||
if (iflag < 0)
|
if (iflag < 0)
|
||||||
goto algo_end;
|
goto algo_end;
|
||||||
@ -72,17 +139,17 @@ int ei_lmdif(
|
|||||||
|
|
||||||
/* calculate the jacobian matrix. */
|
/* calculate the jacobian matrix. */
|
||||||
|
|
||||||
iflag = ei_fdjac2(Functor, x, fvec, fjac, epsfcn);
|
iflag = ei_fdjac2(functor, x, fvec, fjac, epsfcn);
|
||||||
nfev += n;
|
nfev += n;
|
||||||
if (iflag < 0)
|
if (iflag < 0)
|
||||||
break;
|
break;
|
||||||
|
|
||||||
/* if requested, call Functor.f to enable printing of iterates. */
|
/* if requested, call functor.f to enable printing of iterates. */
|
||||||
|
|
||||||
if (nprint > 0) {
|
if (nprint > 0) {
|
||||||
iflag = 0;
|
iflag = 0;
|
||||||
if ((iter - 1) % nprint == 0)
|
if ((iter - 1) % nprint == 0)
|
||||||
iflag = Functor.debug(x, fvec);
|
iflag = functor.debug(x, fvec);
|
||||||
if (iflag < 0)
|
if (iflag < 0)
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
@ -178,7 +245,7 @@ int ei_lmdif(
|
|||||||
|
|
||||||
/* evaluate the function at x + p and calculate its norm. */
|
/* evaluate the function at x + p and calculate its norm. */
|
||||||
|
|
||||||
iflag = Functor.f(wa2, wa4);
|
iflag = functor.f(wa2, wa4);
|
||||||
++nfev;
|
++nfev;
|
||||||
if (iflag < 0)
|
if (iflag < 0)
|
||||||
goto algo_end;
|
goto algo_end;
|
||||||
@ -275,7 +342,7 @@ algo_end:
|
|||||||
if (iflag < 0)
|
if (iflag < 0)
|
||||||
info = iflag;
|
info = iflag;
|
||||||
if (nprint > 0)
|
if (nprint > 0)
|
||||||
iflag = Functor.debug(x, fvec);
|
iflag = functor.debug(x, fvec);
|
||||||
return info;
|
return info;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -1,34 +0,0 @@
|
|||||||
|
|
||||||
template<typename FunctorType, typename Scalar>
|
|
||||||
int ei_lmdif1(
|
|
||||||
const FunctorType &Functor,
|
|
||||||
Matrix< Scalar, Dynamic, 1 > &x,
|
|
||||||
Matrix< Scalar, Dynamic, 1 > &fvec,
|
|
||||||
Scalar tol = ei_sqrt(epsilon<Scalar>())
|
|
||||||
)
|
|
||||||
{
|
|
||||||
const int n = x.size(), m=fvec.size();
|
|
||||||
int info, nfev=0;
|
|
||||||
Matrix< Scalar, Dynamic, Dynamic > fjac(m, n);
|
|
||||||
Matrix< Scalar, Dynamic, 1> diag, qtf;
|
|
||||||
VectorXi ipvt;
|
|
||||||
|
|
||||||
/* check the input parameters for errors. */
|
|
||||||
if (n <= 0 || m < n || tol < 0.) {
|
|
||||||
printf("ei_lmder1 bad args : m,n,tol,...");
|
|
||||||
return 0;
|
|
||||||
}
|
|
||||||
|
|
||||||
info = ei_lmdif(
|
|
||||||
Functor,
|
|
||||||
x, fvec,
|
|
||||||
nfev,
|
|
||||||
fjac, ipvt, qtf, diag,
|
|
||||||
1,
|
|
||||||
100.,
|
|
||||||
(n+1)*200,
|
|
||||||
tol, tol, Scalar(0.), Scalar(0.)
|
|
||||||
);
|
|
||||||
return (info==8)?4:info;
|
|
||||||
}
|
|
||||||
|
|
@ -541,7 +541,9 @@ void testLmdif1()
|
|||||||
x.setConstant(n, 1.);
|
x.setConstant(n, 1.);
|
||||||
|
|
||||||
// do the computation
|
// do the computation
|
||||||
info = ei_lmdif1(lmdif_functor(), x, fvec);
|
lmdif_functor functor;
|
||||||
|
LevenbergMarquardtNumericalDiff<lmdif_functor,double> lm(functor);
|
||||||
|
info = lm.minimize(x, fvec);
|
||||||
|
|
||||||
// check return value
|
// check return value
|
||||||
VERIFY( 1 == info);
|
VERIFY( 1 == info);
|
||||||
@ -569,7 +571,9 @@ void testLmdif()
|
|||||||
x.setConstant(n, 1.);
|
x.setConstant(n, 1.);
|
||||||
|
|
||||||
// do the computation
|
// do the computation
|
||||||
info = ei_lmdif(lmdif_functor(), x, fvec, nfev, fjac, ipvt, qtf, diag);
|
lmdif_functor functor;
|
||||||
|
LevenbergMarquardtNumericalDiff<lmdif_functor,double> lm(functor);
|
||||||
|
info = lm.minimize(x, fvec, nfev, fjac, ipvt, qtf, diag);
|
||||||
|
|
||||||
// check return values
|
// check return values
|
||||||
VERIFY( 1 == info);
|
VERIFY( 1 == info);
|
||||||
|
Loading…
x
Reference in New Issue
Block a user