/* * MINPACK-1 Least Squares Fitting Library * * Original public domain version by B. Garbow, K. Hillstrom, J. More' * (Argonne National Laboratory, MINPACK project, March 1980) * See the file DISCLAIMER for copyright information. * * Tranlation to C Language by S. Moshier (moshier.net) * * Enhancements and packaging by C. Markwardt * (comparable to IDL fitting routine MPFIT * see http://cow.physics.wisc.edu/~craigm/idl/idl.html) */ /* Main mpfit library routines (double precision) $Id: mpfit.c,v 1.17 2009/02/18 23:08:49 craigm Exp $ */ #include #include #include #include #include "mpfit.h" //#include "idl_export.h" /* Forward declarations of functions in this module */ int mp_fdjac2(mp_func funct, int m, int n, int *ifree, int npar, double *x, double *fvec, double *fjac, int ldfjac, double epsfcn, double *wa, void *priv, int *nfev, double *step, double *dstep, int *dside, int *qulimited, double *ulimit, int *ddebug, double *ddrtol, double *ddatol); void mp_qrfac(int m, int n, double *a, int lda, int pivot, int *ipvt, int lipvt, double *rdiag, double *acnorm, double *wa); void mp_qrsolv(int n, double *r, int ldr, int *ipvt, double *diag, double *qtb, double *x, double *sdiag, double *wa); void mp_lmpar(int n, double *r, int ldr, int *ipvt, int *ifree, double *diag, double *qtb, double delta, double *par, double *x, double *sdiag, double *wa1, double *wa2); double mp_enorm(int n, double *x); double mp_dmax1(double a, double b); double mp_dmin1(double a, double b); int mp_min0(int a, int b); int mp_covar(int n, double *r, int ldr, int *ipvt, double tol, double *wa); /* Macro to call user function */ #define mp_call(funct, m, n, x, fvec, dvec, priv) (*(funct))(m,n,x,fvec,dvec,priv) /* Macro to safely allocate memory */ #define mp_malloc(dest,type,size) \ dest = (type *) malloc( sizeof(type)*size ); \ if (dest == 0) { \ info = MP_ERR_MEMORY; \ goto CLEANUP; \ } else { \ int _k; \ for (_k=0; _k<(size); _k++) dest[_k] = 0; \ } /* * ********** * * subroutine mpfit * * the purpose of mpfit is to minimize the sum of the squares of * m nonlinear functions in n variables by a modification of * the levenberg-marquardt algorithm. the user must provide a * subroutine which calculates the functions. the jacobian is * then calculated by a finite-difference approximation. * * mp_funct funct - function to be minimized * int m - number of data points * int npar - number of fit parameters * double *xall - array of n initial parameter values * upon return, contains adjusted parameter values * mp_par *pars - array of npar structures specifying constraints; * or 0 (null pointer) for unconstrained fitting * [ see README and mpfit.h for definition & use of mp_par] * mp_config *config - pointer to structure which specifies the * configuration of mpfit(); or 0 (null pointer) * if the default configuration is to be used. * See README and mpfit.h for definition and use * of config. * void *private - any private user data which is to be passed directly * to funct without modification by mpfit(). * mp_result *result - pointer to structure, which upon return, contains * the results of the fit. The user should zero this * structure. If any of the array values are to be * returned, the user should allocate storage for them * and assign the corresponding pointer in *result. * Upon return, *result will be updated, and * any of the non-null arrays will be filled. * * * FORTRAN DOCUMENTATION BELOW * * * the subroutine statement is * * subroutine lmdif(fcn,m,n,x,fvec,ftol,xtol,gtol,maxfev,epsfcn, * diag,mode,factor,nprint,info,nfev,fjac, * ldfjac,ipvt,qtf,wa1,wa2,wa3,wa4) * * where * * fcn is the name of the user-supplied subroutine which * calculates the functions. fcn must be declared * in an external statement in the user calling * program, and should be written as follows. * * subroutine fcn(m,n,x,fvec,iflag) * integer m,n,iflag * double precision x(n),fvec(m) * ---------- * calculate the functions at x and * return this vector in fvec. * ---------- * return * end * * the value of iflag should not be changed by fcn unless * the user wants to terminate execution of lmdif. * in this case set iflag to a negative integer. * * m is a positive integer input variable set to the number * of functions. * * n is a positive integer input variable set to the number * of variables. n must not exceed m. * * x is an array of length n. on input x must contain * an initial estimate of the solution vector. on output x * contains the final estimate of the solution vector. * * fvec is an output array of length m which contains * the functions evaluated at the output x. * * ftol is a nonnegative input variable. termination * occurs when both the actual and predicted relative * reductions in the sum of squares are at most ftol. * therefore, ftol measures the relative error desired * in the sum of squares. * * xtol is a nonnegative input variable. termination * occurs when the relative error between two consecutive * iterates is at most xtol. therefore, xtol measures the * relative error desired in the approximate solution. * * gtol is a nonnegative input variable. termination * occurs when the cosine of the angle between fvec and * any column of the jacobian is at most gtol in absolute * value. therefore, gtol measures the orthogonality * desired between the function vector and the columns * of the jacobian. * * maxfev is a positive integer input variable. termination * occurs when the number of calls to fcn is at least * maxfev by the end of an iteration. * * epsfcn is an input variable used in determining a suitable * step length for the forward-difference approximation. this * approximation assumes that the relative errors in the * functions are of the order of epsfcn. if epsfcn is less * than the machine precision, it is assumed that the relative * errors in the functions are of the order of the machine * precision. * * diag is an array of length n. if mode = 1 (see * below), diag is internally set. if mode = 2, diag * must contain positive entries that serve as * multiplicative scale factors for the variables. * * mode is an integer input variable. if mode = 1, the * variables will be scaled internally. if mode = 2, * the scaling is specified by the input diag. other * values of mode are equivalent to mode = 1. * * factor is a positive input variable used in determining the * initial step bound. this bound is set to the product of * factor and the euclidean norm of diag*x if nonzero, or else * to factor itself. in most cases factor should lie in the * interval (.1,100.). 100. is a generally recommended value. * * nprint is an integer input variable that enables controlled * printing of iterates if it is positive. in this case, * fcn is called with iflag = 0 at the beginning of the first * iteration and every nprint iterations thereafter and * immediately prior to return, with x and fvec available * for printing. if nprint is not positive, no special calls * of fcn with iflag = 0 are made. * * info is an integer output variable. if the user has * terminated execution, info is set to the (negative) * value of iflag. see description of fcn. otherwise, * info is set as follows. * * info = 0 improper input parameters. * * info = 1 both actual and predicted relative reductions * in the sum of squares are at most ftol. * * info = 2 relative error between two consecutive iterates * is at most xtol. * * info = 3 conditions for info = 1 and info = 2 both hold. * * info = 4 the cosine of the angle between fvec and any * column of the jacobian is at most gtol in * absolute value. * * info = 5 number of calls to fcn has reached or * exceeded maxfev. * * info = 6 ftol is too small. no further reduction in * the sum of squares is possible. * * info = 7 xtol is too small. no further improvement in * the approximate solution x is possible. * * info = 8 gtol is too small. fvec is orthogonal to the * columns of the jacobian to machine precision. * * nfev is an integer output variable set to the number of * calls to fcn. * * fjac is an output m by n array. the upper n by n submatrix * of fjac contains an upper triangular matrix r with * diagonal elements of nonincreasing magnitude such that * * t t t * p *(jac *jac)*p = r *r, * * where p is a permutation matrix and jac is the final * calculated jacobian. column j of p is column ipvt(j) * (see below) of the identity matrix. the lower trapezoidal * part of fjac contains information generated during * the computation of r. * * ldfjac is a positive integer input variable not less than m * which specifies the leading dimension of the array fjac. * * ipvt is an integer output array of length n. ipvt * defines a permutation matrix p such that jac*p = q*r, * where jac is the final calculated jacobian, q is * orthogonal (not stored), and r is upper triangular * with diagonal elements of nonincreasing magnitude. * column j of p is column ipvt(j) of the identity matrix. * * qtf is an output array of length n which contains * the first n elements of the vector (q transpose)*fvec. * * wa1, wa2, and wa3 are work arrays of length n. * * wa4 is a work array of length m. * * subprograms called * * user-supplied ...... fcn * * minpack-supplied ... dpmpar,enorm,fdjac2,lmpar,qrfac * * fortran-supplied ... dabs,dmax1,dmin1,dsqrt,mod * * argonne national laboratory. minpack project. march 1980. * burton s. garbow, kenneth e. hillstrom, jorge j. more * * ********** */ int mpfit(mp_func funct, int m, int npar, double *xall, mp_par *pars, mp_config *config, void *private, mp_result *result) { mp_config conf; int i, j, info, iflag, nfree, npegged, iter; int qanylim = 0, qanypegged = 0; int ij,jj,l; double actred,delta,dirder,fnorm,fnorm1,gnorm, orignorm; double par,pnorm,prered,ratio; double sum,temp,temp1,temp2,temp3,xnorm, alpha; static double one = 1.0; static double p1 = 0.1; static double p5 = 0.5; static double p25 = 0.25; static double p75 = 0.75; static double p0001 = 1.0e-4; static double zero = 0.0; int nfev = 0; double *step = 0, *dstep = 0, *llim = 0, *ulim = 0; int *pfixed = 0, *mpside = 0, *ifree = 0, *qllim = 0, *qulim = 0; int *ddebug = 0; double *ddrtol = 0, *ddatol = 0; double *fvec = 0, *qtf = 0; double *x = 0, *xnew = 0, *fjac = 0, *diag = 0; double *wa1 = 0, *wa2 = 0, *wa3 = 0, *wa4 = 0; int *ipvt = 0; int ldfjac; /* Default configuration */ conf.ftol = 1e-10; conf.xtol = 1e-10; conf.gtol = 1e-10; conf.stepfactor = 100.0; conf.nprint = 1; conf.epsfcn = MP_MACHEP0; conf.maxiter = 200; conf.douserscale = 0; conf.maxfev = 0; conf.covtol = 1e-14; conf.nofinitecheck = 0; if (config) { /* Transfer any user-specified configurations */ if (config->ftol > 0) conf.ftol = config->ftol; if (config->xtol > 0) conf.xtol = config->xtol; if (config->gtol > 0) conf.gtol = config->gtol; if (config->stepfactor > 0) conf.stepfactor = config->stepfactor; if (config->nprint >= 0) conf.nprint = config->nprint; if (config->epsfcn > 0) conf.epsfcn = config->epsfcn; if (config->maxiter > 0) conf.maxiter = config->maxiter; if (config->douserscale != 0) conf.douserscale = config->douserscale; if (config->covtol > 0) conf.covtol = config->covtol; if (config->nofinitecheck > 0) conf.nofinitecheck = config->nofinitecheck; conf.maxfev = config->maxfev; } info = 0; iflag = 0; nfree = 0; npegged = 0; if (funct == 0) { return MP_ERR_FUNC; } if ((m <= 0) || (xall == 0)) { return MP_ERR_NPOINTS; } if (npar <= 0) { return MP_ERR_NFREE; } fnorm = -1.0; fnorm1 = -1.0; xnorm = -1.0; delta = 0.0; /* FIXED parameters? */ mp_malloc(pfixed, int, npar); if (pars) for (i=0; i pars[i].limits[1])) ) { info = MP_ERR_INITBOUNDS; goto CLEANUP; } if ( (pars[i].fixed == 0) && pars[i].limited[0] && pars[i].limited[1] && (pars[i].limits[0] >= pars[i].limits[1])) { info = MP_ERR_BOUNDS; goto CLEANUP; } } mp_malloc(qulim, int, nfree); mp_malloc(qllim, int, nfree); mp_malloc(ulim, double, nfree); mp_malloc(llim, double, nfree); for (i=0; i 0)) { ij = j*ldfjac; for (i=0; i= ulim[j])); int dwa1 = fabs(wa1[j]) > MP_MACHEP0; if (lpegged && (wa1[j] < 0)) wa1[j] = 0; if (upegged && (wa1[j] > 0)) wa1[j] = 0; if (dwa1 && qllim[j] && ((x[j] + wa1[j]) < llim[j])) { alpha = mp_dmin1(alpha, (llim[j]-x[j])/wa1[j]); } if (dwa1 && qulim[j] && ((x[j] + wa1[j]) > ulim[j])) { alpha = mp_dmin1(alpha, (ulim[j]-x[j])/wa1[j]); } } /* Scale the resulting vector, advance to the next position */ for (j=0; j= 0) ? (+1) : (-1); sgnl = (llim[j] >= 0) ? (+1) : (-1); ulim1 = ulim[j]*(1-sgnu*MP_MACHEP0) - ((ulim[j] == 0)?(MP_MACHEP0):0); llim1 = llim[j]*(1+sgnl*MP_MACHEP0) + ((llim[j] == 0)?(MP_MACHEP0):0); if (qulim[j] && (wa2[j] >= ulim1)) { wa2[j] = ulim[j]; } if (qllim[j] && (wa2[j] <= llim1)) { wa2[j] = llim[j]; } } } for (j=0; j= zero) { temp = p5; } else { temp = p5*dirder/(dirder + p5*actred); } if (((p1*fnorm1) >= fnorm) || (temp < p1) ) { temp = p1; } delta = temp*mp_dmin1(delta,pnorm/p1); par = par/temp; } else { if ((par == zero) || (ratio >= p75) ) { delta = pnorm/p5; par = p5*par; } } /* * test for successful iteration. */ if (ratio >= p0001) { /* * successful iteration. update x, fvec, and their norms. */ for (j=0; j 0) && (nfev >= conf.maxfev)) { /* Too many function evaluations */ info = MP_MAXITER; } if (iter >= conf.maxiter) { /* Too many iterations */ info = MP_MAXITER; } if ((fabs(actred) <= MP_MACHEP0) && (prered <= MP_MACHEP0) && (p5*ratio <= one) ) { info = MP_FTOL; } if (delta <= MP_MACHEP0*xnorm) { info = MP_XTOL; } if (gnorm <= MP_MACHEP0) { info = MP_GTOL; } if (info != 0) { goto L300; } /* * end of the inner loop. repeat if iteration unsuccessful. */ if (ratio < p0001) goto L200; /* * end of the outer loop. */ goto OUTER_LOOP; L300: /* * termination, either normal or user imposed. */ if (iflag < 0) { info = iflag; } iflag = 0; for (i=0; i 0) && (info > 0)) { iflag = mp_call(funct, m, npar, xall, fvec, 0, private); nfev += 1; } /* Compute number of pegged parameters */ npegged = 0; if (pars) for (i=0; icovar || result->xerror)) { mp_covar(nfree, fjac, ldfjac, ipvt, conf.covtol, wa2); if (result->covar) { /* Zero the destination covariance array */ for (j=0; j<(npar*npar); j++) result->covar[j] = 0; /* Transfer the covariance array */ for (j=0; jcovar[ifree[j]*npar+ifree[i]] = fjac[j*ldfjac+i]; } } } if (result->xerror) { for (j=0; jxerror[j] = 0; for (j=0; j 0) result->xerror[ifree[j]] = sqrt(cc); } } } if (result) { strcpy(result->version, MPFIT_VERSION); result->bestnorm = mp_dmax1(fnorm,fnorm1); result->bestnorm *= result->bestnorm; result->orignorm = orignorm; result->status = info; result->niter = iter; result->nfev = nfev; result->npar = npar; result->nfree = nfree; result->npegged = npegged; result->nfunc = m; /* Copy residuals if requested */ if (result->resid) { for (j=0; jresid[j] = fvec[j]; } } CLEANUP: if (fvec) free(fvec); if (qtf) free(qtf); if (x) free(x); if (xnew) free(xnew); if (fjac) free(fjac); if (diag) free(diag); if (wa1) free(wa1); if (wa2) free(wa2); if (wa3) free(wa3); if (wa4) free(wa4); if (ipvt) free(ipvt); if (pfixed) free(pfixed); if (step) free(step); if (dstep) free(dstep); if (mpside) free(mpside); if (ddebug) free(ddebug); if (ddrtol) free(ddrtol); if (ddatol) free(ddatol); if (ifree) free(ifree); if (qllim) free(qllim); if (qulim) free(qulim); if (llim) free(llim); if (ulim) free(ulim); return info; } /************************fdjac2.c*************************/ int mp_fdjac2(mp_func funct, int m, int n, int *ifree, int npar, double *x, double *fvec, double *fjac, int ldfjac, double epsfcn, double *wa, void *priv, int *nfev, double *step, double *dstep, int *dside, int *qulimited, double *ulimit, int *ddebug, double *ddrtol, double *ddatol) { /* * ********** * * subroutine fdjac2 * * this subroutine computes a forward-difference approximation * to the m by n jacobian matrix associated with a specified * problem of m functions in n variables. * * the subroutine statement is * * subroutine fdjac2(fcn,m,n,x,fvec,fjac,ldfjac,iflag,epsfcn,wa) * * where * * fcn is the name of the user-supplied subroutine which * calculates the functions. fcn must be declared * in an external statement in the user calling * program, and should be written as follows. * * subroutine fcn(m,n,x,fvec,iflag) * integer m,n,iflag * double precision x(n),fvec(m) * ---------- * calculate the functions at x and * return this vector in fvec. * ---------- * return * end * * the value of iflag should not be changed by fcn unless * the user wants to terminate execution of fdjac2. * in this case set iflag to a negative integer. * * m is a positive integer input variable set to the number * of functions. * * n is a positive integer input variable set to the number * of variables. n must not exceed m. * * x is an input array of length n. * * fvec is an input array of length m which must contain the * functions evaluated at x. * * fjac is an output m by n array which contains the * approximation to the jacobian matrix evaluated at x. * * ldfjac is a positive integer input variable not less than m * which specifies the leading dimension of the array fjac. * * iflag is an integer variable which can be used to terminate * the execution of fdjac2. see description of fcn. * * epsfcn is an input variable used in determining a suitable * step length for the forward-difference approximation. this * approximation assumes that the relative errors in the * functions are of the order of epsfcn. if epsfcn is less * than the machine precision, it is assumed that the relative * errors in the functions are of the order of the machine * precision. * * wa is a work array of length m. * * subprograms called * * user-supplied ...... fcn * * minpack-supplied ... dpmpar * * fortran-supplied ... dabs,dmax1,dsqrt * * argonne national laboratory. minpack project. march 1980. * burton s. garbow, kenneth e. hillstrom, jorge j. more * ********** */ int i,j,ij; int iflag = 0; double eps,h,temp; static double zero = 0.0; double **dvec = 0; int has_analytical_deriv = 0, has_numerical_deriv = 0; int has_debug_deriv = 0; temp = mp_dmax1(epsfcn,MP_MACHEP0); eps = sqrt(temp); ij = 0; dvec = (double **) malloc(sizeof(double **)*npar); if (dvec == 0) return MP_ERR_MEMORY; for (j=0; j 0) h = step[ifree[j]]; if (dstep && dstep[ifree[j]] > 0) h = fabs(dstep[ifree[j]]*temp); if (h == zero) h = eps; /* If negative step requested, or we are against the upper limit */ if ((dside && dsidei == -1) || (dside && dsidei == 0 && qulimited && ulimit && qulimited[j] && (temp > (ulimit[j]-h)))) { h = -h; } x[ifree[j]] = temp + h; iflag = mp_call(funct, m, npar, x, wa, 0, priv); if (nfev) *nfev = *nfev + 1; if (iflag < 0 ) goto DONE; x[ifree[j]] = temp; if (dsidei <= 1) { /* COMPUTE THE ONE-SIDED DERIVATIVE */ if (! debug) { /* Non-debug path for speed */ for (i=0; i da + fabs(fjold)*dr))) { printf(" %10d %10.4g %10.4g %10.4g %10.4g %10.4g\n", i, fvec[i], fjold, fjac[ij], fjold-fjac[ij], (fjold == 0)?(0):((fjold-fjac[ij])/fjold)); } } } } else { /* COMPUTE THE TWO-SIDED DERIVATIVE */ for (i=0; i da + fabs(fjold)*dr))) { printf(" %10d %10.4g %10.4g %10.4g %10.4g %10.4g\n", i, fvec[i], fjold, fjac[ij], fjold-fjac[ij], (fjold == 0)?(0):((fjold-fjac[ij])/fjold)); } } } } } if (has_debug_deriv) { printf("FJAC DEBUG END\n"); } DONE: if (dvec) free(dvec); if (iflag < 0) return iflag; return 0; /* * last card of subroutine fdjac2. */ } /************************qrfac.c*************************/ void mp_qrfac(int m, int n, double *a, int lda, int pivot, int *ipvt, int lipvt, double *rdiag, double *acnorm, double *wa) { /* * ********** * * subroutine qrfac * * this subroutine uses householder transformations with column * pivoting (optional) to compute a qr factorization of the * m by n matrix a. that is, qrfac determines an orthogonal * matrix q, a permutation matrix p, and an upper trapezoidal * matrix r with diagonal elements of nonincreasing magnitude, * such that a*p = q*r. the householder transformation for * column k, k = 1,2,...,min(m,n), is of the form * * t * i - (1/u(k))*u*u * * where u has zeros in the first k-1 positions. the form of * this transformation and the method of pivoting first * appeared in the corresponding linpack subroutine. * * the subroutine statement is * * subroutine qrfac(m,n,a,lda,pivot,ipvt,lipvt,rdiag,acnorm,wa) * * where * * m is a positive integer input variable set to the number * of rows of a. * * n is a positive integer input variable set to the number * of columns of a. * * a is an m by n array. on input a contains the matrix for * which the qr factorization is to be computed. on output * the strict upper trapezoidal part of a contains the strict * upper trapezoidal part of r, and the lower trapezoidal * part of a contains a factored form of q (the non-trivial * elements of the u vectors described above). * * lda is a positive integer input variable not less than m * which specifies the leading dimension of the array a. * * pivot is a logical input variable. if pivot is set true, * then column pivoting is enforced. if pivot is set false, * then no column pivoting is done. * * ipvt is an integer output array of length lipvt. ipvt * defines the permutation matrix p such that a*p = q*r. * column j of p is column ipvt(j) of the identity matrix. * if pivot is false, ipvt is not referenced. * * lipvt is a positive integer input variable. if pivot is false, * then lipvt may be as small as 1. if pivot is true, then * lipvt must be at least n. * * rdiag is an output array of length n which contains the * diagonal elements of r. * * acnorm is an output array of length n which contains the * norms of the corresponding columns of the input matrix a. * if this information is not needed, then acnorm can coincide * with rdiag. * * wa is a work array of length n. if pivot is false, then wa * can coincide with rdiag. * * subprograms called * * minpack-supplied ... dpmpar,enorm * * fortran-supplied ... dmax1,dsqrt,min0 * * argonne national laboratory. minpack project. march 1980. * burton s. garbow, kenneth e. hillstrom, jorge j. more * * ********** */ int i,ij,jj,j,jp1,k,kmax,minmn; double ajnorm,sum,temp; static double zero = 0.0; static double one = 1.0; static double p05 = 0.05; /* * compute the initial column norms and initialize several arrays. */ ij = 0; for (j=0; j rdiag[kmax]) kmax = k; } if (kmax == j) goto L40; ij = m * j; jj = m * kmax; for (i=0; i kp1) { ik = kk + 1; for (i=kp1; i jp1) { ij = jp1 + ldr * j; for (i=jp1; i= 1) { for (k=0; k= 0) { ij = ldr * j; for (i=0; i<=jm1; i++) { wa1[i] -= r[ij]*temp; ij += 1; } } } } for (j=0; j= n) { for (j=0; j= 0) { ij = jj; for (i=0; i<=jm1; i++) { sum += r[ij]*wa1[i]; ij += 1; } } wa1[j] = (wa1[j] - sum)/r[j+ldr*j]; jj += ldr; /* [i+ldr*j] */ } temp = mp_enorm(n,wa1); parl = ((fp/delta)/temp)/temp; } /* * calculate an upper bound, paru, for the zero of the function. */ jj = 0; for (j=0; j zero) parl = mp_dmax1(parl, *par); if (fp < zero) paru = mp_dmin1(paru, *par); /* * compute an improved estimate for par. */ *par = mp_dmax1(parl, *par + parc); /* * end of an iteration. */ goto L150; L220: /* * termination. */ if (iter == 0) *par = zero; /* * last card of subroutine lmpar. */ } /************************enorm.c*************************/ double mp_enorm(int n, double *x) { /* * ********** * * function enorm * * given an n-vector x, this function calculates the * euclidean norm of x. * * the euclidean norm is computed by accumulating the sum of * squares in three different sums. the sums of squares for the * small and large components are scaled so that no overflows * occur. non-destructive underflows are permitted. underflows * and overflows do not occur in the computation of the unscaled * sum of squares for the intermediate components. * the definitions of small, intermediate and large components * depend on two constants, rdwarf and rgiant. the main * restrictions on these constants are that rdwarf**2 not * underflow and rgiant**2 not overflow. the constants * given here are suitable for every known computer. * * the function statement is * * double precision function enorm(n,x) * * where * * n is a positive integer input variable. * * x is an input array of length n. * * subprograms called * * fortran-supplied ... dabs,dsqrt * * argonne national laboratory. minpack project. march 1980. * burton s. garbow, kenneth e. hillstrom, jorge j. more * * ********** */ int i; double agiant,floatn,s1,s2,s3,xabs,x1max,x3max; double ans, temp; double rdwarf = MP_RDWARF; double rgiant = MP_RGIANT; static double zero = 0.0; static double one = 1.0; s1 = zero; s2 = zero; s3 = zero; x1max = zero; x3max = zero; floatn = n; agiant = rgiant/floatn; for (i=0; i rdwarf) && (xabs < agiant)) { /* * sum for intermediate components. */ s2 += xabs*xabs; continue; } if (xabs > rdwarf) { /* * sum for large components. */ if (xabs > x1max) { temp = x1max/xabs; s1 = one + s1*temp*temp; x1max = xabs; } else { temp = xabs/x1max; s1 += temp*temp; } continue; } /* * sum for small components. */ if (xabs > x3max) { temp = x3max/xabs; s3 = one + s3*temp*temp; x3max = xabs; } else { if (xabs != zero) { temp = xabs/x3max; s3 += temp*temp; } } } /* * calculation of norm. */ if (s1 != zero) { temp = s1 + (s2/x1max)/x1max; ans = x1max*sqrt(temp); return(ans); } if (s2 != zero) { if (s2 >= x3max) temp = s2*(one+(x3max/s2)*(x3max*s3)); else temp = x3max*((s2/x3max)+(x3max*s3)); ans = sqrt(temp); } else { ans = x3max*sqrt(s3); } return(ans); /* * last card of function enorm. */ } /************************lmmisc.c*************************/ double mp_dmax1(double a, double b) { if (a >= b) return(a); else return(b); } double mp_dmin1(double a, double b) { if (a <= b) return(a); else return(b); } int mp_min0(int a, int b) { if (a <= b) return(a); else return(b); } /************************covar.c*************************/ /* c ********** c c subroutine covar c c given an m by n matrix a, the problem is to determine c the covariance matrix corresponding to a, defined as c c t c inverse(a *a) . c c this subroutine completes the solution of the problem c if it is provided with the necessary information from the c qr factorization, with column pivoting, of a. that is, if c a*p = q*r, where p is a permutation matrix, q has orthogonal c columns, and r is an upper triangular matrix with diagonal c elements of nonincreasing magnitude, then covar expects c the full upper triangle of r and the permutation matrix p. c the covariance matrix is then computed as c c t t c p*inverse(r *r)*p . c c if a is nearly rank deficient, it may be desirable to compute c the covariance matrix corresponding to the linearly independent c columns of a. to define the numerical rank of a, covar uses c the tolerance tol. if l is the largest integer such that c c abs(r(l,l)) .gt. tol*abs(r(1,1)) , c c then covar computes the covariance matrix corresponding to c the first l columns of r. for k greater than l, column c and row ipvt(k) of the covariance matrix are set to zero. c c the subroutine statement is c c subroutine covar(n,r,ldr,ipvt,tol,wa) c c where c c n is a positive integer input variable set to the order of r. c c r is an n by n array. on input the full upper triangle must c contain the full upper triangle of the matrix r. on output c r contains the square symmetric covariance matrix. c c ldr is a positive integer input variable not less than n c which specifies the leading dimension of the array r. c c ipvt is an integer input array of length n which defines the c permutation matrix p such that a*p = q*r. column j of p c is column ipvt(j) of the identity matrix. c c tol is a nonnegative input variable used to define the c numerical rank of a in the manner described above. c c wa is a work array of length n. c c subprograms called c c fortran-supplied ... dabs c c argonne national laboratory. minpack project. august 1980. c burton s. garbow, kenneth e. hillstrom, jorge j. more c c ********** */ int mp_covar(int n, double *r, int ldr, int *ipvt, double tol, double *wa) { int i, ii, j, jj, k, l; int kk, kj, ji, j0, k0, jj0; int sing; double one = 1.0, temp, tolr, zero = 0.0; /* * form the inverse of r in the full upper triangle of r. */ #if 0 for (j=0; j= 0) { for (k=0; k <= l; k++) { k0 = k*ldr; for (j=0; j l); j0 = j*ldr; jj0 = jj*ldr; for (i=0; i<=j; i++) { ji = j0+i; if (sing) r[ji] = zero; ii = ipvt[i]; if (ii > jj) r[jj0+ii] = r[ji]; if (ii < jj) r[ii*ldr+jj] = r[ji]; } wa[jj] = r[j0+j]; } /* * Symmetrize the covariance matrix in r */ for (j=0; j