/* cs_lusol: solve A*x=b using a sparse LU factorization */ void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { double tol ; CS_INT order ; if (nargout > 1 || nargin < 2 || nargin > 4) { mexErrMsgTxt ("Usage: x = cs_lusol(A,b,order,tol)") ; } order = (nargin < 3) ? 2 : mxGetScalar (pargin [2]) ; order = CS_MAX (order, 0) ; order = CS_MIN (order, 3) ; if (nargin == 2) { tol = 1 ; /* normal partial pivoting */ } else if (nargin == 3) { tol = (order == 1) ? 0.001 : 1 ; /* tol = 0.001 for amd(A+A') */ } else { tol = mxGetScalar (pargin [3]) ; } if (mxIsComplex (pargin [0]) || mxIsComplex (pargin [1])) { #ifndef NCOMPLEX cs_cl *A, Amatrix ; cs_complex_t *x ; A = cs_cl_mex_get_sparse (&Amatrix, 1, pargin [0]) ; /* get A */ x = cs_cl_mex_get_double (A->n, pargin [1]) ; /* x = b */ if (!cs_cl_lusol (order, A, x, tol)) /* x = A\x */ { mexErrMsgTxt ("failed (singular or out of memory)") ; } cs_cl_free (A->x) ; /* complex copy no longer needed */ pargout [0] = cs_cl_mex_put_double (A->n, x) ; /* return x */ #else mexErrMsgTxt ("complex matrices not supported") ; #endif } else { cs_dl *A, Amatrix ; double *x, *b ; A = cs_dl_mex_get_sparse (&Amatrix, 1, 1, pargin [0]) ; /* get A */ b = cs_dl_mex_get_double (A->n, pargin [1]) ; /* get b */ x = cs_dl_mex_put_double (A->n, b, &(pargout [0])) ; /* x = b */ if (!cs_dl_lusol (order, A, x, tol)) /* x = A\x */ { mexErrMsgTxt ("failed (singular or out of memory)") ; } } }
/* cs_qrsol: solve least squares or underdetermined problem */ void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { CS_INT k, order ; if (nargout > 1 || nargin < 2 || nargin > 3) { mexErrMsgTxt ("Usage: x = cs_qrsol(A,b,order)") ; } order = (nargin < 3) ? 3 : mxGetScalar (pargin [2]) ; order = CS_MAX (order, 0) ; order = CS_MIN (order, 3) ; if (mxIsComplex (pargin [0]) || mxIsComplex (pargin [1])) { #ifndef NCOMPLEX cs_cl *A, Amatrix ; cs_complex_t *x, *b ; A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ; /* get A */ b = cs_cl_mex_get_double (A->m, pargin [1]) ; /* get b */ x = cs_dl_calloc (CS_MAX (A->m, A->n), sizeof (cs_complex_t)) ; for (k = 0 ; k < A->m ; k++) x [k] = b [k] ; /* x = b */ cs_free (b) ; if (!cs_cl_qrsol (order, A, x)) /* x = A\x */ { mexErrMsgTxt ("QR solve failed") ; } pargout [0] = cs_cl_mex_put_double (A->n, x) ; /* return x */ #else mexErrMsgTxt ("complex matrices not supported") ; #endif } else { cs_dl *A, Amatrix ; double *x, *b ; A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */ b = cs_dl_mex_get_double (A->m, pargin [1]) ; /* get b */ x = cs_dl_calloc (CS_MAX (A->m, A->n), sizeof (double)) ; /* x = b */ for (k = 0 ; k < A->m ; k++) x [k] = b [k] ; if (!cs_dl_qrsol (order, A, x)) /* x = A\x */ { mexErrMsgTxt ("QR solve failed") ; } cs_dl_mex_put_double (A->n, x, &(pargout [0])) ; /* return x */ cs_free (x) ; } }
/* cs_sqr: symbolic sparse QR factorization */ void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { double s ; cs_dls *S ; cs_dl Amatrix, *A ; CS_INT m, n, order, *p ; if (nargout > 7 || nargin != 1) { mexErrMsgTxt ("Usage: [vnz,rnz,parent,c,leftmost,p,q] = cs_sqr(A)") ; } A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */ m = A->m ; n = A->n ; if (m < n) mexErrMsgTxt ("A must have # rows >= # columns") ; order = (nargout == 7) ? 3 : 0 ; /* determine ordering */ S = cs_dl_sqr (order, A, 1) ; /* symbolic QR ordering & analysis*/ if (!S) mexErrMsgTxt ("cs_sqr failed") ; s = S->lnz ; cs_dl_mex_put_double (1, &s, &(pargout [0])) ; /* return nnz(V) */ s = S->unz ; cs_dl_mex_put_double (1, &s, &(pargout [1])) ; /* return nnz(R) */ pargout [2] = cs_dl_mex_put_int (S->parent, n, 1, 0) ; /* return parent */ pargout [3] = cs_dl_mex_put_int (S->cp, n, 0, 0) ; /* return c */ pargout [4] = cs_dl_mex_put_int (S->leftmost, m, 1, 0) ;/* return leftmost*/ p = cs_dl_pinv (S->pinv, S->m2) ; /* p = pinv' */ pargout [5] = cs_dl_mex_put_int (p, S->m2, 1, 1) ; /* return p */ if (nargout > 6) { pargout [6] = cs_dl_mex_put_int (S->q, n, 1, 0) ; /* return q */ } cs_dl_sfree (S) ; }
/* cs_qr: sparse QR factorization */ void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { CS_INT m, n, order, *p ; if (nargout > 5 || nargin != 1) { mexErrMsgTxt ("Usage: [V,beta,p,R,q] = cs_qr(A)") ; } order = (nargout == 5) ? 3 : 0 ; /* determine ordering */ m = mxGetM (pargin [0]) ; n = mxGetN (pargin [0]) ; if (m < n) mexErrMsgTxt ("A must have # rows >= # columns") ; if (mxIsComplex (pargin [0])) { #ifndef NCOMPLEX cs_cls *S ; cs_cln *N ; cs_cl Amatrix, *A, *D ; A = cs_cl_mex_get_sparse (&Amatrix, 0, pargin [0]) ; /* get A */ S = cs_cl_sqr (order, A, 1) ; /* symbolic QR ordering & analysis*/ N = cs_cl_qr (A, S) ; /* numeric QR factorization */ cs_free (A->x) ; if (!N) mexErrMsgTxt ("qr failed") ; cs_cl_dropzeros (N->L) ; /* drop zeros from V and sort */ D = cs_cl_transpose (N->L, 1) ; cs_cl_spfree (N->L) ; N->L = cs_cl_transpose (D, 1) ; cs_cl_spfree (D) ; cs_cl_dropzeros (N->U) ; /* drop zeros from R and sort */ D = cs_cl_transpose (N->U, 1) ; cs_cl_spfree (N->U) ; N->U = cs_cl_transpose (D, 1) ; cs_cl_spfree (D) ; m = N->L->m ; /* m may be larger now */ p = cs_cl_pinv (S->pinv, m) ; /* p = pinv' */ pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ; /* return V */ cs_dl_mex_put_double (n, N->B, &(pargout [1])) ; /* return beta */ pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ; /* return p */ pargout [3] = cs_cl_mex_put_sparse (&(N->U)) ; /* return R */ pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ; /* return q */ cs_cl_nfree (N) ; cs_cl_sfree (S) ; #else mexErrMsgTxt ("complex matrices not supported") ; #endif } else { cs_dls *S ; cs_dln *N ; cs_dl Amatrix, *A, *D ; A = cs_dl_mex_get_sparse (&Amatrix, 0, 1, pargin [0]) ; /* get A */ S = cs_dl_sqr (order, A, 1) ; /* symbolic QR ordering & analysis*/ N = cs_dl_qr (A, S) ; /* numeric QR factorization */ if (!N) mexErrMsgTxt ("qr failed") ; cs_dl_dropzeros (N->L) ; /* drop zeros from V and sort */ D = cs_dl_transpose (N->L, 1) ; cs_dl_spfree (N->L) ; N->L = cs_dl_transpose (D, 1) ; cs_dl_spfree (D) ; cs_dl_dropzeros (N->U) ; /* drop zeros from R and sort */ D = cs_dl_transpose (N->U, 1) ; cs_dl_spfree (N->U) ; N->U = cs_dl_transpose (D, 1) ; cs_dl_spfree (D) ; m = N->L->m ; /* m may be larger now */ p = cs_dl_pinv (S->pinv, m) ; /* p = pinv' */ pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ; /* return V */ cs_dl_mex_put_double (n, N->B, &(pargout [1])) ; /* return beta */ pargout [2] = cs_dl_mex_put_int (p, m, 1, 1) ; /* return p */ pargout [3] = cs_dl_mex_put_sparse (&(N->U)) ; /* return R */ pargout [4] = cs_dl_mex_put_int (S->q, n, 1, 0) ; /* return q */ cs_dl_nfree (N) ; cs_dl_sfree (S) ; } }