int cholmod_l_drop (double tol, cholmod_sparse * A) { double aij; double *Ax; long long *Ap, *Ai, *Anz; long long packed, i, j, nrow, ncol, p, pend, nz, values; Ap = A->p; ncol = A->ncol; nz = 0; for (j = 0; j < ncol; j++) for (; p < pend; p++) { i = Ai[p]; aij = Ax[p]; if (i <= j && (fabs (aij) > tol || ((aij) != (aij)))) nz++; } Ap[ncol] = nz; cholmod_l_reallocate_sparse (nz, A, 0); }
void mexFunction ( int nargout, mxArray *pargout [ ], int nargin, const mxArray *pargin [ ] ) { void *G ; cholmod_dense *X = NULL ; cholmod_sparse *A = NULL, *Z = NULL ; cholmod_common Common, *cm ; Long *Ap = NULL, *Ai ; double *Ax, *Az = NULL ; char filename [MAXLEN] ; Long nz, k, is_complex = FALSE, nrow = 0, ncol = 0, allzero ; int mtype ; /* ---------------------------------------------------------------------- */ /* start CHOLMOD and set parameters */ /* ---------------------------------------------------------------------- */ cm = &Common ; cholmod_l_start (cm) ; sputil_config (SPUMONI, cm) ; /* ---------------------------------------------------------------------- */ /* get inputs */ /* ---------------------------------------------------------------------- */ if (nargin < 1 || nargin > 2 || nargout > 2) { mexErrMsgTxt ("usage: [A Z] = mread (filename, prefer_binary)") ; } if (!mxIsChar (pargin [0])) { mexErrMsgTxt ("mread requires a filename") ; } mxGetString (pargin [0], filename, MAXLEN) ; sputil_file = fopen (filename, "r") ; if (sputil_file == NULL) { mexErrMsgTxt ("cannot open file") ; } if (nargin > 1) { cm->prefer_binary = (mxGetScalar (pargin [1]) != 0) ; } /* ---------------------------------------------------------------------- */ /* read the matrix, as either a dense or sparse matrix */ /* ---------------------------------------------------------------------- */ G = cholmod_l_read_matrix (sputil_file, 1, &mtype, cm) ; fclose (sputil_file) ; sputil_file = NULL ; if (G == NULL) { mexErrMsgTxt ("could not read file") ; } /* get the specific matrix (A or X), and change to ZOMPLEX if needed */ if (mtype == CHOLMOD_SPARSE) { A = (cholmod_sparse *) G ; nrow = A->nrow ; ncol = A->ncol ; is_complex = (A->xtype == CHOLMOD_COMPLEX) ; Ap = A->p ; Ai = A->i ; if (is_complex) { /* if complex, ensure A is ZOMPLEX */ cholmod_l_sparse_xtype (CHOLMOD_ZOMPLEX, A, cm) ; } Ax = A->x ; Az = A->z ; } else if (mtype == CHOLMOD_DENSE) { X = (cholmod_dense *) G ; nrow = X->nrow ; ncol = X->ncol ; is_complex = (X->xtype == CHOLMOD_COMPLEX) ; if (is_complex) { /* if complex, ensure X is ZOMPLEX */ cholmod_l_dense_xtype (CHOLMOD_ZOMPLEX, X, cm) ; } Ax = X->x ; Az = X->z ; } else { mexErrMsgTxt ("invalid file") ; } /* ---------------------------------------------------------------------- */ /* if requested, extract the zero entries and place them in Z */ /* ---------------------------------------------------------------------- */ if (nargout > 1) { if (mtype == CHOLMOD_SPARSE) { /* A is a sparse real/zomplex double matrix */ Z = sputil_extract_zeros (A, cm) ; } else { /* input is full; just return an empty Z matrix */ Z = cholmod_l_spzeros (nrow, ncol, 0, CHOLMOD_REAL, cm) ; } } /* ---------------------------------------------------------------------- */ /* prune the zero entries from A and set nzmax(A) to nnz(A) */ /* ---------------------------------------------------------------------- */ if (mtype == CHOLMOD_SPARSE) { sputil_drop_zeros (A) ; cholmod_l_reallocate_sparse (cholmod_l_nnz (A, cm), A, cm) ; } /* ---------------------------------------------------------------------- */ /* change a complex matrix to real if its imaginary part is all zero */ /* ---------------------------------------------------------------------- */ if (is_complex) { if (mtype == CHOLMOD_SPARSE) { nz = Ap [ncol] ; } else { nz = nrow * ncol ; } allzero = TRUE ; for (k = 0 ; k < nz ; k++) { if (Az [k] != 0) { allzero = FALSE ; break ; } } if (allzero) { /* discard the all-zero imaginary part */ if (mtype == CHOLMOD_SPARSE) { cholmod_l_sparse_xtype (CHOLMOD_REAL, A, cm) ; } else { cholmod_l_dense_xtype (CHOLMOD_REAL, X, cm) ; } } } /* ---------------------------------------------------------------------- */ /* return results to MATLAB */ /* ---------------------------------------------------------------------- */ if (mtype == CHOLMOD_SPARSE) { pargout [0] = sputil_put_sparse (&A, cm) ; } else { pargout [0] = sputil_put_dense (&X, cm) ; } if (nargout > 1) { pargout [1] = sputil_put_sparse (&Z, cm) ; } /* ---------------------------------------------------------------------- */ /* free workspace */ /* ---------------------------------------------------------------------- */ cholmod_l_finish (cm) ; cholmod_l_print_common (" ", cm) ; }
template <typename Entry> SuiteSparseQR_factorization <Entry> *spqr_1factor ( // inputs, not modified int ordering, // all, except 3:given treated as 0:fixed double tol, // only accept singletons above tol. If tol <= -2, // then use the default tolerance Int bncols, // number of columns of B int keepH, // if TRUE, keep the Householder vectors cholmod_sparse *A, // m-by-n sparse matrix Int ldb, // if dense, the leading dimension of B Int *Bp, // size bncols+1, column pointers of B Int *Bi, // size bnz = Bp [bncols], row indices of B Entry *Bx, // size bnz, numerical values of B // workspace and parameters cholmod_common *cc ) { spqr_symbolic *QRsym ; spqr_numeric <Entry> *QRnum ; SuiteSparseQR_factorization <Entry> *QR ; Int *Yp, *Yi, *Q1fill, *R1p, *R1j, *P1inv, *Ap, *Ai ; Entry *Yx, *R1x, *Ax ; Int noY, anz, a2nz, r1nz, ynz, i, j, k, p, p2, bnz, py, n1rows, n1cols, n2, Bsparse, d, iold, inew, m, n ; cholmod_sparse *Y ; #ifdef TIMING double t0 = spqr_time ( ) ; double t1, t2 ; #endif // ------------------------------------------------------------------------- // get inputs and allocate result // ------------------------------------------------------------------------- m = A->nrow ; n = A->ncol ; Ap = (Int *) A->p ; Ai = (Int *) A->i ; Ax = (Entry *) A->x ; QR = (SuiteSparseQR_factorization <Entry> *) cholmod_l_malloc (1, sizeof (SuiteSparseQR_factorization <Entry>), cc) ; if (cc->status < CHOLMOD_OK) { // out of memory return (NULL) ; } QR->QRsym = NULL ; QR->QRnum = NULL ; QR->R1p = NULL ; QR->R1j = NULL ; QR->R1x = NULL ; QR->P1inv = NULL ; QR->Q1fill = NULL ; QR->Rmap = NULL ; QR->RmapInv = NULL ; QR->HP1inv = NULL ; QR->narows = m ; QR->nacols = n ; QR->bncols = bncols ; QR->n1rows = 0 ; QR->n1cols = 0 ; QR->r1nz = 0 ; r1nz = 0 ; // B is an optional input. It can be sparse or dense Bsparse = (Bp != NULL && Bi != NULL) ; if (Bx == NULL) { // B is not present; force bncols to be zero bncols = 0 ; } // ------------------------------------------------------------------------- // find the default tol, if requested // ------------------------------------------------------------------------- if (tol <= SPQR_DEFAULT_TOL) { tol = spqr_tol <Entry> (A, cc) ; } if (tol < 0) { // no rank detection will be performed QR->allow_tol = FALSE ; tol = EMPTY ; } else { QR->allow_tol = TRUE ; } QR->tol = tol ; // ------------------------------------------------------------------------- // find singletons and construct column pointers for the A part of Y // ------------------------------------------------------------------------- // These return R1p, P1inv, and Y; but they are all NULL if out of memory. // Note that only Y->p is allocated (Y->i and Y->x are dummy placeholders // of one Int and one Entry, each, actually). The entries of Y are // allocated later, below. if (ordering == SPQR_ORDERING_GIVEN) { ordering = SPQR_ORDERING_FIXED ; } if (ordering == SPQR_ORDERING_FIXED) { // fixed ordering: find column singletons without permuting columns Q1fill = NULL ; spqr_1fixed <Entry> (tol, bncols, A, &R1p, &P1inv, &Y, &n1cols, &n1rows, cc) ; } else { // natural or fill-reducing ordering: find column singletons with // column permutations allowed, then permute the pruned submatrix with // a fill-reducing ordering if ordering is not SPQR_ORDERING_NATURAL. spqr_1colamd <Entry> (ordering, tol, bncols, A, &Q1fill, &R1p, &P1inv, &Y, &n1cols, &n1rows, cc) ; ordering = cc->SPQR_istat [7] ; } if (cc->status < CHOLMOD_OK) { // out of memory spqr_freefac (&QR, cc) ; return (NULL) ; } QR->R1p = R1p ; QR->P1inv = P1inv ; QR->Q1fill = Q1fill ; QR->n1rows = n1rows ; QR->n1cols = n1cols ; noY = (Y == NULL) ; // A will be factorized, not Y ASSERT (noY == (n1cols == 0 && bncols == 0)) ; Yp = noY ? NULL : (Int *) Y->p ; anz = Ap [n] ; // nonzeros in A a2nz = noY ? anz : Yp [n-n1cols] ; // nonzeros in S2 n2 = n - n1cols ; // number of columns of S2 // Y is NULL, or of size (m-n1rows)-by-(n-n1cols+bncols) ASSERT (IMPLIES (Y != NULL, ((Int) Y->nrow == m-n1rows))) ; ASSERT (IMPLIES (Y != NULL, ((Int) Y->ncol == n-n1cols+bncols))) ; // Y, if allocated, has no space for any entries yet ynz = 0 ; // ------------------------------------------------------------------------- // construct the column pointers for the B or B2 part of Y // ------------------------------------------------------------------------- if (noY) { // A will be factorized instead of Y. There is no B. C or X can exist // as empty matrices with rows but no columns ASSERT (Yp == NULL) ; ASSERT (R1p == NULL) ; ASSERT (P1inv == NULL) ; ASSERT (n1rows == 0) ; ASSERT (a2nz == Ap [n]) ; ASSERT (bncols == 0) ; } else if (n1cols == 0) { // --------------------------------------------------------------------- // construct the column pointers for the B part of Y = [S B] // --------------------------------------------------------------------- ASSERT (R1p == NULL) ; ASSERT (P1inv == NULL) ; ASSERT (n1rows == 0) ; ASSERT (a2nz == Ap [n]) ; ynz = a2nz ; if (Bsparse) { // B is sparse for (k = 0 ; k < bncols ; k++) { Yp [(n-n1cols)+k] = ynz ; d = Bp [k+1] - Bp [k] ; ynz += d ; } } else { // B is dense Entry *B1 = Bx ; for (k = 0 ; k < bncols ; k++) { // count the nonzero entries in column k of B Yp [(n-n1cols)+k] = ynz ; d = 0 ; for (i = 0 ; i < m ; i++) { if (B1 [i] != (Entry) 0) { d++ ; } } B1 += ldb ; ynz += d ; } } Yp [(n-n1cols)+bncols] = ynz ; } else { // --------------------------------------------------------------------- // construct the column pointers for the B2 part of Y = [S2 B2] // --------------------------------------------------------------------- ynz = a2nz ; if (Bsparse) { // B is sparse for (k = 0 ; k < bncols ; k++) { // count the nonzero entries in column k of B2 Yp [(n-n1cols)+k] = ynz ; d = 0 ; for (p = Bp [k] ; p < Bp [k+1] ; p++) { iold = Bi [p] ; inew = P1inv [iold] ; if (inew >= n1rows) { d++ ; } } ynz += d ; } } else { // B is dense Entry *B1 = Bx ; for (k = 0 ; k < bncols ; k++) { // count the nonzero entries in column k of B2 Yp [(n-n1cols)+k] = ynz ; d = 0 ; for (iold = 0 ; iold < m ; iold++) { inew = P1inv [iold] ; if (inew >= n1rows && B1 [iold] != (Entry) 0) { d++ ; } } B1 += ldb ; ynz += d ; } } Yp [(n-n1cols)+bncols] = ynz ; } // ------------------------------------------------------------------------- // allocate the nonzeros for Y // ------------------------------------------------------------------------- if (noY) { // no singletons found, and B is empty. pass Y=A to QR factorization, // and pass in Q1fill as the "user-provided" ordering ASSERT (Yp == NULL) ; Yi = NULL ; Yx = NULL ; } else { cholmod_l_reallocate_sparse (ynz, Y, cc) ; Yi = (Int *) Y->i ; Yx = (Entry *) Y->x ; } if (cc->status < CHOLMOD_OK) { // out of memory spqr_freefac (&QR, cc) ; cholmod_l_free_sparse (&Y, cc) ; return (NULL) ; } // ------------------------------------------------------------------------- // create the pattern and values of Y and R1 // ------------------------------------------------------------------------- if (noY) { // --------------------------------------------------------------------- // R1 does not exist // --------------------------------------------------------------------- ASSERT (R1p == NULL) ; R1j = NULL ; R1x = NULL ; } else if (n1cols == 0) { // --------------------------------------------------------------------- // R1 does not exist // --------------------------------------------------------------------- ASSERT (R1p == NULL) ; R1j = NULL ; R1x = NULL ; // --------------------------------------------------------------------- // construct the A part of Y = [S B] // --------------------------------------------------------------------- ASSERT (anz == a2nz) ; py = 0 ; for (k = 0 ; k < n ; k++) { j = Q1fill ? Q1fill [k] : k ; ASSERT (py == Yp [k]) ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { Yi [py] = Ai [p] ; Yx [py] = Ax [p] ; py++ ; } } ASSERT (py == anz) ; ASSERT (py == Yp [n]) ; // --------------------------------------------------------------------- // construct the B part of Y = [S B] // --------------------------------------------------------------------- if (Bsparse) { // B is sparse bnz = Bp [bncols] ; for (p = 0 ; p < bnz ; p++) { Yi [py++] = Bi [p] ; } py = anz ; for (p = 0 ; p < bnz ; p++) { Yx [py++] = Bx [p] ; } } else { // B is dense Entry *B1 = Bx ; for (k = 0 ; k < bncols ; k++) { ASSERT (py == Yp [n+k]) ; for (i = 0 ; i < m ; i++) { Entry bij = B1 [i] ; if (bij != (Entry) 0) { Yi [py] = i ; Yx [py] = bij ; py++ ; } } B1 += ldb ; } } ASSERT (py == ynz) ; } else { // --------------------------------------------------------------------- // R1p = cumsum ([0 R1p]) // --------------------------------------------------------------------- r1nz = spqr_cumsum (n1rows, R1p) ; // Int overflow cannot occur PR (("total nonzeros in R1: %ld\n", r1nz)) ; // --------------------------------------------------------------------- // allocate R1 // --------------------------------------------------------------------- R1j = (Int *) cholmod_l_malloc (r1nz, sizeof (Int ), cc) ; R1x = (Entry *) cholmod_l_malloc (r1nz, sizeof (Entry), cc) ; QR->R1j = R1j ; QR->R1x = R1x ; QR->r1nz = r1nz ; if (cc->status < CHOLMOD_OK) { // out of memory spqr_freefac (&QR, cc) ; cholmod_l_free_sparse (&Y, cc) ; return (NULL) ; } // --------------------------------------------------------------------- // scan A and construct R11 // --------------------------------------------------------------------- // At this point, R1p [i] points to the start of row i: // for (Int t = 0 ; t <= n1rows ; t++) Rsave [t] = R1p [t] ; for (k = 0 ; k < n1cols ; k++) { j = Q1fill ? Q1fill [k] : k ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { // row i of A is row inew after singleton permutation i = Ai [p] ; inew = P1inv [i] ; ASSERT (inew < n1rows) ; // A (i,j) is in a singleton row. It becomes R1 (inew,k) p2 = R1p [inew]++ ; ASSERT (p2 < R1p [inew+1]) ; R1j [p2] = k ; R1x [p2] = Ax [p] ; } } // --------------------------------------------------------------------- // scan A and construct R12 and the S2 part of Y = [S2 B2] // --------------------------------------------------------------------- py = 0 ; for ( ; k < n ; k++) { j = Q1fill ? Q1fill [k] : k ; ASSERT (py == Yp [k-n1cols]) ; for (p = Ap [j] ; p < Ap [j+1] ; p++) { // row i of A is row inew after singleton permutation i = Ai [p] ; inew = P1inv [i] ; if (inew < n1rows) { // A (i,j) is in a singleton row. It becomes R1 (inew,k) p2 = R1p [inew]++ ; ASSERT (p2 < R1p [inew+1]) ; R1j [p2] = k ; R1x [p2] = Ax [p] ; } else { // A (i,j) is not in a singleton row. Place it in // Y (inew-n1rows, k-n1cols) Yi [py] = inew - n1rows ; Yx [py] = Ax [p] ; py++ ; } } } ASSERT (py == Yp [n-n1cols]) ; // --------------------------------------------------------------------- // restore the row pointers for R1 // --------------------------------------------------------------------- spqr_shift (n1rows, R1p) ; // the row pointers are back to what they were: // for (Int t = 0 ; t <= n1rows ; t++) ASSERT (Rsave [t] == R1p [t]) ; // --------------------------------------------------------------------- // construct the B2 part of Y = [S2 B2] // --------------------------------------------------------------------- if (Bsparse) { // B is sparse for (k = 0 ; k < bncols ; k++) { // construct the nonzero entries in column k of B2 ASSERT (py == Yp [k+(n-n1cols)]) ; for (p = Bp [k] ; p < Bp [k+1] ; p++) { iold = Bi [p] ; inew = P1inv [iold] ; if (inew >= n1rows) { Yi [py] = inew - n1rows ; Yx [py] = Bx [p] ; py++ ; } } } } else { // B is dense Entry *B1 = Bx ; for (k = 0 ; k < bncols ; k++) { // construct the nonzero entries in column k of B2 ASSERT (py == Yp [k+(n-n1cols)]) ; for (iold = 0 ; iold < m ; iold++) { inew = P1inv [iold] ; if (inew >= n1rows) { Entry bij = B1 [iold] ; if (bij != (Entry) 0) { Yi [py] = inew - n1rows ; Yx [py] = bij ; py++ ; } } } B1 += ldb ; } } ASSERT (py == ynz) ; } // ------------------------------------------------------------------------- // QR factorization of A or Y // ------------------------------------------------------------------------- if (noY) { // factorize A, with fill-reducing ordering already given in Q1fill QRsym = spqr_analyze (A, SPQR_ORDERING_GIVEN, Q1fill, tol >= 0, keepH, cc) ; #ifdef TIMING t1 = spqr_time ( ) ; #endif QRnum = spqr_factorize <Entry> (&A, FALSE, tol, n, QRsym, cc) ; } else { // fill-reducing ordering is already applied to Y; free Y when loaded QRsym = spqr_analyze (Y, SPQR_ORDERING_FIXED, NULL, tol >= 0, keepH, cc) ; #ifdef TIMING t1 = spqr_time ( ) ; #endif QRnum = spqr_factorize <Entry> (&Y, TRUE, tol, n2, QRsym, cc) ; // Y has been freed ASSERT (Y == NULL) ; } // record the actual ordering used (this will have been changed to GIVEN // or FIXED, in spqr_analyze, but change it back to the ordering used by // spqr_1fixed or spqr_1colamd. cc->SPQR_istat [7] = ordering ; QR->QRsym = QRsym ; QR->QRnum = QRnum ; if (cc->status < CHOLMOD_OK) { // out of memory spqr_freefac (&QR, cc) ; return (NULL) ; } cc->SPQR_istat [0] += r1nz ; // nnz (R) // rank estimate of A, including singletons but excluding the columns of // of B, in case [A B] was factorized. QR->rank = n1rows + QRnum->rank1 ; // ------------------------------------------------------------------------- // construct global row permutation if H is kept and singletons exist // ------------------------------------------------------------------------- // If there are no singletons, then HP1inv [0:m-1] and HPinv [0:m-1] would // be identical, so HP1inv is not needed. ASSERT ((n1cols == 0) == (P1inv == NULL)) ; ASSERT (IMPLIES (n1cols == 0, n1rows == 0)) ; if (keepH && n1cols > 0) { // construct the global row permutation. Currently, the row indices // in H reflect the global R. P1inv is the singleton permutation, // where a row index of Y = (P1inv (row of A) - n1rows), and // row of R2 = QRnum->HPinv (row of Y). Combine these two into // HP1inv, where a global row of R = HP1inv (a row of A) Int kk ; Int *HP1inv, *HPinv ; QR->HP1inv = HP1inv = (Int *) cholmod_l_malloc (m, sizeof (Int), cc) ; HPinv = QRnum->HPinv ; if (cc->status < CHOLMOD_OK) { // out of memory spqr_freefac (&QR, cc) ; return (NULL) ; } for (i = 0 ; i < m ; i++) { // i is a row of A, k is a row index after row singletons are // permuted. Then kk is a row index of the global R. k = P1inv ? P1inv [i] : i ; ASSERT (k >= 0 && k < m) ; if (k < n1rows) { kk = k ; } else { // k-n1rows is a row index of Y, the matrix factorized by // the QR factorization kernels (in QRsym and QRnum). // HPinv [k-n1rows] gives a row index of R2, to which n1rows // must be added to give a row of the global R. kk = HPinv [k - n1rows] + n1rows ; } ASSERT (kk >= 0 && kk < m) ; HP1inv [i] = kk ; } } // ------------------------------------------------------------------------- // find the mapping for the squeezed R, if A is rank deficient // ------------------------------------------------------------------------- if (QR->rank < n && !spqr_rmap <Entry> (QR, cc)) { // out of memory spqr_freefac (&QR, cc) ; return (NULL) ; } // ------------------------------------------------------------------------- // output statistics // ------------------------------------------------------------------------- cc->SPQR_istat [4] = QR->rank ; // estimated rank of A cc->SPQR_istat [5] = n1cols ; // number of columns singletons cc->SPQR_istat [6] = n1rows ; // number of singleton rows cc->SPQR_xstat [1] = tol ; // tol used #ifdef TIMING t2 = spqr_time ( ) ; cc->other1 [1] = t1 - t0 ; // analyze time, including singletons cc->other1 [2] = t2 - t1 ; // factorize time #endif return (QR) ; }