/* augmented measurements */ static void LMBLEC_FUNC(LM_REAL *p, LM_REAL *hx, int m, int n, void *adata) { struct LMBLEC_DATA *data=(struct LMBLEC_DATA *)adata; int nn; register int i, j, *typ; register LM_REAL *lb, *ub, *w, tmp; nn=n-m; lb=data->lb; ub=data->ub; w=data->w; typ=data->bctype; (*(data->func))(p, hx, m, nn, data->adata); for(i=nn, j=0; i<n; ++i, ++j){ switch(typ[j]){ case __BC_INTERVAL__: tmp=(LM_CNST(2.0)*p[j]-(lb[j]+ub[j]))/(ub[j]-lb[j]); hx[i]=w[j]*__MAX__(tmp*tmp-LM_CNST(1.0), LM_CNST(0.0)); break; case __BC_LOW__: hx[i]=w[j]*__MAX__(lb[j]-p[j], LM_CNST(0.0)); break; case __BC_HIGH__: hx[i]=w[j]*__MAX__(p[j]-ub[j], LM_CNST(0.0)); break; } } }
/* blocked multiplication of the transpose of the nxm matrix a with itself (i.e. a^T a) * using a block size of bsize. The product is returned in b. * Since a^T a is symmetric, its computation can be speeded up by computing only its * upper triangular part and copying it to the lower part. * * More details on blocking can be found at * http://www-2.cs.cmu.edu/afs/cs/academic/class/15213-f02/www/R07/section_a/Recitation07-SectionA.pdf */ void TRANS_MAT_MAT_MULT(LM_REAL *a, LM_REAL *b, int n, int m) { #ifdef HAVE_LAPACK /* use BLAS matrix multiply */ LM_REAL alpha=CNST(1.0), beta=CNST(0.0); /* Fool BLAS to compute a^T*a avoiding transposing a: a is equivalent to a^T in column major, * therefore BLAS computes a*a^T with a and a*a^T in column major, which is equivalent to * computing a^T*a in row major! */ GEMM("N", "T", &m, &m, &n, &alpha, a, &m, a, &m, &beta, b, &m); #else /* no LAPACK, use blocking-based multiply */ register int i, j, k, jj, kk; register LM_REAL sum, *bim, *akm; const int bsize=__BLOCKSZ__; #define __MIN__(x, y) (((x)<=(y))? (x) : (y)) #define __MAX__(x, y) (((x)>=(y))? (x) : (y)) /* compute upper triangular part using blocking */ for(jj=0; jj<m; jj+=bsize){ for(i=0; i<m; ++i){ bim=b+i*m; for(j=__MAX__(jj, i); j<__MIN__(jj+bsize, m); ++j) bim[j]=0.0; //b[i*m+j]=0.0; } for(kk=0; kk<n; kk+=bsize){ for(i=0; i<m; ++i){ bim=b+i*m; for(j=__MAX__(jj, i); j<__MIN__(jj+bsize, m); ++j){ sum=0.0; for(k=kk; k<__MIN__(kk+bsize, n); ++k){ akm=a+k*m; sum+=akm[i]*akm[j]; //a[k*m+i]*a[k*m+j]; } bim[j]+=sum; //b[i*m+j]+=sum; } } } } /* copy upper triangular part to the lower one */ for(i=0; i<m; ++i) for(j=0; j<i; ++j) b[i*m+j]=b[j*m+i]; #undef __MIN__ #undef __MAX__ #endif /* HAVE_LAPACK */ }
/* check the supplied matlab function and its Jacobian. Returns 1 on error, 0 otherwise */ static int checkFuncAndJacobian(double *p, int m, int n, int havejac, struct mexdata *dat) { mxArray *lhs[1]; register int i; int ret=0; double *mp; mexSetTrapFlag(1); /* handle errors in the MEX-file */ mp=mxGetPr(dat->rhs[0]); for(i=0; i<m; ++i) mp[i]=p[i]; /* attempt to call the supplied func */ i=mexCallMATLAB(1, lhs, dat->nrhs, dat->rhs, dat->fname); if(i){ fprintf(stderr, "levmar: error calling '%s'.\n", dat->fname); ret=1; } else if(!mxIsDouble(lhs[0]) || mxIsComplex(lhs[0]) || !(mxGetM(lhs[0])==1 || mxGetN(lhs[0])==1) || __MAX__(mxGetM(lhs[0]), mxGetN(lhs[0]))!=n){ fprintf(stderr, "levmar: '%s' should produce a real vector with %d elements (got %d).\n", dat->fname, n, __MAX__(mxGetM(lhs[0]), mxGetN(lhs[0]))); ret=1; } /* delete the matrix created by matlab */ mxDestroyArray(lhs[0]); if(havejac){ /* attempt to call the supplied jac */ i=mexCallMATLAB(1, lhs, dat->nrhs, dat->rhs, dat->jacname); if(i){ fprintf(stderr, "levmar: error calling '%s'.\n", dat->jacname); ret=1; } else if(!mxIsDouble(lhs[0]) || mxIsComplex(lhs[0]) || mxGetM(lhs[0])!=n || mxGetN(lhs[0])!=m){ fprintf(stderr, "levmar: '%s' should produce a real %dx%d matrix (got %dx%d).\n", dat->jacname, n, m, mxGetM(lhs[0]), mxGetN(lhs[0])); ret=1; } else if(mxIsSparse(lhs[0])){ fprintf(stderr, "levmar: '%s' should produce a real dense matrix (got a sparse one).\n", dat->jacname); ret=1; } /* delete the matrix created by matlab */ mxDestroyArray(lhs[0]); } mexSetTrapFlag(0); /* on error terminate the MEX-file and return control to the MATLAB prompt */ return ret; }
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *Prhs[]) { register int i; register double *pdbl; mxArray **prhs=(mxArray **)&Prhs[0], *At, *Ct; struct mexdata mdata; int len, status; double *p, *p0, *ret, *x; int m, n, havejac, Arows, Crows, itmax, nopts, mintype, nextra; double opts[LM_OPTS_SZ]={LM_INIT_MU, LM_STOP_THRESH, LM_STOP_THRESH, LM_STOP_THRESH, LM_DIFF_DELTA}; double info[LM_INFO_SZ]; double *lb=NULL, *ub=NULL, *A=NULL, *b=NULL, *wghts=NULL, *C=NULL, *d=NULL, *covar=NULL; /* parse input args; start by checking their number */ if((nrhs<5)) matlabFmtdErrMsgTxt("levmar: at least 5 input arguments required (got %d).", nrhs); if(nlhs>4) matlabFmtdErrMsgTxt("levmar: too many output arguments (max. 4, got %d).", nlhs); else if(nlhs<2) matlabFmtdErrMsgTxt("levmar: too few output arguments (min. 2, got %d).", nlhs); /* note that in order to accommodate optional args, prhs & nrhs are adjusted accordingly below */ /** func **/ /* first argument must be a string , i.e. a char row vector */ if(mxIsChar(prhs[0])!=1) mexErrMsgTxt("levmar: first argument must be a string."); if(mxGetM(prhs[0])!=1) mexErrMsgTxt("levmar: first argument must be a string (i.e. char row vector)."); /* store supplied name */ len=mxGetN(prhs[0])+1; mdata.fname=mxCalloc(len, sizeof(char)); status=mxGetString(prhs[0], mdata.fname, len); if(status!=0) mexErrMsgTxt("levmar: not enough space. String is truncated."); /** jac (optional) **/ /* check whether second argument is a string */ if(mxIsChar(prhs[1])==1){ if(mxGetM(prhs[1])!=1) mexErrMsgTxt("levmar: second argument must be a string (i.e. row vector)."); /* store supplied name */ len=mxGetN(prhs[1])+1; mdata.jacname=mxCalloc(len, sizeof(char)); status=mxGetString(prhs[1], mdata.jacname, len); if(status!=0) mexErrMsgTxt("levmar: not enough space. String is truncated."); havejac=1; ++prhs; --nrhs; } else{ mdata.jacname=NULL; havejac=0; } #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: %s analytic Jacobian\n", havejac? "with" : "no"); #endif /* DEBUG */ /* CHECK if(!mxIsDouble(prhs[1]) || mxIsComplex(prhs[1]) || !(mxGetM(prhs[1])==1 && mxGetN(prhs[1])==1)) */ /** p0 **/ /* the second required argument must be a real row or column vector */ if(!mxIsDouble(prhs[1]) || mxIsComplex(prhs[1]) || !(mxGetM(prhs[1])==1 || mxGetN(prhs[1])==1)) mexErrMsgTxt("levmar: p0 must be a real vector."); p0=mxGetPr(prhs[1]); /* determine if we have a row or column vector and retrieve its * size, i.e. the number of parameters */ if(mxGetM(prhs[1])==1){ m=mxGetN(prhs[1]); mdata.isrow_p0=1; } else{ m=mxGetM(prhs[1]); mdata.isrow_p0=0; } /* copy input parameter vector to avoid destroying it */ p=mxMalloc(m*sizeof(double)); for(i=0; i<m; ++i) p[i]=p0[i]; /** x **/ /* the third required argument must be a real row or column vector */ if(!mxIsDouble(prhs[2]) || mxIsComplex(prhs[2]) || !(mxGetM(prhs[2])==1 || mxGetN(prhs[2])==1)) mexErrMsgTxt("levmar: x must be a real vector."); x=mxGetPr(prhs[2]); n=__MAX__(mxGetM(prhs[2]), mxGetN(prhs[2])); /** itmax **/ /* the fourth required argument must be a scalar */ if(!mxIsDouble(prhs[3]) || mxIsComplex(prhs[3]) || mxGetM(prhs[3])!=1 || mxGetN(prhs[3])!=1) mexErrMsgTxt("levmar: itmax must be a scalar."); itmax=(int)mxGetScalar(prhs[3]); /** opts **/ /* the fifth required argument must be a real row or column vector */ if(!mxIsDouble(prhs[4]) || mxIsComplex(prhs[4]) || (!(mxGetM(prhs[4])==1 || mxGetN(prhs[4])==1) && !(mxGetM(prhs[4])==0 && mxGetN(prhs[4])==0))) mexErrMsgTxt("levmar: opts must be a real vector."); pdbl=mxGetPr(prhs[4]); nopts=__MAX__(mxGetM(prhs[4]), mxGetN(prhs[4])); if(nopts!=0){ /* if opts==[], nothing needs to be done and the defaults are used */ if(nopts>LM_OPTS_SZ) matlabFmtdErrMsgTxt("levmar: opts must have at most %d elements, got %d.", LM_OPTS_SZ, nopts); else if(nopts<((havejac)? LM_OPTS_SZ-1 : LM_OPTS_SZ)) matlabFmtdWarnMsgTxt("levmar: only the %d first elements of opts specified, remaining set to defaults.", nopts); for(i=0; i<nopts; ++i) opts[i]=pdbl[i]; } #ifdef DEBUG else{ fflush(stderr); fprintf(stderr, "LEVMAR: empty options vector, using defaults\n"); } #endif /* DEBUG */ /** mintype (optional) **/ /* check whether sixth argument is a string */ if(nrhs>=6 && mxIsChar(prhs[5])==1 && mxGetM(prhs[5])==1){ char *minhowto; /* examine supplied name */ len=mxGetN(prhs[5])+1; minhowto=mxCalloc(len, sizeof(char)); status=mxGetString(prhs[5], minhowto, len); if(status!=0) mexErrMsgTxt("levmar: not enough space. String is truncated."); for(i=0; minhowto[i]; ++i) minhowto[i]=tolower(minhowto[i]); if(!strncmp(minhowto, "unc", 3)) mintype=MIN_UNCONSTRAINED; else if(!strncmp(minhowto, "bc", 2)) mintype=MIN_CONSTRAINED_BC; else if(!strncmp(minhowto, "lec", 3)) mintype=MIN_CONSTRAINED_LEC; else if(!strncmp(minhowto, "blec", 4)) mintype=MIN_CONSTRAINED_BLEC; else if(!strncmp(minhowto, "bleic", 5)) mintype=MIN_CONSTRAINED_BLEIC; else if(!strncmp(minhowto, "blic", 4)) mintype=MIN_CONSTRAINED_BLIC; else if(!strncmp(minhowto, "leic", 4)) mintype=MIN_CONSTRAINED_LEIC; else if(!strncmp(minhowto, "lic", 3)) mintype=MIN_CONSTRAINED_BLIC; else matlabFmtdErrMsgTxt("levmar: unknown minimization type '%s'.", minhowto); mxFree(minhowto); ++prhs; --nrhs; } else mintype=MIN_UNCONSTRAINED; if(mintype==MIN_UNCONSTRAINED) goto extraargs; /* arguments below this point are optional and their presence depends * upon the minimization type determined above */ /** lb, ub **/ if(nrhs>=7 && (mintype==MIN_CONSTRAINED_BC || mintype==MIN_CONSTRAINED_BLEC || mintype==MIN_CONSTRAINED_BLIC || mintype==MIN_CONSTRAINED_BLEIC)){ /* check if the next two arguments are real row or column vectors */ if(mxIsDouble(prhs[5]) && !mxIsComplex(prhs[5]) && (mxGetM(prhs[5])==1 || mxGetN(prhs[5])==1)){ if(mxIsDouble(prhs[6]) && !mxIsComplex(prhs[6]) && (mxGetM(prhs[6])==1 || mxGetN(prhs[6])==1)){ if((i=__MAX__(mxGetM(prhs[5]), mxGetN(prhs[5])))!=m) matlabFmtdErrMsgTxt("levmar: lb must have %d elements, got %d.", m, i); if((i=__MAX__(mxGetM(prhs[6]), mxGetN(prhs[6])))!=m) matlabFmtdErrMsgTxt("levmar: ub must have %d elements, got %d.", m, i); lb=mxGetPr(prhs[5]); ub=mxGetPr(prhs[6]); prhs+=2; nrhs-=2; } } } /** A, b **/ if(nrhs>=7 && (mintype==MIN_CONSTRAINED_LEC || mintype==MIN_CONSTRAINED_BLEC || mintype==MIN_CONSTRAINED_LEIC || mintype==MIN_CONSTRAINED_BLEIC)){ /* check if the next two arguments are a real matrix and a real row or column vector */ if(mxIsDouble(prhs[5]) && !mxIsComplex(prhs[5]) && mxGetM(prhs[5])>=1 && mxGetN(prhs[5])>=1){ if(mxIsDouble(prhs[6]) && !mxIsComplex(prhs[6]) && (mxGetM(prhs[6])==1 || mxGetN(prhs[6])==1)){ if((i=mxGetN(prhs[5]))!=m) matlabFmtdErrMsgTxt("levmar: A must have %d columns, got %d.", m, i); if((i=__MAX__(mxGetM(prhs[6]), mxGetN(prhs[6])))!=(Arows=mxGetM(prhs[5]))) matlabFmtdErrMsgTxt("levmar: b must have %d elements, got %d.", Arows, i); At=prhs[5]; b=mxGetPr(prhs[6]); A=getTranspose(At); prhs+=2; nrhs-=2; } } } /* wghts */ /* check if we have a weights vector */ if(nrhs>=6 && mintype==MIN_CONSTRAINED_BLEC){ /* only check if we have seen both box & linear constraints */ if(mxIsDouble(prhs[5]) && !mxIsComplex(prhs[5]) && (mxGetM(prhs[5])==1 || mxGetN(prhs[5])==1)){ if(__MAX__(mxGetM(prhs[5]), mxGetN(prhs[5]))==m){ wghts=mxGetPr(prhs[5]); ++prhs; --nrhs; } } } /** C, d **/ if(nrhs>=7 && (mintype==MIN_CONSTRAINED_BLEIC || mintype==MIN_CONSTRAINED_BLIC || mintype==MIN_CONSTRAINED_LEIC || mintype==MIN_CONSTRAINED_LIC)){ /* check if the next two arguments are a real matrix and a real row or column vector */ if(mxIsDouble(prhs[5]) && !mxIsComplex(prhs[5]) && mxGetM(prhs[5])>=1 && mxGetN(prhs[5])>=1){ if(mxIsDouble(prhs[6]) && !mxIsComplex(prhs[6]) && (mxGetM(prhs[6])==1 || mxGetN(prhs[6])==1)){ if((i=mxGetN(prhs[5]))!=m) matlabFmtdErrMsgTxt("levmar: C must have %d columns, got %d.", m, i); if((i=__MAX__(mxGetM(prhs[6]), mxGetN(prhs[6])))!=(Crows=mxGetM(prhs[5]))) matlabFmtdErrMsgTxt("levmar: d must have %d elements, got %d.", Crows, i); Ct=prhs[5]; d=mxGetPr(prhs[6]); C=getTranspose(Ct); prhs+=2; nrhs-=2; } } } /* arguments below this point are assumed to be extra arguments passed * to every invocation of the fitting function and its Jacobian */ extraargs: /* handle any extra args and allocate memory for * passing the current parameter estimate to matlab */ nextra=nrhs-5; mdata.nrhs=nextra+1; mdata.rhs=(mxArray **)mxMalloc(mdata.nrhs*sizeof(mxArray *)); for(i=0; i<nextra; ++i) mdata.rhs[i+1]=(mxArray *)prhs[nrhs-nextra+i]; /* discard 'const' modifier */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: %d extra args\n", nextra); #endif /* DEBUG */ if(mdata.isrow_p0){ /* row vector */ mdata.rhs[0]=mxCreateDoubleMatrix(1, m, mxREAL); /* mxSetM(mdata.rhs[0], 1); mxSetN(mdata.rhs[0], m); */ } else{ /* column vector */ mdata.rhs[0]=mxCreateDoubleMatrix(m, 1, mxREAL); /* mxSetM(mdata.rhs[0], m); mxSetN(mdata.rhs[0], 1); */ } /* ensure that the supplied function & Jacobian are as expected */ if(checkFuncAndJacobian(p, m, n, havejac, &mdata)){ status=LM_ERROR; goto cleanup; } if(nlhs>3) /* covariance output required */ covar=mxMalloc(m*m*sizeof(double)); /* invoke levmar */ switch(mintype){ case MIN_UNCONSTRAINED: /* no constraints */ if(havejac) status=dlevmar_der(func, jacfunc, p, x, m, n, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_dif(func, p, x, m, n, itmax, opts, info, NULL, covar, (void *)&mdata); #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_der()/dlevmar_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_BC: /* box constraints */ if(havejac) status=dlevmar_bc_der(func, jacfunc, p, x, m, n, lb, ub, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_bc_dif(func, p, x, m, n, lb, ub, itmax, opts, info, NULL, covar, (void *)&mdata); #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_bc_der()/dlevmar_bc_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_LEC: /* linear equation constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_lec_der(func, jacfunc, p, x, m, n, A, b, Arows, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_lec_dif(func, p, x, m, n, A, b, Arows, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no linear constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_lec_der()/dlevmar_lec_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_BLEC: /* box & linear equation constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_blec_der(func, jacfunc, p, x, m, n, lb, ub, A, b, Arows, wghts, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_blec_dif(func, p, x, m, n, lb, ub, A, b, Arows, wghts, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no box & linear constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_blec_der()/dlevmar_blec_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_BLEIC: /* box, linear equation & inequalities constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_bleic_der(func, jacfunc, p, x, m, n, lb, ub, A, b, Arows, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_bleic_dif(func, p, x, m, n, lb, ub, A, b, Arows, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no box, linear equation & inequality constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_bleic_der()/dlevmar_bleic_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_BLIC: /* box, linear inequalities constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_bleic_der(func, jacfunc, p, x, m, n, lb, ub, NULL, NULL, 0, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_bleic_dif(func, p, x, m, n, lb, ub, NULL, NULL, 0, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no box & linear inequality constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_blic_der()/dlevmar_blic_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_LEIC: /* linear equation & inequalities constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_bleic_der(func, jacfunc, p, x, m, n, NULL, NULL, A, b, Arows, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_bleic_dif(func, p, x, m, n, NULL, NULL, A, b, Arows, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no linear equation & inequality constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_leic_der()/dlevmar_leic_dif()\n"); #endif /* DEBUG */ break; case MIN_CONSTRAINED_LIC: /* linear inequalities constraints */ #ifdef HAVE_LAPACK if(havejac) status=dlevmar_bleic_der(func, jacfunc, p, x, m, n, NULL, NULL, NULL, NULL, 0, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); else status=dlevmar_bleic_dif(func, p, x, m, n, NULL, NULL, NULL, NULL, 0, C, d, Crows, itmax, opts, info, NULL, covar, (void *)&mdata); #else mexErrMsgTxt("levmar: no linear equation & inequality constraints support, HAVE_LAPACK was not defined during MEX-file compilation."); #endif /* HAVE_LAPACK */ #ifdef DEBUG fflush(stderr); fprintf(stderr, "LEVMAR: calling dlevmar_lic_der()/dlevmar_lic_dif()\n"); #endif /* DEBUG */ break; default: mexErrMsgTxt("levmar: unexpected internal error."); } #ifdef DEBUG fflush(stderr); printf("LEVMAR: minimization returned %d in %g iter, reason %g\n\tSolution: ", status, info[5], info[6]); for(i=0; i<m; ++i) printf("%.7g ", p[i]); printf("\n\n\tMinimization info:\n\t"); for(i=0; i<LM_INFO_SZ; ++i) printf("%g ", info[i]); printf("\n"); #endif /* DEBUG */ /* copy back return results */ /** ret **/ plhs[0]=mxCreateDoubleMatrix(1, 1, mxREAL); ret=mxGetPr(plhs[0]); ret[0]=(double)status; /** popt **/ plhs[1]=(mdata.isrow_p0==1)? mxCreateDoubleMatrix(1, m, mxREAL) : mxCreateDoubleMatrix(m, 1, mxREAL); pdbl=mxGetPr(plhs[1]); for(i=0; i<m; ++i) pdbl[i]=p[i]; /** info **/ if(nlhs>2){ plhs[2]=mxCreateDoubleMatrix(1, LM_INFO_SZ, mxREAL); pdbl=mxGetPr(plhs[2]); for(i=0; i<LM_INFO_SZ; ++i) pdbl[i]=info[i]; } /** covar **/ if(nlhs>3){ plhs[3]=mxCreateDoubleMatrix(m, m, mxREAL); pdbl=mxGetPr(plhs[3]); for(i=0; i<m*m; ++i) /* covariance matrices are symmetric, thus no need to transpose! */ pdbl[i]=covar[i]; } cleanup: /* cleanup */ mxDestroyArray(mdata.rhs[0]); if(A) mxFree(A); if(C) mxFree(C); mxFree(mdata.fname); if(havejac) mxFree(mdata.jacname); mxFree(p); mxFree(mdata.rhs); if(covar) mxFree(covar); if(status==LM_ERROR) mexWarnMsgTxt("levmar: optimization returned with an error!"); }