/* free workspace for demo3 */ static cs_long_t done3 (cs_long_t ok, cs_cls *S, cs_cln *N, cs_complex_t *y, cs_cl *W, cs_cl *E, cs_long_t *p) { cs_cl_sfree (S) ; cs_cl_nfree (N) ; cs_cl_free (y) ; cs_cl_spfree (W) ; cs_cl_spfree (E) ; cs_cl_free (p) ; return (ok) ; }
/* free a problem */ problem *free_problem (problem *Prob) { if (!Prob) return (NULL) ; cs_cl_spfree (Prob->A) ; if (Prob->sym) cs_cl_spfree (Prob->C) ; cs_cl_free (Prob->b) ; cs_cl_free (Prob->x) ; cs_cl_free (Prob->resid) ; return (cs_cl_free (Prob)) ; }
/* 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)") ; } } }
/* return a complex sparse matrix to MATLAB */ mxArray *cs_cl_mex_put_sparse (cs_cl **Ahandle) { cs_cl *A ; double *x, *z ; mxArray *Amatlab ; CS_INT k ; A = *Ahandle ; if (!A) mexErrMsgTxt ("failed") ; Amatlab = mxCreateSparse (0, 0, 0, mxCOMPLEX) ; mxSetM (Amatlab, A->m) ; mxSetN (Amatlab, A->n) ; mxSetNzmax (Amatlab, A->nzmax) ; cs_cl_free (mxGetJc (Amatlab)) ; cs_cl_free (mxGetIr (Amatlab)) ; cs_cl_free (mxGetPr (Amatlab)) ; cs_cl_free (mxGetPi (Amatlab)) ; mxSetJc (Amatlab, (void *) (A->p)) ; /* assign A->p pointer to MATLAB A */ mxSetIr (Amatlab, (void *) (A->i)) ; x = cs_dl_malloc (A->nzmax, sizeof (double)) ; z = cs_dl_malloc (A->nzmax, sizeof (double)) ; for (k = 0 ; k < A->nzmax ; k++) { x [k] = creal (A->x [k]) ; /* copy and split numerical values */ z [k] = cimag (A->x [k]) ; } cs_cl_free (A->x) ; /* free copy of complex values */ mxSetPr (Amatlab, x) ; mxSetPi (Amatlab, z) ; cs_cl_free (A) ; /* frees A struct only, not A->p, etc */ *Ahandle = NULL ; return (Amatlab) ; }
/* copy a complex vector back to MATLAB and free it */ mxArray *cs_cl_mex_put_double (CS_INT n, cs_complex_t *b) { double *x, *z ; mxArray *X ; CS_INT k ; X = mxCreateDoubleMatrix (n, 1, mxCOMPLEX) ; /* create x */ x = mxGetPr (X) ; z = mxGetPi (X) ; for (k = 0 ; k < n ; k++) { x [k] = creal (b [k]) ; /* copy x = b */ z [k] = cimag (b [k]) ; } cs_cl_free (b) ; return (X) ; }
/* cs_lu: sparse LU factorization, with optional fill-reducing ordering */ void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { CS_INT n, order, *p ; double tol ; if (nargout > 4 || nargin > 3 || nargin < 1) { mexErrMsgTxt ("Usage: [L,U,p,q] = cs_lu (A,tol)") ; } if (nargin == 2) /* determine tol and ordering */ { tol = mxGetScalar (pargin [1]) ; order = (nargout == 4) ? 1 : 0 ; /* amd (A+A'), or natural */ } else { tol = 1 ; order = (nargout == 4) ? 2 : 0 ; /* amd(S'*S) w/dense rows or I */ } if (mxIsComplex (pargin [0])) { #ifndef NCOMPLEX cs_cls *S ; cs_cln *N ; cs_cl Amatrix, *A, *D ; A = cs_cl_mex_get_sparse (&Amatrix, 1, pargin [0]) ; /* get A */ n = A->n ; S = cs_cl_sqr (order, A, 0) ; /* symbolic ordering, no QR bound */ N = cs_cl_lu (A, S, tol) ; /* numeric factorization */ if (!N) mexErrMsgTxt ("cs_lu failed (singular, or out of memory)") ; cs_cl_free (A->x) ; /* complex copy no longer needed */ cs_cl_dropzeros (N->L) ; /* drop zeros from L and sort it */ 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 U and sort it */ D = cs_cl_transpose (N->U, 1) ; cs_cl_spfree (N->U) ; N->U = cs_cl_transpose (D, 1) ; cs_cl_spfree (D) ; p = cs_cl_pinv (N->pinv, n) ; /* p=pinv' */ pargout [0] = cs_cl_mex_put_sparse (&(N->L)) ; /* return L */ pargout [1] = cs_cl_mex_put_sparse (&(N->U)) ; /* return U */ pargout [2] = cs_dl_mex_put_int (p, n, 1, 1) ; /* return p */ /* return Q */ if (nargout == 4) pargout [3] = cs_dl_mex_put_int (S->q, n, 1, 0) ; 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, 1, 1, pargin [0]) ; /* get A */ n = A->n ; S = cs_dl_sqr (order, A, 0) ; /* symbolic ordering, no QR bound */ N = cs_dl_lu (A, S, tol) ; /* numeric factorization */ if (!N) mexErrMsgTxt ("cs_lu failed (singular, or out of memory)") ; cs_dl_dropzeros (N->L) ; /* drop zeros from L and sort it */ 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 U and sort it */ D = cs_dl_transpose (N->U, 1) ; cs_dl_spfree (N->U) ; N->U = cs_dl_transpose (D, 1) ; cs_dl_spfree (D) ; p = cs_dl_pinv (N->pinv, n) ; /* p=pinv' */ pargout [0] = cs_dl_mex_put_sparse (&(N->L)) ; /* return L */ pargout [1] = cs_dl_mex_put_sparse (&(N->U)) ; /* return U */ pargout [2] = cs_dl_mex_put_int (p, n, 1, 1) ; /* return p */ /* return Q */ if (nargout == 4) pargout [3] = cs_dl_mex_put_int (S->q, n, 1, 0) ; cs_dl_nfree (N) ; cs_dl_sfree (S) ; } }