/* 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) ; }
mxArray *mxCreateNumericMatrixE(mwSize m, mwSize n, mxClassID classid, mxComplexity ComplexFlag) { mwSize sz = m*n*sizeof(double); mxArray *a = mxCreateNumericMatrix(1, 1, classid, ComplexFlag); mxSetM(a,m); mxSetN(a,n); mxSetPr(a, mxRealloc(mxGetPr(a),sz)); if(ComplexFlag == mxCOMPLEX) { mxSetPi(a, mxRealloc(mxGetPi(a),sz)); } return a; }
static int put_values ( Long nz, mxArray *A, double *Ax, // complex case: size 2*nz and freed on return, // real case: size nz, not freed on return. Long is_complex, cholmod_common *cc ) { Long imag_all_zero = TRUE ; if (is_complex) { // A is complex, stored in interleaved form; split it for MATLAB Long k ; double z, *Ax2, *Az2 ; mxFree (mxGetPi (A)) ; Ax2 = (double *) cholmod_l_malloc (nz, sizeof (double), cc) ; Az2 = (double *) cholmod_l_malloc (nz, sizeof (double), cc) ; for (k = 0 ; k < nz ; k++) { Ax2 [k] = Ax [2*k] ; z = Ax [2*k+1] ; if (z != 0) { imag_all_zero = FALSE ; } Az2 [k] = z ; } mxSetPr (A, Ax2) ; if (imag_all_zero) { // free the imaginary part, converting A to real cholmod_l_free (nz, sizeof (double), Az2, cc) ; Az2 = NULL ; } mxSetPi (A, Az2) ; // NOTE: the input Ax is freed cholmod_l_free (nz, sizeof (Complex), Ax, cc) ; } else { // A is real; just set Ax and return (do not free Ax) mxSetPr (A, Ax) ; } return (TRUE) ; }
mxArray *sfmult_yalloc // return Y ( Int m, Int n, int Ycomplex // true if Y is complex ) { // (TO DO): guard against integer overflow mxArray *Y = mxCreateDoubleMatrix (0, 0, Ycomplex ? mxCOMPLEX : mxREAL) ; MXFREE (mxGetPr (Y)) ; mxSetPr (Y, mxMalloc (MAX (m*n, 1) * sizeof (double))) ; if (Ycomplex) { MXFREE (mxGetPi (Y)) ; mxSetPi (Y, mxMalloc (MAX (m*n, 1) * sizeof (double))) ; } mxSetM (Y, m) ; mxSetN (Y, n) ; return (Y) ; }
mxArray *mxCreateNumericArrayE(mwSize ndim, const mwSize *dims, mxClassID classid, mxComplexity ComplexFlag) { mxArray *a; mwSize i; mwSize *dims1 = mxMalloc(ndim*sizeof(mwSize)); mwSize sz = 1; for(i=0;i<ndim;i++) { sz *= dims[i]; dims1[i] = 1; } a = mxCreateNumericArray(ndim,dims1,classid,ComplexFlag); sz *= mxGetElementSize(a); mxSetDimensions(a, dims, ndim); mxFree(dims1); mxSetData(a, mxRealloc(mxGetData(a), sz)); if(ComplexFlag == mxCOMPLEX) { mxSetPi(a, mxRealloc(mxGetPi(a),sz)); } return a; }
/* ************************************************************ PROCEDURE mexFunction - Entry for Matlab [lab,q] = eigK(x,K) Computes spectral coefficients of x w.r.t. K REMARK If this function is used internally by SeDuMi, then complex numbers are stored in a single real vector. To make it invokable from the Matlab command-line by the user, we also allow Matlab complex vector x. ************************************************************ */ void mexFunction(const int nlhs, mxArray *plhs[], const int nrhs, const mxArray *prhs[]) { mxArray *output_array[3], *Xk, *hXk; coneK cK; int k, nk, nksqr, lendiag,i,ii,nkp1, lenfull; double *lab,*q,*qpi,*labk,*xwork,*xpiwork; const double *x,*xpi; /* ------------------------------------------------------------ Check for proper number of arguments ------------------------------------------------------------ */ mxAssert(nrhs >= NPARIN, "eigK requires more input arguments"); mxAssert(nlhs <= NPAROUT, "eigK produces less output arguments"); /* ------------------------------------------------------------ Disassemble cone K structure ------------------------------------------------------------ */ conepars(K_IN, &cK); /* ------------------------------------------------------------ Compute statistics based on cone K structure ------------------------------------------------------------ */ lendiag = cK.lpN + 2 * (cK.lorN + cK.rconeN) + cK.rLen + cK.hLen; lenfull = cK.lpN + cK.qDim + cK.rDim + cK.hDim; if(cK.rconeN > 0) for(i = 0; i < cK.rconeN; i++) lenfull += cK.rconeNL[i]; /* ------------------------------------------------------------ Get input vector x ------------------------------------------------------------ */ mxAssert(mxGetM(X_IN) * mxGetN(X_IN) == lenfull, "Size mismatch x"); mxAssert(!mxIsSparse(X_IN), "x must be full (not sparse)."); x = mxGetPr(X_IN); if(mxIsComplex(X_IN)) xpi = mxGetPi(X_IN) + cK.lpN; /* ------------------------------------------------------------ Allocate output LAB(diag), eigvec Q(full for psd) ------------------------------------------------------------ */ LAB_OUT = mxCreateDoubleMatrix(lendiag, 1, mxREAL); lab = mxGetPr(LAB_OUT); if(nlhs > 1){ if(mxIsComplex(X_IN)){ Q_OUT = mxCreateDoubleMatrix(cK.rDim, 1, mxCOMPLEX); qpi = mxGetPi(Q_OUT); } else Q_OUT = mxCreateDoubleMatrix(cK.rDim + cK.hDim, 1, mxREAL); q = mxGetPr(Q_OUT); } /* ------------------------------------------------------------ Allocate working arrays: ------------------------------------------------------------ */ Xk = mxCreateDoubleMatrix(0,0,mxREAL); hXk = mxCreateDoubleMatrix(0,0,mxCOMPLEX); if(mxIsComplex(X_IN)){ xwork = (double *) mxCalloc(MAX(1,2 * SQR(cK.rMaxn)), sizeof(double)); xpiwork = xwork + SQR(cK.rMaxn); } else xwork =(double *) mxCalloc(MAX(1,SQR(cK.rMaxn)+2*SQR(cK.hMaxn)), sizeof(double)); /* ------------------------------------------------------------ The actual job is done here:. ------------------------------------------------------------ */ if(cK.lpN){ /* ------------------------------------------------------------ LP: lab = x ------------------------------------------------------------ */ memcpy(lab, x, cK.lpN * sizeof(double)); lab += cK.lpN; x += cK.lpN; } /* ------------------------------------------------------------ CONSIDER FIRST MATLAB-REAL-TYPE: ------------------------------------------------------------ */ if(!mxIsComplex(X_IN)){ /* Not Matlab-type complex */ /* ------------------------------------------------------------ LORENTZ: (I) lab = qeig(x) ------------------------------------------------------------ */ for(k = 0; k < cK.lorN; k++){ nk = cK.lorNL[k]; qeig(lab,x,nk); lab += 2; x += nk; } /* ------------------------------------------------------------ RCONE: LAB = eig(X) (Lorentz-Rcone's are not used internally) ------------------------------------------------------------ */ for(k = 0; k < cK.rconeN; k++){ nk = cK.rconeNL[k]; rconeeig(lab,x[0],x[1],realssqr(x+2,nk-2)); lab += 2; x += nk; } /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ if(nlhs < 2){ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(1, output_array, 1, &Xk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); /* ------------------------------------------------------------ With mexCallMATLAB, we invoked the mexFunction "eig", which allocates a matrix struct *output_array[0], AND a block for the float data of that matrix. ==> mxDestroyArray() does not only free the float data, it also releases the matrix struct (and this is what we want). ------------------------------------------------------------ */ mxDestroyArray(output_array[0]); lab += nk; x += SQR(nk); } /* ------------------------------------------------------------ WARNING: Matlab's eig doesn't recognize Hermitian, hence VERY slow ------------------------------------------------------------ */ for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork + nksqr); mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); mxDestroyArray(output_array[0]); lab += nk; x += 2 * nksqr; } } else{ /* ------------------------------------------------------------ SDP: (II) (Q,LAB) = eig(X) ------------------------------------------------------------ */ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(2, output_array, 1, &Xk, "eig"); nksqr = SQR(nk); /* copy Q-matrix */ memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += nksqr; q += nksqr; } for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork+nksqr); mexCallMATLAB(2, output_array, 1, &hXk, "eig"); memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); q += nksqr; if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(q, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += 2 * nksqr; q += nksqr; } } /* [lab,q] = eigK */ } /* !iscomplex */ else{ /* is MATLAB type complex */ /* ------------------------------------------------------------ LORENTZ: (I) lab = qeig(x) ------------------------------------------------------------ */ for(k = 0; k < cK.lorN; k++){ nk = cK.lorNL[k]; cxqeig(lab,x,xpi,nk); lab += 2; x += nk; xpi += nk; } /* ------------------------------------------------------------ RCONE: LAB = eig(X) (Lorentz-Rcone's are not used internally) ------------------------------------------------------------ */ for(k = 0; k < cK.rconeN; k++){ nk = cK.rconeNL[k]; rconeeig(lab,x[0],x[1], realssqr(x+2,nk-2) + realssqr(xpi+2,nk-2)); lab += 2; x += nk; xpi += nk; } /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ for(k = 0; k < cK.sdpN; k++){ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xpiwork,xpi,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xpiwork); if(nlhs < 2){ mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); } else{ mexCallMATLAB(2, output_array, 1, &hXk, "eig"); memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(qpi, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[1]); q += nksqr; qpi += nksqr; } mxDestroyArray(output_array[0]); lab += nk; x += nksqr; xpi += nksqr; } } /* iscomplex */ /* ------------------------------------------------------------ Release PSD-working arrays. ------------------------------------------------------------ */ mxSetM(Xk,0); mxSetN(Xk,0); mxSetPr(Xk, (double *) NULL); mxDestroyArray(Xk); mxSetM(hXk,0); mxSetN(hXk,0); mxSetPr(hXk, (double *) NULL); mxSetPi(hXk, (double *) NULL); mxDestroyArray(hXk); mxFree(xwork); }
mxArray *ssmult_transpose // returns C = A' or A.' ( const mxArray *A, int conj // compute A' if true, compute A.' if false ) { Int *Cp, *Ci, *Ap, *Ai, *W ; double *Cx, *Cz, *Ax, *Az ; // (TO DO): do single too mxArray *C ; Int p, pend, q, i, j, n, m, anz, cnz ; int C_is_complex ; //-------------------------------------------------------------------------- // get inputs //-------------------------------------------------------------------------- m = mxGetM (A) ; n = mxGetN (A) ; Ap = mxGetJc (A) ; Ai = mxGetIr (A) ; Ax = mxGetPr (A) ; Az = mxGetPi (A) ; anz = Ap [n] ; C_is_complex = mxIsComplex (A) ; //-------------------------------------------------------------------------- // allocate C but do not initialize it //-------------------------------------------------------------------------- cnz = MAX (anz, 1) ; C = mxCreateSparse (0, 0, 0, C_is_complex ? mxCOMPLEX : mxREAL) ; MXFREE (mxGetJc (C)) ; MXFREE (mxGetIr (C)) ; MXFREE (mxGetPr (C)) ; MXFREE (mxGetPi (C)) ; Cp = mxMalloc ((m+1) * sizeof (Int)) ; Ci = mxMalloc (MAX (cnz,1) * sizeof (Int)) ; Cx = mxMalloc (MAX (cnz,1) * sizeof (double)) ; Cz = C_is_complex ? mxMalloc (MAX (cnz,1) * sizeof (double)) : NULL ; mxSetJc (C, Cp) ; mxSetIr (C, Ci) ; mxSetPr (C, Cx) ; mxSetPi (C, Cz) ; mxSetNzmax (C, cnz) ; mxSetM (C, n) ; mxSetN (C, m) ; //-------------------------------------------------------------------------- // allocate workspace //-------------------------------------------------------------------------- W = mxCalloc (MAX (m,1), sizeof (Int)) ; //-------------------------------------------------------------------------- // compute row counts //-------------------------------------------------------------------------- for (p = 0 ; p < anz ; p++) { W [Ai [p]]++ ; } //-------------------------------------------------------------------------- // compute column pointers of C and copy back into W //-------------------------------------------------------------------------- for (p = 0, i = 0 ; i < m ; i++) { Cp [i] = p ; p += W [i] ; W [i] = Cp [i] ; } Cp [m] = p ; //-------------------------------------------------------------------------- // C = A' //-------------------------------------------------------------------------- p = 0 ; if (!C_is_complex) { // C = A' (real case) for (j = 0 ; j < n ; j++) { pend = Ap [j+1] ; for ( ; p < pend ; p++) { q = W [Ai [p]]++ ; // find position for C(j,i) Ci [q] = j ; // place A(i,j) as entry C(j,i) Cx [q] = Ax [p] ; } } } else if (conj) { // C = A' (complex conjugate) for (j = 0 ; j < n ; j++) { pend = Ap [j+1] ; for ( ; p < pend ; p++) { q = W [Ai [p]]++ ; // find position for C(j,i) Ci [q] = j ; // place A(i,j) as entry C(j,i) Cx [q] = Ax [p] ; Cz [q] = -Az [p] ; } } } else { // C = A.' (complex case) for (j = 0 ; j < n ; j++) { pend = Ap [j+1] ; for ( ; p < pend ; p++) { q = W [Ai [p]]++ ; // find position for C(j,i) Ci [q] = j ; // place A(i,j) as entry C(j,i) Cx [q] = Ax [p] ; Cz [q] = Az [p] ; } } } //-------------------------------------------------------------------------- // free workspace and return result //-------------------------------------------------------------------------- MXFREE (W) ; return (C) ; }
/* ************************************************************ PROCEDURE mexFunction - Entry for Matlab ************************************************************ */ void mexFunction(const int nlhs, mxArray *plhs[], const int nrhs, const mxArray *prhs[]) { mxArray *output_array[3], *Xk, *hXk; coneK cK; int k, nk, nkp1, nksqr, lendiag,lenud, lenfull, i,ii; double *lab,*q,*labk,*xwork; const double *x; /* ------------------------------------------------------------ Check for proper number of arguments ------------------------------------------------------------ */ mxAssert(nrhs >= NPARIN, "psdeig requires more input arguments"); mxAssert(nlhs <= NPAROUT, "psdeig produces less output arguments"); /* ------------------------------------------------------------ Disassemble cone K structure ------------------------------------------------------------ */ conepars(K_IN, &cK); /* ------------------------------------------------------------ Compute statistics based on cone K structure ------------------------------------------------------------ */ lendiag = cK.rLen + cK.hLen; lenud = cK.rDim + cK.hDim; lenfull = cK.lpN + cK.qDim + lenud; /* ------------------------------------------------------------ Get input vector x ------------------------------------------------------------ */ mxAssert(!mxIsSparse(X_IN), "x must be full (not sparse)."); x = mxGetPr(X_IN); if(mxGetM(X_IN) * mxGetN(X_IN) != lenud){ mxAssert(mxGetM(X_IN) * mxGetN(X_IN) == lenfull, "Size mismatch x"); x += cK.lpN + cK.qDim; /* point to PSD part */ } /* ------------------------------------------------------------ Allocate output LAB(diag), eigvec Q(full for psd) ------------------------------------------------------------ */ LAB_OUT = mxCreateDoubleMatrix(lendiag, 1, mxREAL); lab = mxGetPr(LAB_OUT); if(nlhs > 1){ Q_OUT = mxCreateDoubleMatrix(lenud, 1, mxREAL); q = mxGetPr(Q_OUT); } /* ------------------------------------------------------------ Allocate working arrays: ------------------------------------------------------------ */ Xk = mxCreateDoubleMatrix(0,0,mxREAL); hXk = mxCreateDoubleMatrix(0,0,mxCOMPLEX); xwork =(double *) mxCalloc(MAX(1,SQR(cK.rMaxn)+2*SQR(cK.hMaxn)), sizeof(double)); /* ------------------------------------------------------------ PSD: (I) LAB = eig(X) ------------------------------------------------------------ */ if(nlhs < 2){ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(1, output_array, 1, &Xk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); /* ------------------------------------------------------------ With mexCallMATLAB, we invoked the mexFunction "eig", which allocates a matrix struct *output_array[0], AND a block for the float data of that matrix. ==> mxDestroyArray() does not only free the float data, it also releases the matrix struct (and this is what we want). ------------------------------------------------------------ */ mxDestroyArray(output_array[0]); lab += nk; x += SQR(nk); } /* ------------------------------------------------------------ WARNING: Matlab's eig doesn't recognize Hermitian, hence VERY slow ------------------------------------------------------------ */ for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork + nksqr); mexCallMATLAB(1, output_array, 1, &hXk, "eig"); memcpy(lab, mxGetPr(output_array[0]), nk * sizeof(double)); mxDestroyArray(output_array[0]); lab += nk; x += 2 * nksqr; } } /* nlhs < 2 */ else{ /* ------------------------------------------------------------ SDP: (II) (Q,LAB) = eig(X) ------------------------------------------------------------ */ for(k=0; k < cK.rsdpN; k++){ /* real symmetric */ nk = cK.sdpNL[k]; symproj(xwork,x,nk); /* make it symmetric */ mxSetM(Xk, nk); mxSetN(Xk, nk); mxSetPr(Xk, xwork); mexCallMATLAB(2, output_array, 1, &Xk, "eig"); nksqr = SQR(nk); /* copy Q-matrix */ memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); lab += nk; x += nksqr; q += nksqr; } for(; k < cK.sdpN; k++){ /* complex Hermitian */ nk = cK.sdpNL[k]; nksqr = SQR(nk); symproj(xwork,x,nk); /* make it Hermitian */ skewproj(xwork + nksqr,x+nksqr,nk); mxSetM(hXk, nk); mxSetN(hXk, nk); mxSetPr(hXk, xwork); mxSetPi(hXk, xwork+nksqr); #ifdef USE_SVD mexCallMATLAB(3, output_array, 1, &hXk, "svd"); #else mexCallMATLAB(2, output_array, 1, &hXk, "eig"); #endif memcpy(q, mxGetPr(output_array[0]), nksqr * sizeof(double)); q += nksqr; if(mxIsComplex(output_array[0])) /* if any imaginary part */ memcpy(q, mxGetPi(output_array[0]), nksqr * sizeof(double)); nkp1 = nk + 1; /* copy diag(Lab) */ labk = mxGetPr(output_array[1]); for(i = 0, ii = 0; i < nk; i++, ii += nkp1) lab[i] = labk[ii]; mxDestroyArray(output_array[0]); mxDestroyArray(output_array[1]); #ifdef USE_SVD mxDestroyArray(output_array[2]); #endif lab += nk; x += 2 * nksqr; q += nksqr; } } /* [lab,q] = eigK */ /* ------------------------------------------------------------ Release PSD-working arrays. ------------------------------------------------------------ */ mxSetM(Xk,0); mxSetN(Xk,0); mxSetPr(Xk, (double *) NULL); mxDestroyArray(Xk); mxSetM(hXk,0); mxSetN(hXk,0); mxSetPr(hXk, (double *) NULL); mxSetPi(hXk, (double *) NULL); mxDestroyArray(hXk); mxFree(xwork); }
void mexFunction(int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[]) { double *xr = mxGetPr(prhs[0]); /* Input field vector, real part */ double *xi = mxGetPi(prhs[0]); /* Input field vector, imag part */ double *h1r = mxGetPr(prhs[1]); /* 1st output filter, real part */ double *h1i = mxGetPi(prhs[1]); /* 1st output filter, imag part */ double *h2r = mxGetPr(prhs[2]); /* 2nd output filter, real part */ double *h2i = mxGetPi(prhs[2]); /* 2nd output filter, imag part */ double Ntap = mxGetScalar(prhs[3]); /* Length of filters, real scalar*/ double mu = mxGetScalar(prhs[4]); /* Algo constant, real scalar */ double *R = mxGetPr(prhs[5]); /* Convergence radii, real vector*/ double sps = mxGetScalar(prhs[6]); /* Samples x symbol, 1 or 2 */ int Mdim = mxGetM(prhs[0]); /* Length of input field vector */ int Npol = mxGetN(prhs[0]); /* Should be 2 */ int Mfilter1 = mxGetM(prhs[1]); /* Should be == Ntap */ int Nfilter1 = mxGetN(prhs[1]); /* Should be 2 */ int Mfilter2 = mxGetM(prhs[2]); /* Should be == Ntap */ int Nfilter2 = mxGetN(prhs[2]); /* Should be 2 */ bool dontskip; /* Will be better to add some checks and error: * Ntap == Mfilter1 == Mfilter2 == ODD INTEGER * Npol == Nfilter1 == Nfilter2 == 2 */ if ((int)Ntap % 2 == 0) mexErrMsgTxt("Ntaps should be an ODD INTEGER."); /* Chechking the value of samples per symbol */ if ((int)sps == 1) dontskip = 1; else { if ((int)sps == 2) { dontskip = 0; } else { mexErrMsgTxt("Samples x symbol should be either 1 or 2."); } } /* Chech if the vector that are supposed to be complex are really complex */ /* If X, H1 or H2 are purely real numbers (which is often the case for H1 */ /* and H2) then matlab does not allocates the memory for the imaginary */ /* part and the program will segfault if it tries to access it */ if (xi==NULL) { /* If necessary allocates memory for the imaginary part of X */ xi = mxCalloc(Mdim*Npol, sizeof(double)); mxSetPi( (mxArray *)prhs[0], xi); } if (h1i==NULL) { /* If necessary allocates memory for the imaginary part of H1 */ h1i = mxCalloc(Ntap*Npol,sizeof(double)); mxSetPi( (mxArray *)prhs[1], h1i); } if (h2i==NULL) { /* If necessary allocates memory for the imaginary part of H2 */ h2i = mxCalloc(Ntap*Npol,sizeof(double)); mxSetPi( (mxArray *)prhs[2], h2i); } /* Next line allocates memory for a (Mdim-Ntap+1) x Npol complex vector */ plhs[0] = mxCreateDoubleMatrix(Mdim-Ntap+1,Npol,mxCOMPLEX); double *yr = mxGetPr(plhs[0]); /* Output field vector, real part */ double *yi = mxGetPi(plhs[0]); /* Output field vector, imag part */ cmafilter(xr,xi,Mdim,h1r,h1i,h2r,h2i,Ntap,mu,R,yr,yi,dontskip); plhs[1] = mxCreateDoubleMatrix(1,1,mxREAL); double *y = mxGetPr(plhs[1]); *y = 0; plhs[2] = mxCreateDoubleMatrix(1,1,mxREAL); double *z = mxGetPr(plhs[2]); *z = 0; return; }
/* the input problem could be reduced to a conventional 2D matrix product.. so here we use the fast MATLAB version */ mxArray* MATLAB_mm(MxInfo zinfo, const MxInfo xinfo, const MxInfo yinfo, const MxInfo xrest, const MxInfo yrest, const MxInfo xmacc, const MxInfo ymacc){ mxArray *Xmx, *Ymx, *Zmx, *args[2]; /* call matlab to do the real work -- more reliable than dgemm */ /* first create a new matrix with the right size to use in the matlab call*/ /* create and empty array */ Xmx= mxCreateNumericMatrix(0,0,xinfo.dtype,(xinfo.ip==0)?mxREAL:mxCOMPLEX); Ymx= mxCreateNumericMatrix(0,0,yinfo.dtype,(yinfo.ip==0)?mxREAL:mxCOMPLEX); /* and populate it with data */ mxSetPr(Xmx,xinfo.rp); if ( xinfo.ip ) mxSetPi(Xmx,xinfo.ip); mxSetPr(Ymx,yinfo.rp); if ( yinfo.ip ) mxSetPi(Ymx,yinfo.ip); /* Set trap so errors return here so we can clean up correctly */ mexSetTrapFlag(1); /* now do the calls to matlab to get the result */ args[0] = Xmx; args[1]= Ymx; Zmx=0; /* no accumulated dims -- just outer product, or just inner product */ if ( xmacc.nd == 1 && xmacc.sz[0]==1 ) { /* set X as its vector version, implicitly transpose Y and then prod */ mxSetM(Xmx,MAX(sz(xrest,0),1)); mxSetN(Xmx,1); mxSetM(Ymx,1); mxSetN(Ymx,MAX(sz(yrest,1),1)); mexCallMATLAB(1, &Zmx, 2, args, "*"); } else if ( xmacc.stride[0] == 1 && ymacc.stride[0] == 1 ){ /* | x | == macc/rest * macc/rest*/ mxSetM(Xmx,xmacc.sz[0]); mxSetN(Xmx,xinfo.numel/xmacc.sz[0]); mxSetM(Ymx,ymacc.sz[0]); mxSetN(Ymx,yinfo.numel/ymacc.sz[0]); /* transpose X and call matlab */ if( yinfo.numel+zinfo.numel>xinfo.numel*.75 ){/*cheaper to transpose X*/ mxArray *XmxT; mexCallMATLAB(1, &XmxT, 1, &Xmx, ".\'");/* N.B. this copies X!!! */ args[0]=XmxT; mexCallMATLAB(1, &Zmx, 2, args, "*"); mxDestroyArray(XmxT); } else { /*cheaper to transpose y and z */ mxArray *YmxT; mexCallMATLAB(1, &YmxT, 1, &Ymx, ".\'");/* N.B. this copies Y!!! */ args[0]=YmxT; args[1]=Xmx; mexCallMATLAB(1, &Zmx, 2, args, "*"); mxDestroyArray(YmxT); if( ymacc.sz[0]!=yinfo.numel ) {/* transpose result if necessary */ mxArray *ZmxT; mexCallMATLAB(1, &ZmxT, 1, &Zmx, ".\'");/* N.B. this copies Z!!! */ mxDestroyArray(Zmx); Zmx=ZmxT; } } } else if ( xmacc.stride[0] == 1 && ymacc.stride[0] > 1 ){ /* | x _ == macc/rest * rest/macc */ mxSetM(Xmx,xmacc.sz[0]); mxSetN(Xmx,xinfo.numel/xmacc.sz[0]); mxSetM(Ymx,ymacc.stride[0]); mxSetN(Ymx,yinfo.numel/ymacc.stride[0]); if( yinfo.numel+xinfo.numel<zinfo.numel*.75 ){/* cheaper to transpose Z */ /* reverse order of multiply, call matlab and transpose Z */ mxArray *ZmxT; args[0]=Ymx; args[1]=Xmx; mexCallMATLAB(1, &Zmx, 2, args, "*"); mexCallMATLAB(1, &ZmxT, 1, &Zmx, ".\'");/* N.B. this copies Z!!! */ mxDestroyArray(Zmx); Zmx=ZmxT; } else { /* cheaper to transpose twice */ mxArray *XmxT, *YmxT; mexCallMATLAB(1, &XmxT, 1, &Xmx, ".\'");/* N.B. this copies X!!! */ mexCallMATLAB(1, &YmxT, 1, &Ymx, ".\'");/* N.B. this copies Y!!! */ args[0]=XmxT; args[1]=YmxT; mexCallMATLAB(1, &Zmx, 2, args, "*"); mxDestroyArray(XmxT); mxDestroyArray(YmxT); } } else if ( xmacc.stride[0] > 1 && ymacc.stride[0] == 1 ){ /* _ x | == rest/macc * macc/rest */ mxSetM(Xmx,xmacc.stride[0]); mxSetN(Xmx,xinfo.numel/xmacc.stride[0]); mxSetM(Ymx,ymacc.sz[0]); mxSetN(Ymx,yinfo.numel/ymacc.sz[0]); mexCallMATLAB(1, &Zmx, 2, args, "*"); } else if ( xmacc.stride[0] > 1 && ymacc.stride[0] > 1 ){ /* _ x _ == rest/macc * rest/macc */ mxSetM(Xmx,xmacc.stride[0]); mxSetN(Xmx,xinfo.numel/xmacc.stride[0]); mxSetM(Ymx,ymacc.stride[0]); mxSetN(Ymx,yinfo.numel/ymacc.stride[0]); if( xinfo.numel+zinfo.numel>yinfo.numel*.75 ){/* cheaper to transpose Y */ /* transpose Y and call matlab */ mxArray *YmxT; mexCallMATLAB(1, &YmxT, 1, &Ymx, ".\'");/* N.B. this copies Y!!! */ args[1]=YmxT; mexCallMATLAB(1, &Zmx, 2, args, "*"); mxDestroyArray(YmxT); } else { /* cheaper to transpose X and Z */ mxArray *XmxT; mexCallMATLAB(1, &XmxT, 1, &Xmx, ".\'");/* N.B. this copies X!!! */ args[0]=Ymx; args[1]=XmxT; mexCallMATLAB(1, &Zmx, 2, args, "*"); mxDestroyArray(XmxT); if( xmacc.sz[0]!=xinfo.numel ) {/* transpose result if necessary */ mxArray *ZmxT; mexCallMATLAB(1, &ZmxT, 1, &Zmx, ".\'");/* N.B. this copies Z!!! */ mxDestroyArray(Zmx); Zmx=ZmxT; } } } else { ERROR("tprod: somethings gone horibbly wrong!"); } /* Set trap so errors return here so we can clean up correctly */ mexSetTrapFlag(0); /* set the tempory X and Y matrices back to empty without ref to data to stop matlab "helpfully" double freeing them? */ mxSetM(Xmx,0);mxSetN(Xmx,0);mxSetPr(Xmx,0);mxSetPi(Xmx,0); mxDestroyArray(Xmx); mxSetM(Ymx,0);mxSetN(Ymx,0);mxSetPr(Ymx,0);mxSetPi(Ymx,0); mxDestroyArray(Ymx); return Zmx; }
int flip(mxArray *X, const mxArray *DL, const mxArray *A, const mxArray *dDL, const mxArray *s, mwSize offset, mwSize T, const mxArray *h, const mxArray *wws, const mxArray *wVs) { const mwSize *sz = mxGetDimensions(dDL); mwSize D = sz[0], M = sz[1], p = sz[3], Ndt = sz[4]; mwSize N = mxGetM(DL); mwSize Tdt = ceil((double) N / (double) Ndt); mwSize lenh = mxGetM(h); mwSize pad = (lenh - 1) / 2; // find maximum double *DLpr = mxGetPr(DL); double max = 0, d; int iMax = -1, jMax = -1; for (mwSize j = 0; j != M; ++j) { for (mwSize i = offset + 1; i != offset + T + 1; ++i) { d = DLpr[N * j + i]; if (d >= max) { max = d; iMax = i; jMax = j; } } } double Xij; double *Xpr = mxGetPr(X), *Xpi = mxGetPi(X); double *hpr = mxGetPr(h), *Apr = mxGetPr(A); mwIndex *Xir = mxGetIr(X), *Xjc = mxGetJc(X); double sgn; if (max > 0) { // find location in sparse array double a = 0, r = 0; int sub; mwSize l = Xjc[jMax]; while (l != Xjc[jMax + 1] && Xir[l] <= iMax) { if (Xir[l] == iMax) { a = Xpr[l]; r = Xpi[l]; break; } else { ++l; } } if (a == 0) { // add spike - subsample // determine amplitude and subsample shift sgn = 1; max = 0; for (mwSize j = 0; j != p; ++j) { double m = 0; for (mwSize i = 0; i != lenh; ++i) { m = m + DLpr[jMax * N + iMax - pad + i] * hpr[j * lenh + i]; } if (m > max) { sub = j; max = m; } } for (mwSize i = 0; i != lenh; ++i) { a = a + Apr[N * jMax + iMax - pad + i] * hpr[sub * lenh + i]; } r = (double) (sub - (int) p / 2) / (double) p; // CHECK // grow sparse array if necessary mwSize nzmax = mxGetNzmax(X); if (Xjc[M] == nzmax) { // grow X nzmax *= 2; mxSetNzmax(X, nzmax); Xir = (mwIndex*) mxRealloc(Xir, nzmax * sizeof(*Xir)); mxSetIr(X, Xir); Xpr = (double*) mxRealloc(Xpr, nzmax * sizeof(*Xpr)); mxSetPr(X, Xpr); Xpi = (double*) mxRealloc(Xpi, nzmax * sizeof(*Xpi)); mxSetPi(X, Xpi); } // add values to sparse array // real: amplitude // imag: subsample (> 0 => shift right, < 0 => shift left) for (mwSize j = jMax; j != M; ++j) { ++Xjc[j + 1]; } for (mwSize i = Xjc[M] - 1; i > l; --i) { Xir[i] = Xir[i - 1]; Xpr[i] = Xpr[i - 1]; Xpi[i] = Xpi[i - 1]; } Xir[l] = iMax; Xpr[l] = a; Xpi[l] = r; } else { // remove spike sgn = -1; for (mwSize j = jMax + 1; j != M + 1; ++j) { --Xjc[j]; } for (; l != Xjc[M]; ++l) { Xir[l] = Xir[l + 1]; Xpr[l] = Xpr[l + 1]; Xpi[l] = Xpi[l + 1]; } } // update change in posterior double DLij = DLpr[N * jMax + iMax]; sub = (int) p / 2 - round(r * (double) p); mwSize t = iMax / Tdt; mwSize start = D * M * (jMax + M * (sub + p * t)); mwSize ii; double dA; double *dDLpr = mxGetPr(dDL), *wwspr = mxGetPr(wws), *wVspr = mxGetPr(wVs), *spr = mxGetPr(s); for (mwSize j = 0; j != M; ++j) { for (mwSize i = 0; i != D; ++i) { dA = dDLpr[start + D * j + i] * a * sgn / wwspr[Ndt * j + t]; ii = N * j + iMax + spr[i]; Apr[ii] = Apr[ii] - dA; DLpr[ii] = DLpr[ii] - dA * (wVspr[ii] + a * dDLpr[start + D * j + i]); DLpr[N * jMax + iMax] = -DLij; } } } else { iMax = -1; } return iMax; }
void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { Long *Ap, *Ai, *Zp, *Zi ; double *Ax, *Az, *Zx ; Long p, j, build_upper, zero_handling, nrow, ncol, mkind, skind, asize, znz, status ; char filename [LEN+1], title [73], key [9], mtype [4] ; /* ---------------------------------------------------------------------- */ /* check inputs */ /* ---------------------------------------------------------------------- */ if (nargin != 1 || nargout > 5 || !mxIsChar (pargin [0])) { mexErrMsgTxt ("Usage: [A Z title key mtype] = RBread (filename)") ; } /* ---------------------------------------------------------------------- */ /* get filename */ /* ---------------------------------------------------------------------- */ if (mxGetString (pargin [0], filename, LEN) != 0) { mexErrMsgTxt ("filename too long") ; } /* ---------------------------------------------------------------------- */ /* read the matrix */ /* ---------------------------------------------------------------------- */ build_upper = TRUE ; /* always build upper tri. part */ zero_handling = (nargout > 1) ? 2 : 1 ; /* prune or extract zeros */ status = RBread (filename, build_upper, zero_handling, title, key, mtype, &nrow, &ncol, &mkind, &skind, &asize, &znz, &Ap, &Ai, &Ax, &Az, &Zp, &Zi) ; if (status != RBIO_OK) { RBerror (status) ; mexErrMsgTxt ("error reading file") ; } /* ---------------------------------------------------------------------- */ /* return A to MATLAB */ /* ---------------------------------------------------------------------- */ pargout [0] = mxCreateSparse (0, 0, 0, (mkind == 2) ? mxCOMPLEX : mxREAL) ; mxFree (mxGetJc (pargout [0])) ; mxFree (mxGetIr (pargout [0])) ; mxFree (mxGetPr (pargout [0])) ; if (mkind == 2) mxFree (mxGetPi (pargout [0])) ; mxSetM (pargout [0], nrow) ; mxSetN (pargout [0], ncol) ; mxSetNzmax (pargout [0], asize) ; mxSetJc (pargout [0], (mwIndex *) Ap) ; mxSetIr (pargout [0], (mwIndex *) Ai) ; mxSetPr (pargout [0], Ax) ; if (mkind == 2) mxSetPi (pargout [0], Az) ; /* ---------------------------------------------------------------------- */ /* return Z to MATLAB */ /* ---------------------------------------------------------------------- */ if (nargout > 1) { Zx = (double *) SuiteSparse_malloc (znz, sizeof (double)) ; for (p = 0 ; p < znz ; p++) { Zx [p] = 1 ; } pargout [1] = mxCreateSparse (0, 0, 0, mxREAL) ; mxFree (mxGetJc (pargout [1])) ; mxFree (mxGetIr (pargout [1])) ; mxFree (mxGetPr (pargout [1])) ; mxSetM (pargout [1], nrow) ; mxSetN (pargout [1], ncol) ; mxSetNzmax (pargout [1], MAX (znz,1)) ; mxSetJc (pargout [1], (mwIndex *) Zp) ; mxSetIr (pargout [1], (mwIndex *) Zi) ; mxSetPr (pargout [1], Zx) ; } /* ---------------------------------------------------------------------- */ /* return title */ /* ---------------------------------------------------------------------- */ if (nargout > 2) { pargout [2] = mxCreateString (title) ; } /* ---------------------------------------------------------------------- */ /* return key */ /* ---------------------------------------------------------------------- */ if (nargout > 3) { pargout [3] = mxCreateString (key) ; } /* ---------------------------------------------------------------------- */ /* return mtype */ /* ---------------------------------------------------------------------- */ if (nargout > 4) { pargout [4] = mxCreateString (mtype) ; } }
mxArray *ssmult_dot /* returns C = A'*B */ ( const mxArray *A, const mxArray *B, int ac, /* if true: conj(A) if false: A. ignored if A real */ int bc, /* if true: conj(B) if false: B. ignored if B real */ int cc /* if true: conj(C) if false: C. ignored if C real */ ) { double cx, cz, ax, az, bx, bz ; mxArray *C ; double *Ax, *Az, *Bx, *Bz, *Cx, *Cz ; Int *Ap, *Ai, *Bp, *Bi, *Cp, *Ci ; Int m, n, k, cnzmax, i, j, p, paend, pbend, ai, bi, cnz, pa, pb, zallzero, A_is_complex, B_is_complex, C_is_complex ; /* ---------------------------------------------------------------------- */ /* get inputs */ /* ---------------------------------------------------------------------- */ m = mxGetM (A) ; n = mxGetN (A) ; k = mxGetN (B) ; if (m != mxGetM (B)) ssmult_invalid (ERROR_DIMENSIONS) ; Ap = mxGetJc (A) ; Ai = mxGetIr (A) ; Ax = mxGetPr (A) ; Az = mxGetPi (A) ; A_is_complex = mxIsComplex (A) ; Bp = mxGetJc (B) ; Bi = mxGetIr (B) ; Bx = mxGetPr (B) ; Bz = mxGetPi (B) ; B_is_complex = mxIsComplex (B) ; /* ---------------------------------------------------------------------- */ /* allocate C as an n-by-k full matrix but do not initialize it */ /* ---------------------------------------------------------------------- */ /* NOTE: integer overflow cannot occur here, because this function is not called unless O(n*k) is less than O(m+nnz(A)). The test is done in the caller, not here. */ cnzmax = n*k ; cnzmax = MAX (cnzmax, 1) ; Cx = mxMalloc (cnzmax * sizeof (double)) ; C_is_complex = A_is_complex || B_is_complex ; Cz = C_is_complex ? mxMalloc (cnzmax * sizeof (double)) : NULL ; /* ---------------------------------------------------------------------- */ /* C = A'*B using sparse dot products */ /* ---------------------------------------------------------------------- */ /* NOTE: this method REQUIRES the columns of A and B to be sorted on input. That is, the row indices in each column must appear in ascending order. This is the standard in all versions of MATLAB to date, and likely will be for some time. However, if MATLAB were to use unsorted sparse matrices in the future (a lazy sort) then a test should be included in ssmult to not use ssmult_dot if A or B are unsorted, or they should be sorted on input. */ cnz = 0 ; for (j = 0 ; j < k ; j++) { for (i = 0 ; i < n ; i++) { /* compute C (i,j) = A (:,i)' * B (:,j) */ pa = Ap [i] ; paend = Ap [i+1] ; pb = Bp [j] ; pbend = Bp [j+1] ; if (pa == paend /* nnz (A (:,i)) == 0 */ || pb == pbend /* nnz (B (:,j)) == 0 */ || Ai [paend-1] < Bi [pb] /* max(find(A(:,i)))<min(find(B(:,j))) */ || Ai [pa] > Bi [pbend-1]) /* min(find(A(:,i)))>max(find(B(:,j))) */ { Cx [i+j*n] = 0 ; /* no work to do */ if (C_is_complex) { Cz [i+j*n] = 0 ; } continue ; } cx = 0 ; cz = 0 ; while (pa < paend && pb < pbend) { /* The dot product looks like the merge in mergesort, except */ /* no "clean-up" phase is need when one list is exhausted. */ ai = Ai [pa] ; bi = Bi [pb] ; if (ai == bi) { /* c += A (ai,i) * B (ai,j), and "consume" both entries */ if (!C_is_complex) { cx += Ax [pa] * Bx [pb] ; } else { /* complex case */ ax = Ax [pa] ; bx = Bx [pb] ; az = Az ? (ac ? (-Az [pa]) : Az [pa]) : 0.0 ; bz = Bz ? (bc ? (-Bz [pb]) : Bz [pb]) : 0.0 ; cx += ax * bx - az * bz ; cz += az * bx + ax * bz ; } pa++ ; pb++ ; } else if (ai < bi) { /* consume A(ai,i) and discard it, since B(ai,j) is zero */ pa++ ; } else { /* consume B(bi,j) and discard it, since A(ai,i) is zero */ pb++ ; } } Cx [i+j*n] = cx ; if (C_is_complex) { Cz [i+j*n] = cz ; } } /* count the number of nonzeros in C(:,j) */ for (i = 0 ; i < n ; i++) { /* This could be done above, except for the gcc compiler bug when cx is an 80-bit nonzero in register above, but becomes zero here when stored into memory. We need the latter, to correctly handle the case when cx underflows to zero in 64-bit floating-point. Do not attempt to "optimize" this code by doing this test above, unless the gcc compiler bug is fixed (as of gcc version 4.1.0). */ if (Cx [i+j*n] != 0 || (C_is_complex && Cz [i+j*n] != 0)) { cnz++ ; } } } /* ---------------------------------------------------------------------- */ /* convert C to real if the imaginary part is all zero */ /* ---------------------------------------------------------------------- */ if (C_is_complex) { zallzero = 1 ; for (p = 0 ; zallzero && p < cnzmax ; p++) { if (Cz [p] != 0) { zallzero = 0 ; } } if (zallzero) { /* the imaginary part of C is all zero */ C_is_complex = 0 ; mxFree (Cz) ; Cz = NULL ; } } /* ---------------------------------------------------------------------- */ /* allocate integer part of C but do not initialize it */ /* ---------------------------------------------------------------------- */ cnz = MAX (cnz, 1) ; C = mxCreateSparse (0, 0, 0, C_is_complex ? mxCOMPLEX : mxREAL) ; mxFree (mxGetJc (C)) ; mxFree (mxGetIr (C)) ; mxFree (mxGetPr (C)) ; mxFree (mxGetPi (C)) ; Cp = mxMalloc ((k + 1) * sizeof (Int)) ; Ci = mxMalloc (cnz * sizeof (Int)) ; mxSetJc (C, Cp) ; mxSetIr (C, Ci) ; mxSetM (C, n) ; mxSetN (C, k) ; /* ---------------------------------------------------------------------- */ /* C = sparse (C). Note that this is done in-place. */ /* ---------------------------------------------------------------------- */ p = 0 ; for (j = 0 ; j < k ; j++) { Cp [j] = p ; for (i = 0 ; i < n ; i++) { cx = Cx [i+j*n] ; cz = (C_is_complex ? Cz [i+j*n] : 0) ; if (cx != 0 || cz != 0) { Ci [p] = i ; Cx [p] = cx ; if (C_is_complex) Cz [p] = (cc ? (-cz) : cz) ; p++ ; } } } Cp [k] = p ; /* ---------------------------------------------------------------------- */ /* reduce the size of Cx and Cz and return result */ /* ---------------------------------------------------------------------- */ if (cnz < cnzmax) { Cx = mxRealloc (Cx, cnz * sizeof (double)) ; if (C_is_complex) { Cz = mxRealloc (Cz, cnz * sizeof (double)) ; } } mxSetNzmax (C, cnz) ; mxSetPr (C, Cx) ; if (C_is_complex) { mxSetPi (C, Cz) ; } return (C) ; }
void mexFunction( int nlhs, mxArray *plhs[], int nrhs, const mxArray *prhs[] ) { /* Declare variable */ mwSize m,n; mwSize nzmax; mwIndex *irs,*jcs,j,k; int cmplx,isfull; double *pr,*pi,*si,*sr; double percent_sparse; /* Check for proper number of input and output arguments */ if (nrhs != 1) { mexErrMsgIdAndTxt( "MATLAB:fulltosparse:invalidNumInputs", "One input argument required."); } if(nlhs > 1){ mexErrMsgIdAndTxt( "MATLAB:fulltosparse:maxlhs", "Too many output arguments."); } /* Check data type of input argument */ if (!(mxIsDouble(prhs[0]))){ mexErrMsgIdAndTxt( "MATLAB:fulltosparse:inputNotDouble", "Input argument must be of type double."); } if (mxGetNumberOfDimensions(prhs[0]) != 2){ mexErrMsgIdAndTxt( "MATLAB:fulltosparse:inputNot2D", "Input argument must be two dimensional\n"); } /* Get the size and pointers to input data */ m = mxGetM(prhs[0]); n = mxGetN(prhs[0]); pr = mxGetPr(prhs[0]); pi = mxGetPi(prhs[0]); cmplx = (pi==NULL ? 0 : 1); /* Allocate space for sparse matrix * NOTE: Assume at most 20% of the data is sparse. Use ceil * to cause it to round up. */ percent_sparse = 0.2; nzmax=(mwSize)ceil((double)m*(double)n*percent_sparse); plhs[0] = mxCreateSparse(m,n,nzmax,cmplx); sr = mxGetPr(plhs[0]); si = mxGetPi(plhs[0]); irs = mxGetIr(plhs[0]); jcs = mxGetJc(plhs[0]); /* Copy nonzeros */ k = 0; isfull=0; for (j=0; (j<n); j++) { mwSize i; jcs[j] = k; for (i=0; (i<m ); i++) { if (IsNonZero(pr[i]) || (cmplx && IsNonZero(pi[i]))) { /* Check to see if non-zero element will fit in * allocated output array. If not, increase percent_sparse * by 10%, recalculate nzmax, and augment the sparse array */ if (k>=nzmax){ mwSize oldnzmax = nzmax; percent_sparse += 0.1; nzmax = (mwSize)ceil((double)m*(double)n*percent_sparse); /* make sure nzmax increases atleast by 1 */ if (oldnzmax == nzmax) nzmax++; mxSetNzmax(plhs[0], nzmax); mxSetPr(plhs[0], mxRealloc(sr, nzmax*sizeof(double))); if(si != NULL) mxSetPi(plhs[0], mxRealloc(si, nzmax*sizeof(double))); mxSetIr(plhs[0], mxRealloc(irs, nzmax*sizeof(mwIndex))); sr = mxGetPr(plhs[0]); si = mxGetPi(plhs[0]); irs = mxGetIr(plhs[0]); } sr[k] = pr[i]; if (cmplx){ si[k]=pi[i]; } irs[k] = i; k++; } } pr += m; pi += m; } jcs[n] = k; }