26static int mp_fdjac2(mp_func funct,
27 int m,
int n,
int *ifree,
int npar,
double *x,
double *fvec,
28 double *fjac,
double epsfcn,
29 double *wa,
void *priv,
int *
nfev,
30 double *step,
double *dstep,
int *dside,
31 int *qulimited,
double *ulimit,
32 int *ddebug,
double *ddrtol,
double *ddatol);
33static void mp_qrfac(
int m,
int n,
double *a,
35 double *rdiag,
double *acnorm,
double *wa);
36static void mp_qrsolv(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
37 double *qtb,
double *x,
double *sdiag,
double *wa);
38static void mp_lmpar(
int n,
double *r,
int ldr,
int *ipvt,
int *ifree,
double *diag,
39 double *qtb,
double delta,
double *par,
double *x,
40 double *sdiag,
double *wa1,
double *wa2);
41static double mp_enorm(
int n,
double *x);
42static double mp_dmax1(
double a,
double b);
43static double mp_dmin1(
double a,
double b);
44static int mp_min0(
int a,
int b);
45static int mp_covar(
int n,
double *r,
int ldr,
int *ipvt,
double tol,
double *wa);
48#define mp_call(funct, m, n, x, fvec, dvec, priv) (*(funct))(m,n,x,fvec,dvec,priv)
51#define mp_malloc(dest,type,size) \
52 dest = (type *) malloc( sizeof(type)*size ); \
54 info = MP_ERR_MEMORY; \
58 for (_k=0; _k<(size); _k++) dest[_k] = 0; \
271int mpfit(mp_func funct,
int m,
int npar,
272 double *xall, mp_par *pars, mp_config *config,
void *private_data,
276 int i, j, info, iflag, nfree, npegged, iter;
280 double actred,delta,dirder,fnorm,fnorm1,gnorm, orignorm;
281 double par,pnorm,prered,ratio;
282 double sum,temp,temp1,temp2,temp3,xnorm, alpha;
283 static double one = 1.0;
284 static double p1 = 0.1;
285 static double p5 = 0.5;
286 static double p25 = 0.25;
287 static double p75 = 0.75;
288 static double p0001 = 1.0e-4;
289 static double zero = 0.0;
292 double *step = 0, *dstep = 0, *llim = 0, *ulim = 0;
293 int *pfixed = 0, *mpside = 0, *ifree = 0, *qllim = 0, *qulim = 0;
295 double *ddrtol = 0, *ddatol = 0;
297 double *fvec = 0, *qtf = 0;
298 double *x = 0, *xnew = 0, *fjac = 0, *diag = 0;
299 double *wa1 = 0, *wa2 = 0, *wa3 = 0, *wa4 = 0;
308 conf.stepfactor = 100.0;
310 conf.epsfcn = MP_MACHEP0;
312 conf.douserscale = 0;
315 conf.nofinitecheck = 0;
319 if (config->ftol > 0) conf.ftol = config->ftol;
320 if (config->xtol > 0) conf.xtol = config->xtol;
321 if (config->gtol > 0) conf.gtol = config->gtol;
322 if (config->stepfactor > 0) conf.stepfactor = config->stepfactor;
323 if (config->nprint >= 0) conf.nprint = config->nprint;
324 if (config->epsfcn > 0) conf.epsfcn = config->epsfcn;
325 if (config->maxiter > 0) conf.maxiter = config->maxiter;
326 if (config->douserscale != 0) conf.douserscale = config->douserscale;
327 if (config->covtol > 0) conf.covtol = config->covtol;
328 if (config->nofinitecheck > 0) conf.nofinitecheck = config->nofinitecheck;
329 conf.maxfev = config->maxfev;
341 if ((m <= 0) || (xall == 0)) {
342 return MP_ERR_NPOINTS;
355 mp_malloc(pfixed,
int, npar);
356 if (pars)
for (i=0; i<npar; i++) {
357 pfixed[i] = (pars[i].fixed)?1:0;
361 mp_malloc(step,
double, npar);
362 mp_malloc(dstep,
double, npar);
363 mp_malloc(mpside,
int, npar);
364 mp_malloc(ddebug,
int, npar);
365 mp_malloc(ddrtol,
double, npar);
366 mp_malloc(ddatol,
double, npar);
367 if (pars)
for (i=0; i<npar; i++) {
368 step[i] = pars[i].step;
369 dstep[i] = pars[i].relstep;
370 mpside[i] = pars[i].side;
371 ddebug[i] = pars[i].deriv_debug;
372 ddrtol[i] = pars[i].deriv_reltol;
373 ddatol[i] = pars[i].deriv_abstol;
378 mp_malloc(ifree,
int, npar);
379 for (i=0, j=0; i<npar; i++) {
380 if (pfixed[i] == 0) {
391 for (i=0; i<npar; i++) {
392 if ( (pars[i].limited[0] && (xall[i] < pars[i].limits[0])) ||
393 (pars[i].limited[1] && (xall[i] > pars[i].limits[1])) ) {
394 info = MP_ERR_INITBOUNDS;
397 if ( (pars[i].fixed == 0) && pars[i].limited[0] && pars[i].limited[1] &&
398 (pars[i].limits[0] >= pars[i].limits[1])) {
399 info = MP_ERR_BOUNDS;
404 mp_malloc(qulim,
int, nfree);
405 mp_malloc(qllim,
int, nfree);
406 mp_malloc(ulim,
double, nfree);
407 mp_malloc(llim,
double, nfree);
409 for (i=0; i<nfree; i++) {
410 qllim[i] = pars[ifree[i]].limited[0];
411 qulim[i] = pars[ifree[i]].limited[1];
412 llim[i] = pars[ifree[i]].limits[0];
413 ulim[i] = pars[ifree[i]].limits[1];
414 if (qllim[i] || qulim[i]) qanylim = 1;
419 if ((npar <= 0) || (conf.ftol <= 0) || (conf.xtol <= 0) ||
420 (conf.gtol <= 0) || (conf.maxiter < 0) ||
421 (conf.stepfactor <= 0)) {
433 mp_malloc(fvec,
double, m);
434 mp_malloc(qtf,
double, nfree);
435 mp_malloc(x,
double, nfree);
436 mp_malloc(xnew,
double, npar);
437 mp_malloc(fjac,
double, m*nfree);
439 mp_malloc(diag,
double, npar);
440 mp_malloc(wa1,
double, npar);
441 mp_malloc(wa2,
double, npar);
442 mp_malloc(wa3,
double, npar);
443 mp_malloc(wa4,
double, m);
444 mp_malloc(ipvt,
int, npar);
447 iflag = mp_call(funct, m, npar, xall, fvec, 0, private_data);
453 fnorm = mp_enorm(m, fvec);
454 orignorm = fnorm*fnorm;
457 for (i=0; i<npar; i++) {
462 for (i=0; i<nfree; i++) {
463 x[i] = xall[ifree[i]];
470 for (i=0; i<nfree; i++) {
476 for (i=0; i<nfree; i++) {
477 xnew[ifree[i]] = x[i];
483 iflag = mp_fdjac2(funct, m, nfree, ifree, npar, xnew, fvec, fjac,
484 conf.epsfcn, wa4, private_data, &
nfev,
485 step, dstep, mpside, qulim, ulim,
486 ddebug, ddrtol, ddatol);
494 for (j=0; j<nfree; j++) {
495 int lpegged = (qllim[j] && (x[j] == llim[j]));
496 int upegged = (qulim[j] && (x[j] == ulim[j]));
499 if (lpegged || upegged) {
502 for (i=0; i<m; i++, ij++) {
503 sum += fvec[i] * fjac[ij];
506 if (lpegged && (sum > 0)) {
508 for (i=0; i<m; i++, ij++) fjac[ij] = 0;
510 if (upegged && (sum < 0)) {
512 for (i=0; i<m; i++, ij++) fjac[ij] = 0;
518 mp_qrfac(m,nfree,fjac,1,ipvt,wa1,wa2,wa3);
525 if (conf.douserscale == 0) {
526 for (j=0; j<nfree; j++) {
527 diag[ifree[j]] = wa2[j];
528 if (wa2[j] == zero ) {
529 diag[ifree[j]] = one;
538 for (j=0; j<nfree; j++ ) {
539 wa3[j] = diag[ifree[j]] * x[j];
542 xnorm = mp_enorm(nfree, wa3);
543 delta = conf.stepfactor*xnorm;
544 if (delta == zero) delta = conf.stepfactor;
551 for (i=0; i<m; i++ ) {
556 for (j=0; j<nfree; j++ ) {
561 for (i=j; i<m; i++ ) {
562 sum += fjac[ij] * wa4[i];
567 for (i=j; i<m; i++ ) {
568 wa4[i] += fjac[ij] * temp;
581 if (conf.nofinitecheck) {
585 int off = 0, nonfinite = 0;
587 for (j=0; j<nfree; j++) {
588 for (i=0; i<nfree; i++) {
589 if (mpfinite(fjac[off+i]) == 0) nonfinite = 1;
607 for (j=0; j<nfree; j++ ) {
609 if (wa2[l] != zero) {
612 for (i=0; i<=j; i++ ) {
613 sum += fjac[ij]*(qtf[i]/fnorm);
616 gnorm = mp_dmax1(gnorm,fabs(sum/wa2[l]));
625 if (gnorm <= conf.gtol)
627 if (info != 0)
goto L300;
628 if (conf.maxiter == 0)
goto L300;
633 if (conf.douserscale == 0) {
634 for (j=0; j<nfree; j++ ) {
635 diag[ifree[j]] = mp_dmax1(diag[ifree[j]],wa2[j]);
646 mp_lmpar(nfree,fjac,ldfjac,ipvt,ifree,diag,qtf,delta,&par,wa1,wa2,wa3,wa4);
650 for (j=0; j<nfree; j++ ) {
657 for (j=0; j<nfree; j++ ) {
658 wa2[j] = x[j] + wa1[j];
666 for (j=0; j<nfree; j++) {
667 int lpegged = (qllim[j] && (x[j] <= llim[j]));
668 int upegged = (qulim[j] && (x[j] >= ulim[j]));
669 int dwa1 = fabs(wa1[j]) > MP_MACHEP0;
671 if (lpegged && (wa1[j] < 0)) wa1[j] = 0;
672 if (upegged && (wa1[j] > 0)) wa1[j] = 0;
674 if (dwa1 && qllim[j] && ((x[j] + wa1[j]) < llim[j])) {
675 alpha = mp_dmin1(alpha, (llim[j]-x[j])/wa1[j]);
677 if (dwa1 && qulim[j] && ((x[j] + wa1[j]) > ulim[j])) {
678 alpha = mp_dmin1(alpha, (ulim[j]-x[j])/wa1[j]);
683 for (j=0; j<nfree; j++) {
687 wa1[j] = wa1[j] * alpha;
688 wa2[j] = x[j] + wa1[j];
693 sgnu = (ulim[j] >= 0) ? (+1) : (-1);
694 sgnl = (llim[j] >= 0) ? (+1) : (-1);
695 ulim1 = ulim[j]*(1-sgnu*MP_MACHEP0) - ((ulim[j] == 0)?(MP_MACHEP0):0);
696 llim1 = llim[j]*(1+sgnl*MP_MACHEP0) + ((llim[j] == 0)?(MP_MACHEP0):0);
698 if (qulim[j] && (wa2[j] >= ulim1)) {
701 if (qllim[j] && (wa2[j] <= llim1)) {
708 for (j=0; j<nfree; j++ ) {
709 wa3[j] = diag[ifree[j]]*wa1[j];
712 pnorm = mp_enorm(nfree,wa3);
718 delta = mp_dmin1(delta,pnorm);
724 for (i=0; i<nfree; i++) {
725 xnew[ifree[i]] = wa2[i];
728 iflag = mp_call(funct, m, npar, xnew, wa4, 0, private_data);
730 if (iflag < 0)
goto L300;
732 fnorm1 = mp_enorm(m,wa4);
738 if ((p1*fnorm1) < fnorm) {
740 actred = one - temp * temp;
748 for (j=0; j<nfree; j++ ) {
753 for (i=0; i<=j; i++ ) {
754 wa3[i] += fjac[ij]*temp;
764 temp1 = mp_enorm(nfree,wa3)*alpha/fnorm;
765 temp2 = (sqrt(alpha*par)*pnorm)/fnorm;
766 prered = temp1*temp1 + (temp2*temp2)/p5;
767 dirder = -(temp1*temp1 + temp2*temp2);
774 if (prered != zero) {
775 ratio = actred/prered;
783 if (actred >= zero) {
786 temp = p5*dirder/(dirder + p5*actred);
788 if (((p1*fnorm1) >= fnorm)
792 delta = temp*mp_dmin1(delta,pnorm/p1);
795 if ((par == zero) || (ratio >= p75) ) {
804 if (ratio >= p0001) {
809 for (j=0; j<nfree; j++ ) {
811 wa2[j] = diag[ifree[j]]*x[j];
813 for (i=0; i<m; i++ ) {
816 xnorm = mp_enorm(nfree,wa2);
824 if ((fabs(actred) <= conf.ftol) && (prered <= conf.ftol) &&
825 (p5*ratio <= one) ) {
828 if (delta <= conf.xtol*xnorm) {
831 if ((fabs(actred) <= conf.ftol) && (prered <= conf.ftol) && (p5*ratio <= one)
842 if ((conf.maxfev > 0) && (
nfev >= conf.maxfev)) {
846 if (iter >= conf.maxiter) {
850 if ((fabs(actred) <= MP_MACHEP0) && (prered <= MP_MACHEP0) && (p5*ratio <= one) ) {
853 if (delta <= MP_MACHEP0*xnorm) {
856 if (gnorm <= MP_MACHEP0) {
866 if (ratio < p0001)
goto L200;
881 for (i=0; i<nfree; i++) {
882 xall[ifree[i]] = x[i];
885 if ((conf.nprint > 0) && (info > 0)) {
886 iflag = mp_call(funct, m, npar, xall, fvec, 0, private_data);
892 if (pars)
for (i=0; i<npar; i++) {
893 if ((pars[i].limited[0] && (pars[i].limits[0] == xall[i])) ||
894 (pars[i].limited[1] && (pars[i].limits[1] == xall[i]))) {
900 if (result && (result->covar || result->xerror)) {
901 mp_covar(nfree, fjac, ldfjac, ipvt, conf.covtol, wa2);
905 for (j=0; j<(npar*npar); j++) result->covar[j] = 0;
908 for (j=0; j<nfree; j++) {
909 for (i=0; i<nfree; i++) {
910 result->covar[ifree[j]*npar+ifree[i]] = fjac[j*ldfjac+i];
915 if (result->xerror) {
916 for (j=0; j<npar; j++) result->xerror[j] = 0;
918 for (j=0; j<nfree; j++) {
919 double cc = fjac[j*ldfjac+j];
920 if (cc > 0) result->xerror[ifree[j]] = sqrt(cc);
926 strcpy(result->version, MPFIT_VERSION);
927 result->bestnorm = mp_dmax1(fnorm,fnorm1);
928 result->bestnorm *= result->bestnorm;
929 result->orignorm = orignorm;
930 result->status = info;
931 result->niter = iter;
934 result->nfree = nfree;
935 result->npegged = npegged;
940 for (j=0; j<m; j++) result->resid[j] = fvec[j];
946 if (fvec) free(fvec);
949 if (xnew) free(xnew);
950 if (fjac) free(fjac);
951 if (diag) free(diag);
956 if (ipvt) free(ipvt);
957 if (pfixed) free(pfixed);
958 if (step) free(step);
959 if (dstep) free(dstep);
960 if (mpside) free(mpside);
961 if (ddebug) free(ddebug);
962 if (ddrtol) free(ddrtol);
963 if (ddatol) free(ddatol);
964 if (ifree) free(ifree);
965 if (qllim) free(qllim);
966 if (qulim) free(qulim);
967 if (llim) free(llim);
968 if (ulim) free(ulim);
978int mp_fdjac2(mp_func funct,
979 int m,
int n,
int *ifree,
int npar,
double *x,
double *fvec,
980 double *fjac,
double epsfcn,
981 double *wa,
void *priv,
int *
nfev,
982 double *step,
double *dstep,
int *dside,
983 int *qulimited,
double *ulimit,
984 int *ddebug,
double *ddrtol,
double *ddatol)
1066 static double zero = 0.0;
1068 int has_analytical_deriv = 0, has_numerical_deriv = 0;
1069 int has_debug_deriv = 0;
1071 temp = mp_dmax1(epsfcn,MP_MACHEP0);
1076 dvec = (
double **) malloc(
sizeof(
double **)*npar);
1077 if (dvec == 0)
return MP_ERR_MEMORY;
1078 for (j=0; j<npar; j++) dvec[j] = 0;
1081 for (j=0; j<(n*m); j++) fjac[j] = 0;
1085 for (j=0; j<n; j++) {
1086 if (dside && dside[ifree[j]] == 3 && ddebug[ifree[j]] == 0) {
1088 dvec[ifree[j]] = fjac + j*m;
1089 has_analytical_deriv = 1;
1090 }
else if (dside && ddebug[ifree[j]] == 1) {
1092 dvec[ifree[j]] = fjac + j*m;
1093 has_analytical_deriv = 1;
1094 has_numerical_deriv = 1;
1095 has_debug_deriv = 1;
1097 has_numerical_deriv = 1;
1103 if (has_analytical_deriv) {
1104 iflag = mp_call(funct, m, npar, x, wa, dvec, priv);
1106 if (iflag < 0 )
goto DONE;
1109 if (has_debug_deriv) {
1110 printf(
"FJAC DEBUG BEGIN\n");
1111 printf(
"# %10s %10s %10s %10s %10s %10s\n",
1112 "IPNT",
"FUNC",
"DERIV_U",
"DERIV_N",
"DIFF_ABS",
"DIFF_REL");
1116 if (has_numerical_deriv)
for (j=0; j<n; j++) {
1117 int dsidei = (dside)?(dside[ifree[j]]):(0);
1118 int debug = ddebug[ifree[j]];
1119 double dr = ddrtol[ifree[j]], da = ddatol[ifree[j]];
1123 printf(
"FJAC PARM %d\n", ifree[j]);
1127 if (dside && dsidei == 3)
continue;
1130 h = eps * fabs(temp);
1131 if (step && step[ifree[j]] > 0) h = step[ifree[j]];
1132 if (dstep && dstep[ifree[j]] > 0) h = fabs(dstep[ifree[j]]*temp);
1133 if (h == zero) h = eps;
1136 if ((dside && dsidei == -1) ||
1137 (dside && dsidei == 0 &&
1138 qulimited && ulimit && qulimited[j] &&
1139 (temp > (ulimit[j]-h)))) {
1143 x[ifree[j]] = temp + h;
1144 iflag = mp_call(funct, m, npar, x, wa, 0, priv);
1146 if (iflag < 0 )
goto DONE;
1153 for (i=0; i<m; i++, ij++) {
1154 fjac[ij] = (wa[i] - fvec[i])/h;
1158 for (i=0; i<m; i++, ij++) {
1159 double fjold = fjac[ij];
1160 fjac[ij] = (wa[i] - fvec[i])/h;
1161 if ((da == 0 && dr == 0 && (fjold != 0 || fjac[ij] != 0)) ||
1162 ((da != 0 || dr != 0) && (fabs(fjold-fjac[ij]) > da + fabs(fjold)*dr))) {
1163 printf(
" %10d %10.4g %10.4g %10.4g %10.4g %10.4g\n",
1164 i, fvec[i], fjold, fjac[ij], fjold-fjac[ij],
1165 (fjold == 0)?(0):((fjold-fjac[ij])/fjold));
1172 for (i=0; i<m; i++, ij++) {
1177 x[ifree[j]] = temp - h;
1178 iflag = mp_call(funct, m, npar, x, wa, 0, priv);
1180 if (iflag < 0 )
goto DONE;
1186 for (i=0; i<m; i++, ij++) {
1187 fjac[ij] = (fjac[ij] - wa[i])/(2*h);
1190 for (i=0; i<m; i++, ij++) {
1191 double fjold = fjac[ij];
1192 fjac[ij] = (fjac[ij] - wa[i])/(2*h);
1193 if ((da == 0 && dr == 0 && (fjold != 0 || fjac[ij] != 0)) ||
1194 ((da != 0 || dr != 0) && (fabs(fjold-fjac[ij]) > da + fabs(fjold)*dr))) {
1195 printf(
" %10d %10.4g %10.4g %10.4g %10.4g %10.4g\n",
1196 i, fvec[i], fjold, fjac[ij], fjold-fjac[ij],
1197 (fjold == 0)?(0):((fjold-fjac[ij])/fjold));
1205 if (has_debug_deriv) {
1206 printf(
"FJAC DEBUG END\n");
1210 if (dvec) free(dvec);
1211 if (iflag < 0)
return iflag;
1221void mp_qrfac(
int m,
int n,
double *a,
1222 int pivot,
int *ipvt,
1223 double *rdiag,
double *acnorm,
double *wa)
1302 int i,ij,jj,j,jp1,k,kmax,minmn;
1303 double ajnorm,sum,temp;
1304 static double zero = 0.0;
1305 static double one = 1.0;
1306 static double p05 = 0.05;
1315 for (j=0; j<n; j++) {
1316 acnorm[j] = mp_enorm(m,&a[ij]);
1317 rdiag[j] = acnorm[j];
1326 minmn = mp_min0(m,n);
1327 for (j=0; j<minmn; j++) {
1336 if (rdiag[k] > rdiag[kmax])
1352 rdiag[kmax] = rdiag[j];
1355 ipvt[j] = ipvt[kmax];
1364 ajnorm = mp_enorm(m-j,&a[jj]);
1383 for (k=jp1; k<n; k++)
1394 temp = sum/a[j+m*j];
1399 a[ij] -= temp*a[jj];
1403 if ((pivot != 0) && (rdiag[k] != zero))
1405 temp = a[j+m*k]/rdiag[k];
1406 temp = mp_dmax1( zero, one-temp*temp );
1407 rdiag[k] *= sqrt(temp);
1408 temp = rdiag[k]/wa[k];
1409 if ((p05*temp*temp) <= MP_MACHEP0)
1411 rdiag[k] = mp_enorm(m-j-1,&a[jp1+m*k]);
1429void mp_qrsolv(
int n,
double *r,
int ldr,
int *ipvt,
double *diag,
1430 double *qtb,
double *x,
double *sdiag,
double *wa)
1511 int i,ij,ik,kk,j,jp1,k,kp1,l,nsing;
1512 double cosx,cotan,qtbpj,sinx,sum,tanx,temp;
1513 static double zero = 0.0;
1514 static double p25 = 0.25;
1515 static double p5 = 0.5;
1522 for (j=0; j<n; j++) {
1539 for (j=0; j<n; j++) {
1545 if (diag[l] == zero)
1562 if (sdiag[k] == zero)
1565 if (fabs(r[kk]) < fabs(sdiag[k]))
1567 cotan = r[kk]/sdiag[k];
1568 sinx = p5/sqrt(p25+p25*cotan*cotan);
1573 tanx = sdiag[k]/r[kk];
1574 cosx = p5/sqrt(p25+p25*tanx*tanx);
1581 r[kk] = cosx*r[kk] + sinx*sdiag[k];
1582 temp = cosx*wa[k] + sinx*qtbpj;
1583 qtbpj = -sinx*wa[k] + cosx*qtbpj;
1592 for (i=kp1; i<n; i++)
1594 temp = cosx*r[ik] + sinx*sdiag[i];
1595 sdiag[i] = -sinx*r[ik] + cosx*sdiag[i];
1615 for (j=0; j<n; j++) {
1616 if ((sdiag[j] == zero) && (nsing == n))
1624 for (k=0; k<nsing; k++) {
1631 for (i=jp1; i<nsing; i++)
1637 wa[j] = (wa[j] - sum)/sdiag[j];
1643 for (j=0; j<n; j++) {
1655void mp_lmpar(
int n,
double *r,
int ldr,
int *ipvt,
int *ifree,
double *diag,
1656 double *qtb,
double delta,
double *par,
double *x,
1657 double *sdiag,
double *wa1,
double *wa2)
1754 int i,iter,ij,jj,j,jm1,jp1,k,l,nsing;
1755 double dxnorm,fp,gnorm,parc,parl,paru;
1757 static double zero = 0.0;
1759 static double p1 = 0.1;
1760 static double p001 = 0.001;
1768 for (j=0; j<n; j++) {
1770 if ((r[jj] == zero) && (nsing == n))
1778 for (k=0; k<nsing; k++)
1781 wa1[j] = wa1[j]/r[j+ldr*j];
1787 for (i=0; i<=jm1; i++)
1789 wa1[i] -= r[ij]*temp;
1796 for (j=0; j<n; j++) {
1807 wa2[j] = diag[ifree[j]]*x[j];
1808 dxnorm = mp_enorm(n,wa2);
1809 fp = dxnorm - delta;
1810 if (fp <= p1*delta) {
1823 wa1[j] = diag[ifree[l]]*(wa2[l]/dxnorm);
1833 for (i=0; i<=jm1; i++)
1835 sum += r[ij]*wa1[i];
1839 wa1[j] = (wa1[j] - sum)/r[j+ldr*j];
1842 temp = mp_enorm(n,wa1);
1843 parl = ((fp/delta)/temp)/temp;
1849 for (j=0; j<n; j++) {
1852 for (i=0; i<=j; i++)
1854 sum += r[ij]*qtb[i];
1858 wa1[j] = sum/diag[ifree[l]];
1861 gnorm = mp_enorm(n,wa1);
1864 paru = MP_DWARF/mp_dmin1(delta,p1);
1869 *par = mp_dmax1( *par,parl);
1870 *par = mp_dmin1( *par,paru);
1872 *par = gnorm/dxnorm;
1883 *par = mp_dmax1(MP_DWARF,p001*paru);
1884 temp = sqrt( *par );
1886 wa1[j] = temp*diag[ifree[j]];
1887 mp_qrsolv(n,r,ldr,ipvt,wa1,qtb,x,sdiag,wa2);
1889 wa2[j] = diag[ifree[j]]*x[j];
1890 dxnorm = mp_enorm(n,wa2);
1892 fp = dxnorm - delta;
1898 if ((fabs(fp) <= p1*delta)
1899 || ((parl == zero) && (fp <= temp) && (temp < zero))
1905 for (j=0; j<n; j++) {
1907 wa1[j] = diag[ifree[l]]*(wa2[l]/dxnorm);
1910 for (j=0; j<n; j++) {
1911 wa1[j] = wa1[j]/sdiag[j];
1917 for (i=jp1; i<n; i++)
1919 wa1[i] -= r[ij]*temp;
1925 temp = mp_enorm(n,wa1);
1926 parc = ((fp/delta)/temp)/temp;
1931 parl = mp_dmax1(parl, *par);
1933 paru = mp_dmin1(paru, *par);
1937 *par = mp_dmax1(parl, *par + parc);
1958double mp_enorm(
int n,
double *x)
2000 double agiant,floatn,s1,s2,s3,xabs,x1max,x3max;
2002 double rdwarf = MP_RDWARF;
2003 double rgiant = MP_RGIANT;
2004 static double zero = 0.0;
2005 static double one = 1.0;
2013 agiant = rgiant/floatn;
2015 for (i=0; i<n; i++) {
2017 if ((xabs > rdwarf) && (xabs < agiant))
2034 s1 = one + s1*temp*temp;
2050 s3 = one + s3*temp*temp;
2066 temp = s1 + (s2/x1max)/x1max;
2067 ans = x1max*sqrt(temp);
2072 temp = s2*(one+(x3max/s2)*(x3max*s3));
2074 temp = x3max*((s2/x3max)+(x3max*s3));
2079 ans = x3max*sqrt(s3);
2090double mp_dmax1(
double a,
double b)
2099double mp_dmin1(
double a,
double b)
2108int mp_min0(
int a,
int b)
2186int mp_covar(
int n,
double *r,
int ldr,
int *ipvt,
double tol,
double *wa)
2188 int i, ii, j, jj, k, l;
2189 int kk, kj, ji, j0, k0, jj0;
2191 double one = 1.0, temp, tolr, zero = 0.0;
2198 for (j=0; j<n; j++) {
2199 for (i=0; i<n; i++) {
2200 printf(
"%f ", r[j*ldr+i]);
2206 tolr = tol*fabs(r[0]);
2208 for (k=0; k<n; k++) {
2210 if (fabs(r[kk]) <= tolr)
break;
2213 for (j=0; j<k; j++) {
2215 temp = r[kk] * r[kj];
2218 k0 = k*ldr; j0 = j*ldr;
2219 for (i=0; i<=j; i++) {
2220 r[k0+i] += (-temp*r[j0+i]);
2232 for (k=0; k <= l; k++) {
2235 for (j=0; j<k; j++) {
2239 for (i=0; i<=j; i++) {
2240 r[j0+i] += temp*r[k0+i];
2245 for (i=0; i<=k; i++) {
2255 for (j=0; j<n; j++) {
2260 for (i=0; i<=j; i++) {
2263 if (sing) r[ji] = zero;
2265 if (ii > jj) r[jj0+ii] = r[ji];
2266 if (ii < jj) r[ii*ldr+jj] = r[ji];
2274 for (j=0; j<n; j++) {
2276 for (i=0; i<j; i++) {
2277 r[j0+i] = r[i*ldr+j];
2283 for (j=0; j<n; j++) {
2284 for (i=0; i<n; i++) {
2285 printf(
"%f ", r[j*ldr+i]);