mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-12 19:59:05 +08:00
porting lmstr1 to eigen
This commit is contained in:
parent
980c40f72c
commit
de7d14b2b3
@ -130,7 +130,6 @@ Scalar ei_enorm ( int n, const Scalar *x ){
|
|||||||
#include "src/NonLinear/lmder.h"
|
#include "src/NonLinear/lmder.h"
|
||||||
#include "src/NonLinear/hybrd1.h"
|
#include "src/NonLinear/hybrd1.h"
|
||||||
#include "src/NonLinear/hybrd.h"
|
#include "src/NonLinear/hybrd.h"
|
||||||
#include "src/NonLinear/lmstr1.h"
|
|
||||||
#include "src/NonLinear/lmstr.h"
|
#include "src/NonLinear/lmstr.h"
|
||||||
#include "src/NonLinear/lmdif1.h"
|
#include "src/NonLinear/lmdif1.h"
|
||||||
#include "src/NonLinear/lmdif.h"
|
#include "src/NonLinear/lmdif.h"
|
||||||
@ -139,6 +138,7 @@ Scalar ei_enorm ( int n, const Scalar *x ){
|
|||||||
#include "src/NonLinear/chkder.h"
|
#include "src/NonLinear/chkder.h"
|
||||||
#include "src/NonLinear/MathFunctions.h"
|
#include "src/NonLinear/MathFunctions.h"
|
||||||
#include "src/NonLinear/lmder1.h"
|
#include "src/NonLinear/lmder1.h"
|
||||||
|
#include "src/NonLinear/lmstr1.h"
|
||||||
|
|
||||||
//@}
|
//@}
|
||||||
|
|
||||||
|
@ -149,30 +149,6 @@ int ei_hybrj(
|
|||||||
);
|
);
|
||||||
}
|
}
|
||||||
|
|
||||||
template<typename Functor, typename Scalar>
|
|
||||||
int ei_lmstr1(
|
|
||||||
Matrix< Scalar, Dynamic, 1 > &x,
|
|
||||||
Matrix< Scalar, Dynamic, 1 > &fvec,
|
|
||||||
VectorXi &ipvt,
|
|
||||||
Scalar tol = ei_sqrt(epsilon<Scalar>())
|
|
||||||
)
|
|
||||||
{
|
|
||||||
int lwa = 5*x.size()+fvec.size();
|
|
||||||
int ldfjac = fvec.size();
|
|
||||||
Matrix< Scalar, Dynamic, 1 > wa(lwa);
|
|
||||||
Matrix< Scalar, Dynamic, Dynamic > fjac(ldfjac, x.size());
|
|
||||||
|
|
||||||
ipvt.resize(x.size());
|
|
||||||
return lmstr1_template<Scalar>(
|
|
||||||
Functor::f, 0,
|
|
||||||
fvec.size(), x.size(), x.data(), fvec.data(),
|
|
||||||
fjac.data() , ldfjac,
|
|
||||||
tol,
|
|
||||||
ipvt.data(),
|
|
||||||
wa.data(), lwa
|
|
||||||
);
|
|
||||||
}
|
|
||||||
|
|
||||||
template<typename Functor, typename Scalar>
|
template<typename Functor, typename Scalar>
|
||||||
int ei_lmstr(
|
int ei_lmstr(
|
||||||
Matrix< Scalar, Dynamic, 1 > &x,
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
|
@ -11,15 +11,14 @@ int ei_lmder1(
|
|||||||
int info, nfev, njev;
|
int info, nfev, njev;
|
||||||
Matrix< Scalar, Dynamic, Dynamic > fjac(m, n);
|
Matrix< Scalar, Dynamic, Dynamic > fjac(m, n);
|
||||||
Matrix< Scalar, Dynamic, 1> diag;
|
Matrix< Scalar, Dynamic, 1> diag;
|
||||||
ipvt.resize(n);
|
|
||||||
|
|
||||||
|
/* check the input parameters for errors. */
|
||||||
/* check the input parameters for errors. */
|
|
||||||
if (n <= 0 || m < n || tol < 0.) {
|
if (n <= 0 || m < n || tol < 0.) {
|
||||||
printf("ei_lmder1 bad args : m,n,tol,...");
|
printf("ei_lmder1 bad args : m,n,tol,...");
|
||||||
return 0;
|
return 0;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
ipvt.resize(n);
|
||||||
info = ei_lmder<Functor,Scalar>(
|
info = ei_lmder<Functor,Scalar>(
|
||||||
x, fvec,
|
x, fvec,
|
||||||
nfev, njev,
|
nfev, njev,
|
||||||
|
@ -1,61 +1,33 @@
|
|||||||
|
|
||||||
template<typename Scalar>
|
template<typename Functor, typename Scalar>
|
||||||
int lmstr1_template(minpack_funcderstr_mn fcn, void *p, int m, int n, Scalar *x,
|
int ei_lmstr1(
|
||||||
Scalar *fvec, Scalar *fjac, int ldfjac, Scalar tol,
|
Matrix< Scalar, Dynamic, 1 > &x,
|
||||||
int *ipvt, Scalar *wa, int lwa)
|
Matrix< Scalar, Dynamic, 1 > &fvec,
|
||||||
|
VectorXi &ipvt,
|
||||||
|
Scalar tol = ei_sqrt(epsilon<Scalar>())
|
||||||
|
)
|
||||||
{
|
{
|
||||||
/* Initialized data */
|
const int n = x.size(), m=fvec.size();
|
||||||
|
int info, nfev, njev;
|
||||||
|
Matrix< Scalar, Dynamic, Dynamic > fjac(m, n);
|
||||||
|
Matrix< Scalar, Dynamic, 1> diag;
|
||||||
|
|
||||||
const Scalar factor = 100.;
|
/* check the input parameters for errors. */
|
||||||
|
if (n <= 0 || m < n || tol < 0.) {
|
||||||
/* System generated locals */
|
printf("ei_lmstr1 bad args : m,n,tol,...");
|
||||||
int fjac_dim1, fjac_offset;
|
return 0;
|
||||||
|
|
||||||
/* Local variables */
|
|
||||||
int mode, nfev, njev;
|
|
||||||
Scalar ftol, gtol, xtol;
|
|
||||||
int maxfev, nprint;
|
|
||||||
int info;
|
|
||||||
|
|
||||||
/* Parameter adjustments */
|
|
||||||
--fvec;
|
|
||||||
--ipvt;
|
|
||||||
--x;
|
|
||||||
fjac_dim1 = ldfjac;
|
|
||||||
fjac_offset = 1 + fjac_dim1 * 1;
|
|
||||||
fjac -= fjac_offset;
|
|
||||||
--wa;
|
|
||||||
|
|
||||||
/* Function Body */
|
|
||||||
info = 0;
|
|
||||||
|
|
||||||
/* check the input parameters for errors. */
|
|
||||||
|
|
||||||
if (n <= 0 || m < n || ldfjac < n || tol < 0. || lwa < n * 5 +
|
|
||||||
m) {
|
|
||||||
/* goto L10; */
|
|
||||||
return info;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
/* call lmstr. */
|
ipvt.resize(n);
|
||||||
|
info = ei_lmstr<Functor,Scalar>(
|
||||||
maxfev = (n + 1) * 100;
|
x, fvec,
|
||||||
ftol = tol;
|
nfev, njev,
|
||||||
xtol = tol;
|
fjac, ipvt, diag,
|
||||||
gtol = 0.;
|
1,
|
||||||
mode = 1;
|
100.,
|
||||||
nprint = 0;
|
(n+1)*100,
|
||||||
info = lmstr(fcn, p, m, n, &x[1], &fvec[1], &fjac[fjac_offset], ldfjac,
|
tol, tol, Scalar(0.)
|
||||||
ftol, xtol, gtol, maxfev, &wa[1], mode, factor, nprint,
|
);
|
||||||
&nfev, &njev, &ipvt[1], &wa[n + 1], &wa[(n << 1) + 1], &
|
return (info==8)?4:info;
|
||||||
wa[n * 3 + 1], &wa[(n << 2) + 1], &wa[n * 5 + 1]);
|
}
|
||||||
if (info == 8) {
|
|
||||||
info = 4;
|
|
||||||
}
|
|
||||||
/* L10: */
|
|
||||||
return info;
|
|
||||||
|
|
||||||
/* last card of subroutine lmstr1. */
|
|
||||||
|
|
||||||
} /* lmstr1_ */
|
|
||||||
|
|
||||||
|
Loading…
x
Reference in New Issue
Block a user