43 #include <visp/vpLevenbergMarquartd.h>
44 #include <visp/vpMath.h>
51 #define SIGN(x) ((x) < 0 ? -1 : 1)
52 #define SWAP(a,b,c) {(c) = (a); (a) = (b); (b) = (c);}
53 #define MIJ(m,i,j,s) ((m) + ((long) (i) * (long) (s)) + (long) (j))
84 double enorm (
const double *x,
int n)
86 const double rdwarf = 3.834e-20;
87 const double rgiant = 1.304e19;
90 double agiant, floatn;
91 double norm_eucl = 0.0;
92 double s1 = 0.0, s2 = 0.0, s3 = 0.0;
94 double x1max = 0.0 , x3max = 0.0;
97 agiant = rgiant / floatn;
99 for (i = 0; i < n; i++)
101 xabs = std::fabs(x[i]);
102 if ((xabs > rdwarf) && (xabs < agiant))
110 else if (xabs <= rdwarf)
118 if (xabs > std::numeric_limits<double>::epsilon())
119 s3 += (xabs / x3max) * (xabs / x3max);
123 s3 = 1.0 + s3 * (x3max / xabs) * (x3max / xabs);
135 s1 += (xabs / x1max) * (xabs / x1max);
139 s1 = 1.0 + s1 * (x1max / xabs) * (x1max / xabs);
149 if (std::fabs(s1) <= std::numeric_limits<double>::epsilon())
152 if (std::fabs(s2) <= std::numeric_limits<double>::epsilon())
153 norm_eucl = x3max * sqrt(s3);
154 else if (s2 >= x3max)
155 norm_eucl = sqrt (s2 * (1.0 + ( x3max / s2) * (x3max * s3)));
157 norm_eucl = sqrt (x3max * ((s2 / x3max) + (x3max * s3)));
160 norm_eucl = x1max * sqrt (s1 + (s2 / x1max) / x1max);
236 int lmpar(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
237 double *delta,
double *par,
double *x,
double *sdiag,
double *wa1,
double *wa2)
239 const double tol1 = 0.1, tol001 = 0.001;
245 double dxnorm, fp, gnorm, parc;
248 double dwarf = DBL_MIN;
257 for (
int i = 0; i < n; i++)
260 double *pt = MIJ(r, i, i, ldr);
262 if (std::fabs(*pt) <= std::numeric_limits<double>::epsilon() && nsing == n)
271 for (
int k = 0; k < nsing; k++)
273 int i = nsing - 1 - k;
274 wa1[i] /= *MIJ(r, i, i, ldr);
280 for (
unsigned int j = 0; j <= (
unsigned int)jm1; j++)
281 wa1[j] -= *MIJ(r, i, j, ldr) * temp;
286 for (
int j = 0; j < n; j++)
299 for (
int i = 0; i < n; i++)
300 wa2[i] = diag[i] * x[i];
302 dxnorm = enorm(wa2, n);
304 fp = dxnorm - *delta;
306 if (fp > tol1 * (*delta))
317 for (
int i = 0; i < n; i++)
320 wa1[i] = diag[l] * (wa2[l] / dxnorm);
323 for (
int i = 0; i < n; i++)
331 for (
unsigned int j = 0; j <= (
unsigned int)im1; j++)
332 sum += (*MIJ(r, i, j, ldr) * wa1[j]);
334 wa1[i] = (wa1[i] - sum) / *MIJ(r, i, i, ldr);
337 temp = enorm(wa1, n);
338 parl = ((fp / *delta) / temp) / temp;
345 for (
int i = 0; i < n; i++)
349 for (
int j = 0; j <= i; j++)
350 sum += *MIJ(r, i, j, ldr) * qtb[j];
353 wa1[i] = sum / diag[l];
356 gnorm = enorm(wa1, n);
357 paru = gnorm / *delta;
360 if (std::fabs(paru) <= std::numeric_limits<double>::epsilon())
372 if (std::fabs(*par) <= std::numeric_limits<double>::epsilon() )
373 *par = gnorm / dxnorm;
387 if (std::fabs(*par) <= std::numeric_limits<double>::epsilon() )
392 for (
int i = 0; i < n; i++)
393 wa1[i] = temp * diag[i];
395 qrsolv(n, r, ldr, ipvt, wa1, qtb, x, sdiag, wa2);
397 for (
int i = 0; i < n; i++)
398 wa2[i] = diag[i] * x[i];
400 dxnorm = enorm(wa2, n);
402 fp = dxnorm - *delta;
412 if ((std::fabs(fp) <= tol1 * (*delta)) || ((std::fabs(parl) <= std::numeric_limits<double>::epsilon()) && (fp <= temp)
413 && (temp < 0.0)) || (iter == 10))
428 for (
int i = 0; i < n; i++)
431 wa1[i] = diag[l] * (wa2[l] / dxnorm);
434 for (
unsigned int i = 0; i < (
unsigned int)n; i++)
436 wa1[i] = wa1[i] / sdiag[i];
439 if ( (
unsigned int) n >= jp1)
441 for (
unsigned int j = jp1; j < (
unsigned int)n; j++)
442 wa1[j] -= (*MIJ(r, i, j, ldr) * temp);
446 temp = enorm(wa1, n);
447 parc = ((fp / *delta) / temp) / temp;
489 double pythag (
double a,
double b)
491 double pyth, p, r, s, t, u;
496 if (std::fabs(p) <= std::numeric_limits<double>::epsilon() )
502 r = (std::min(std::fabs(a), std::fabs(b)) / p) * (std::min(std::fabs(a), std::fabs(b)) / p);
506 while (std::fabs(t - 4.0) < std::fabs(
vpMath::maximum(t,4.0)) * std::numeric_limits<double>::epsilon() )
511 r *= (s / u) * (s / u);
566 int qrfac(
int m,
int n,
double *a,
int lda,
int *pivot,
int *ipvt,
567 int ,
double *rdiag,
double *acnorm,
double *wa)
569 const double tolerance = 0.05;
571 int i, j, ip1, k, kmax, minmn;
572 double ajnorm, epsmch;
573 double sum, temp, tmp;
578 epsmch = DBL_EPSILON;
584 for (i = 0; i < m; i++)
586 acnorm[i] = enorm(MIJ(a, i, 0, lda), n);
587 rdiag[i] = acnorm[i];
597 for (i = 0; i < minmn; i++)
606 for (k = i; k < m; k++)
608 if (rdiag[k] > rdiag[kmax])
614 for (j = 0; j < n; j++)
615 SWAP(*MIJ(a, i, j, lda),
616 *MIJ(a, kmax, j, lda), tmp);
618 rdiag[kmax] = rdiag[i];
621 SWAP( ipvt[i], ipvt[kmax], k);
629 ajnorm = enorm(MIJ(a, i, i, lda), n - i);
632 if (std::fabs(ajnorm) > std::numeric_limits<double>::epsilon() )
634 if (*MIJ(a, i, i, lda) < 0.0)
637 for (j = i; j < n; j++)
638 *MIJ(a, i, j, lda) /= ajnorm;
639 *MIJ(a, i, i, lda) += 1.0;
649 for (k = ip1; k < m; k++)
652 for (j = i; j < n; j++)
653 sum += *MIJ(a, i, j, lda) * *MIJ(a, k, j, lda);
655 temp = sum / *MIJ(a, i, i, lda);
657 for (j = i; j < n; j++)
658 *MIJ(a, k, j, lda) -= temp * *MIJ(a, i, j, lda);
661 if (pivot && (std::fabs(rdiag[k]) > std::numeric_limits<double>::epsilon()) )
663 temp = *MIJ (a, k, i, lda) / rdiag[k];
666 if (tolerance * (rdiag[k] / wa[k]) * (rdiag[k] / wa[k]) <= epsmch)
668 rdiag[k] = enorm(MIJ(a, k, ip1, lda), (n -1 - (
int) i));
742 int qrsolv (
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
double *qtb,
743 double *x,
double *sdiag,
double *wa)
745 int i, j, jp1, k, kp1, l;
747 double cosi, cotg, qtbpj, sinu, sum, tg, temp;
754 for (i = 0; i < n; i++)
756 for (j = i; j < n; j++)
757 *MIJ(r, i, j, ldr) = *MIJ(r, j, i, ldr);
759 x[i] = *MIJ(r, i, i, ldr);
768 for (i = 0; i < n; i++)
778 if (std::fabs(diag[l]) > std::numeric_limits<double>::epsilon())
780 for (k = i; k < n; k++)
794 for (k = i; k < n; k++)
803 if (std::fabs(sdiag[k]) > std::numeric_limits<double>::epsilon())
805 if (std::fabs(*MIJ(r, k, k, ldr)) >= std::fabs(sdiag[k]))
807 tg = sdiag[k] / *MIJ(r, k, k, ldr);
808 cosi = 0.5 / sqrt(0.25 + 0.25 * (tg * tg));
813 cotg = *MIJ(r, k, k, ldr) / sdiag[k];
814 sinu = 0.5 / sqrt(0.25 + 0.25 * (cotg * cotg));
823 *MIJ(r, k, k, ldr) = cosi * *MIJ(r, k, k, ldr) + sinu * sdiag[k];
824 temp = cosi * wa[k] + sinu * qtbpj;
825 qtbpj = -sinu * wa[k] + cosi * qtbpj;
837 for (j = kp1; j < n; j++)
839 temp = cosi * *MIJ(r, k, j, ldr) +
841 sdiag[j] = - sinu * *MIJ(r, k, j, ldr) +
843 *MIJ(r, k, j, ldr) = temp;
854 sdiag[i] = *MIJ(r, i, i, ldr);
855 *MIJ(r, i, i, ldr) = x[i];
864 for (i = 0; i < n; i++)
867 if ( (std::fabs(sdiag[i]) <= std::numeric_limits<double>::epsilon()) && nsing == n)
876 for (k = 0; k < nsing; k++)
884 for (j = jp1; j < nsing; j++)
885 sum += *MIJ(r, i, j, ldr) * wa[j];
887 wa[i] = (wa[i] - sum) / sdiag[i];
894 for (j = 0; j < n; j++)
1014 int lmder (
void (*ptr_fcn)(
int m,
int n,
double *xc,
double *fvecc,
1015 double *jac,
int ldfjac,
int iflag),
int m,
int n,
double *x,
1016 double *fvec,
double *fjac,
int ldfjac,
double ftol,
double xtol,
1017 double gtol,
unsigned int maxfev,
double *diag,
int mode,
1018 const double factor,
int nprint,
int *info,
unsigned int *nfev,
1019 int *njev,
int *ipvt,
double *qtf,
double *wa1,
double *wa2,
1020 double *wa3,
double *wa4)
1022 const double tol1 = 0.1, tol5 = 0.5, tol25 = 0.25, tol75 = 0.75, tol0001 = 0.0001;
1027 double actred, delta, dirder, epsmch, fnorm, fnorm1, gnorm;
1028 double ratio = DBL_EPSILON;
1029 double par, pnorm, prered;
1030 double sum, temp, temp1, temp2, xnorm = 0.0;
1033 epsmch = DBL_EPSILON;
1058 if ((n <= 0) || (m < n) || (ldfjac < m) || (ftol < 0.0) || (xtol < 0.0)
1059 || (gtol < 0.0) || (maxfev == 0) || (factor <= 0.0))
1070 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1077 for (j = 0; j < n; j++)
1090 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1103 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1119 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1124 fnorm = enorm(fvec, m);
1137 while (count < (
int)maxfev)
1146 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1161 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1172 if ((iter-1) % nprint == 0)
1173 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1186 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1195 qrfac(n, m, fjac, ldfjac, &oncol, ipvt, n, wa1, wa2, wa3);
1206 for (j = 0; j < n; j++)
1210 if (std::fabs(wa2[j]) <= std::numeric_limits<double>::epsilon())
1221 for (j = 0; j < n; j++)
1222 wa3[j] = diag[j] * x[j];
1224 xnorm = enorm (wa3, n);
1225 delta = factor * xnorm;
1228 if (std::fabs(delta) <= std::numeric_limits<double>::epsilon())
1236 for (i = 0; i < m; i++)
1239 for (i = 0; i < n; i++)
1241 double *pt = MIJ(fjac, i, i, ldfjac);
1243 if (std::fabs(*pt) > std::numeric_limits<double>::epsilon() )
1247 for (j = i; j < m; j++)
1248 sum += *MIJ(fjac, i, j, ldfjac) * wa4[j];
1250 temp = - sum / *MIJ(fjac, i, i, ldfjac);
1252 for (j = i; j < m; j++)
1253 wa4[j] += *MIJ(fjac, i, j, ldfjac) * temp;
1257 *MIJ(fjac, i, i, ldfjac) = wa1[i];
1268 if (std::fabs(fnorm) > std::numeric_limits<double>::epsilon())
1270 for (i = 0; i < n; i++)
1274 if (std::fabs(wa2[l]) > std::numeric_limits<double>::epsilon())
1277 for (j = 0; j <= i; j++)
1278 sum += *MIJ(fjac, i, j, ldfjac) * (qtf[j] / fnorm);
1303 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1314 for (j = 0; j < n; j++)
1322 while (ratio < tol0001)
1328 lmpar(n, fjac, ldfjac, ipvt, diag, qtf, &delta, &par, wa1,
1335 for (j = 0; j < n; j++)
1338 wa2[j] = x[j] + wa1[j];
1339 wa3[j] = diag[j] * wa1[j];
1342 pnorm = enorm(wa3, n);
1357 (*ptr_fcn)(m, n, wa2, wa4, fjac, ldfjac, iflag);
1372 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1377 fnorm1 = enorm(wa4, m);
1385 if ((tol1 * fnorm1) < fnorm)
1386 actred = 1.0 - ((fnorm1 / fnorm) * (fnorm1 / fnorm));
1393 for (i = 0; i < n; i++)
1398 for (j = 0; j <= i; j++)
1399 wa3[j] += *MIJ(fjac, i, j, ldfjac) * temp;
1402 temp1 = enorm(wa3, n) / fnorm;
1403 temp2 = (sqrt(par) * pnorm) / fnorm;
1404 prered = (temp1 * temp1) + (temp2 * temp2) / tol5;
1405 dirder = - ((temp1 * temp1) + (temp2 * temp2));
1414 if (std::fabs(prered) > std::numeric_limits<double>::epsilon())
1415 ratio = actred / prered;
1424 if ((std::fabs(par) <= std::numeric_limits<double>::epsilon()) || (ratio <= tol75))
1426 delta = pnorm / tol5;
1436 temp = tol5 * dirder / (dirder + tol5 * actred);
1438 if ((tol1 * fnorm1 >= fnorm) || (temp < tol1))
1448 if (ratio >= tol0001)
1455 for (j = 0; j < n; j++)
1458 wa2[j] = diag[j] * x[j];
1461 for (i = 0; i < m; i++)
1464 xnorm = enorm(wa2, n);
1473 if ((std::fabs(actred) <= ftol) && (prered <= ftol) && (tol5 * ratio <= 1.0))
1476 if (delta <= xtol * xnorm)
1479 if ((std::fabs(actred) <= ftol) && (prered <= ftol) && (tol5 * ratio <= 1.0)
1494 (*ptr_fcn)(m,n,x,fvec,fjac,ldfjac, iflag);
1503 if (*nfev >= maxfev)
1506 if ((std::fabs(actred) <= epsmch) && (prered <= epsmch)
1507 && (tol5 * ratio <= 1.0))
1510 if (delta <= epsmch * xnorm)
1513 if (gnorm <= epsmch)
1527 (*ptr_fcn)(m, n, x, fvec, fjac, ldfjac, iflag);
1618 int lmder1 (
void (*ptr_fcn)(
int m,
int n,
double *xc,
double *fvecc,
1619 double *jac,
int ldfjac,
int iflag),
1620 int m,
int n,
double *x,
double *fvec,
double *fjac,
1621 int ldfjac,
double tol,
int *info,
int *ipvt,
int lwa,
double *wa)
1623 const double factor = 100.0;
1624 unsigned int maxfev, nfev;
1625 int mode, njev, nprint;
1626 double ftol, gtol, xtol;
1633 if ( (m < n) || (ldfjac < m) || (tol < 0.0) ||
1634 (lwa < (5 * n + m)) )
1636 printf(
"%d %d %d %d \n", (m < n), (ldfjac < m), (tol < 0.0) , (lwa < (5 * n + m))) ;
1642 maxfev = (
unsigned int)(100 * (n + 1));
1649 lmder (ptr_fcn , m, n, x, fvec, fjac, ldfjac, ftol, xtol, gtol, maxfev, wa,
1650 mode, factor, nprint, info, &nfev, &njev, ipvt, &wa[n], &wa[2 * n],
1651 &wa[3 * n], &wa[4 * n], &wa[5 * n]);
static Type maximum(const Type &a, const Type &b)
static Type minimum(const Type &a, const Type &b)