template <typename Entry> Long spqr_front ( // input, not modified Long m, // F is m-by-n with leading dimension m Long n, Long npiv, // number of pivot columns double tol, // a column is flagged as dead if its norm is <= tol Long ntol, // apply tol only to first ntol pivot columns Long fchunk, // block size for compact WY Householder reflections, // treated as 1 if fchunk <= 1 // input/output Entry *F, // frontal matrix F of size m-by-n Long *Stair, // size n, entries F (Stair[k]:m-1, k) are all zero, // for each k = 0:n-1, and remain zero on output. char *Rdead, // size npiv; all zero on input. If k is dead, // Rdead [k] is set to 1 // output, not defined on input Entry *Tau, // size n, Householder coefficients // workspace, undefined on input and output Entry *W, // size b*n, where b = min (fchunk,n,m) // input/output double *wscale, double *wssq, cholmod_common *cc // for cc->hypotenuse function ) { Entry tau ; double wk ; Entry *V ; Long k, t, g, g1, nv, k1, k2, i, t0, vzeros, mleft, nleft, vsize, minchunk, rank ; // NOTE: inputs are not checked for NULL (except if debugging enabled) ASSERT (F != NULL) ; ASSERT (Stair != NULL) ; ASSERT (Rdead != NULL) ; ASSERT (Tau != NULL) ; ASSERT (W != NULL) ; ASSERT (m >= 0 && n >= 0) ; npiv = MAX (0, npiv) ; // npiv must be between 0 and n npiv = MIN (n, npiv) ; g1 = 0 ; // row index of first queued-up Householder k1 = 0 ; // pending Householders are in F (g1:t, k1:k2-1) k2 = 0 ; V = F ; // Householder vectors start here g = 0 ; // number of good Householders found nv = 0 ; // number of Householder reflections queued up vzeros = 0 ; // number of explicit zeros in queued-up H's t = 0 ; // staircase of current column fchunk = MAX (fchunk, 1) ; minchunk = MAX (MINCHUNK, fchunk/MINCHUNK_RATIO) ; rank = MIN (m,npiv) ; // F (rank,npiv) is the first entry in C. rank // is also the number of rows in the R block, // and the number of good pivot columns found. ntol = MIN (ntol, npiv) ; // Note ntol can be negative, which means to // not use tol at all. PR (("Front %ld by %ld with %ld pivots\n", m, n, npiv)) ; for (k = 0 ; k < n ; k++) { // --------------------------------------------------------------------- // reduce the kth column of F to eliminate all but "diagonal" F (g,k) // --------------------------------------------------------------------- // get the staircase for the kth column, and operate on the "diagonal" // F (g,k); eliminate F (g+1:t-1, k) to zero t0 = t ; // t0 = staircase of column k-1 t = Stair [k] ; // t = staircase of this column k PR (("k %ld g %ld m %ld n %ld npiv %ld\n", k, g, m, n, npiv)) ; if (g >= m) { // F (g,k) is outside the matrix, so we're done. If this happens // when k < npiv, then there is no contribution block. PR (("hit the wall, npiv: %ld k %ld rank %ld\n", npiv, k, rank)) ; for ( ; k < npiv ; k++) { Rdead [k] = 1 ; Stair [k] = 0 ; // remaining pivot columns all dead Tau [k] = 0 ; } for ( ; k < n ; k++) { Stair [k] = m ; // non-pivotal columns Tau [k] = 0 ; } ASSERT (nv == 0) ; // there can be no pending updates return (rank) ; } // if t < g+1, then this column is all zero; fix staircase so that R is // always inside the staircase. t = MAX (g+1,t) ; Stair [k] = t ; // --------------------------------------------------------------------- // If t just grew a lot since the last t, apply H now to all of F // --------------------------------------------------------------------- // vzeros is the number of zero entries in V after including the next // Householder vector. If it would exceed 50% of the size of V, // apply the pending Householder reflections now, but only if // enough vectors have accumulated. vzeros += nv * (t - t0) ; if (nv >= minchunk) { vsize = (nv*(nv+1))/2 + nv*(t-g1-nv) ; if (vzeros > MAX (16, vsize/2)) { // apply pending block of Householder reflections PR (("(1) apply k1 %ld k2 %ld\n", k1, k2)) ; spqr_larftb ( 0, // method 0: Left, Transpose t0-g1, n-k2, nv, m, m, V, // F (g1:t-1, k1:k1+nv-1) &Tau [k1], // Tau (k1:k-1) &F [INDEX (g1,k2,m)], // F (g1:t-1, k2:n-1) W, cc) ; // size nv*nv + nv*(n-k2) nv = 0 ; // clear queued-up Householder reflections vzeros = 0 ; } } // --------------------------------------------------------------------- // find a Householder reflection that reduces column k // --------------------------------------------------------------------- tau = spqr_private_house (t-g, &F [INDEX (g,k,m)], cc) ; // --------------------------------------------------------------------- // check to see if the kth column is OK // --------------------------------------------------------------------- if (k < ntol && (wk = spqr_abs (F [INDEX (g,k,m)], cc)) <= tol) { // ----------------------------------------------------------------- // norm (F (g:t-1, k)) is too tiny; the kth pivot column is dead // ----------------------------------------------------------------- // keep track of the 2-norm of w, the dead column 2-norms if (wk != 0) { // see also LAPACK's dnrm2 function if ((*wscale) == 0) { // this is the nonzero first entry in w (*wssq) = 1 ; } if ((*wscale) < wk) { double rr = (*wscale) / wk ; (*wssq) = 1 + (*wssq) * rr * rr ; (*wscale) = wk ; } else { double rr = wk / (*wscale) ; (*wssq) += rr * rr ; } } // zero out F (g:m-1,k) and flag it as dead for (i = g ; i < m ; i++) { // This is not strictly necessary. On output, entries outside // the staircase are ignored. F [INDEX (i,k,m)] = 0 ; } Stair [k] = 0 ; Tau [k] = 0 ; Rdead [k] = 1 ; if (nv > 0) { // apply pending block of Householder reflections PR (("(2) apply k1 %ld k2 %ld\n", k1, k2)) ; spqr_larftb ( 0, // method 0: Left, Transpose t0-g1, n-k2, nv, m, m, V, // F (g1:t-1, k1:k1+nv-1) &Tau [k1], // Tau (k1:k-1) &F [INDEX (g1,k2,m)], // F (g1:t-1, k2:n-1) W, cc) ; // size nv*nv + nv*(n-k2) nv = 0 ; // clear queued-up Householder reflections vzeros = 0 ; } } else { // ----------------------------------------------------------------- // one more good pivot column found // ----------------------------------------------------------------- Tau [k] = tau ; // save the Householder coefficient if (nv == 0) { // start the queue of pending Householder updates, and define // the current panel as k1:k2-1 g1 = g ; // first row of V k1 = k ; // first column of V k2 = MIN (n, k+fchunk) ; // k2-1 is last col in panel V = &F [INDEX (g1,k1,m)] ; // pending V starts here // check for switch to unblocked code mleft = m-g1 ; // number of rows left nleft = n-k1 ; // number of columns left if (mleft * (nleft-(fchunk+4)) < SMALL || mleft <= fchunk/2 || fchunk <= 1) { // remaining matrix is small; switch to unblocked code by // including the rest of the matrix in the panel. Always // use unblocked code if fchunk <= 1. k2 = n ; } } nv++ ; // one more pending update; V is F (g1:t-1, k1:k1+nv-1) // ----------------------------------------------------------------- // keep track of "pure" flops, for performance testing only // ----------------------------------------------------------------- // The Householder vector is of length t-g, including the unit // diagonal, and takes 3*(t-g) flops to compute. It will be // applied as a block, but compute the "pure" flops by assuming // that this single Householder vector is computed and then applied // just by itself to the rest of the frontal matrix (columns // k+1:n-1, or n-k-1 columns). Applying the Householder reflection // to just one column takes 4*(t-g) flops. This computation only // works if TBB is disabled, merely because it uses a global // variable to keep track of the flop count. If TBB is used, this // computation may result in a race condition; it is disabled in // that case. FLOP_COUNT ((t-g) * (3 + 4 * (n-k-1))) ; // ----------------------------------------------------------------- // apply the kth Householder reflection to the current panel // ----------------------------------------------------------------- // F (g:t-1, k+1:k2-1) -= v * tau * v' * F (g:t-1, k+1:k2-1), where // v is stored in F (g:t-1,k). This applies just one reflection // to the current panel. PR (("apply 1: k %ld\n", k)) ; spqr_private_apply1 (t-g, k2-k-1, m, &F [INDEX (g,k,m)], tau, &F [INDEX (g,k+1,m)], W, cc) ; g++ ; // one more pivot found // ----------------------------------------------------------------- // apply the Householder reflections if end of panel reached // ----------------------------------------------------------------- if (k == k2-1 || g == m) // or if last pivot is found { // apply pending block of Householder reflections PR (("(3) apply k1 %ld k2 %ld\n", k1, k2)) ; spqr_larftb ( 0, // method 0: Left, Transpose t-g1, n-k2, nv, m, m, V, // F (g1:t-1, k1:k1+nv-1) &Tau [k1], // Tau (k1:k2-1) &F [INDEX (g1,k2,m)], // F (g1:t-1, k2:n-1) W, cc) ; // size nv*nv + nv*(n-k2) nv = 0 ; // clear queued-up Householder reflections vzeros = 0 ; } } // --------------------------------------------------------------------- // determine the rank of the pivot columns // --------------------------------------------------------------------- if (k == npiv-1) { // the rank is the number of good columns found in the first // npiv columns. It is also the number of rows in the R block. // F (rank,npiv) is first entry in the C block. rank = g ; PR (("rank of Front pivcols: %ld\n", rank)) ; } } if (CHECK_BLAS_INT && !cc->blas_ok) { // This cannot occur if the BLAS_INT and the Long are the same integer. // In that case, CHECK_BLAS_INT is FALSE at compile-time, and the // compiler will then remove this as dead code. ERROR (CHOLMOD_INVALID, "problem too large for the BLAS") ; return (0) ; } return (rank) ; }
template <typename Entry> int spqr_1colamd // TRUE if OK, FALSE otherwise ( // inputs, not modified int ordering, // all available, except 0:fixed and 3:given // treated as 1:natural double tol, // only accept singletons above tol Long bncols, // number of columns of B cholmod_sparse *A, // m-by-n sparse matrix // outputs, neither allocated nor defined on input Long **p_Q1fill, // size n+bncols, fill-reducing // or natural ordering Long **p_R1p, // size n1rows+1, R1p [k] = # of nonzeros in kth // row of R1. NULL if n1cols == 0. Long **p_P1inv, // size m, singleton row inverse permutation. // If row i of A is the kth singleton row, then // P1inv [i] = k. NULL if n1cols is zero. cholmod_sparse **p_Y, // on output, only the first n-n1cols+1 entries of // Y->p are defined (if Y is not NULL), where // Y = [A B] or Y = [A2 B2]. If B is empty and // there are no column singletons, Y is NULL Long *p_n1cols, // number of column singletons found Long *p_n1rows, // number of corresponding rows found // workspace and parameters cholmod_common *cc ) { Long *Q1fill, *Degree, *Qrows, *W, *Winv, *ATp, *ATj, *R1p, *P1inv, *Yp, *Ap, *Ai, *Work ; Entry *Ax ; Long p, d, j, i, k, n1cols, n1rows, row, col, pend, n2rows, n2cols = EMPTY, nz2, kk, p2, col2, ynz, fill_reducing_ordering, m, n, xtype, worksize ; cholmod_sparse *AT, *Y ; // ------------------------------------------------------------------------- // get inputs // ------------------------------------------------------------------------- xtype = spqr_type <Entry> ( ) ; m = A->nrow ; n = A->ncol ; Ap = (Long *) A->p ; Ai = (Long *) A->i ; Ax = (Entry *) A->x ; // set outputs to NULL in case of early return *p_Q1fill = NULL ; *p_R1p = NULL ; *p_P1inv = NULL ; *p_Y = NULL ; *p_n1cols = EMPTY ; *p_n1rows = EMPTY ; // ------------------------------------------------------------------------- // allocate result Q1fill (Y, R1p, P1inv allocated later) // ------------------------------------------------------------------------- Q1fill = (Long *) cholmod_l_malloc (n+bncols, sizeof (Long), cc) ; // ------------------------------------------------------------------------- // allocate workspace // ------------------------------------------------------------------------- fill_reducing_ordering = ! ((ordering == SPQR_ORDERING_FIXED) || (ordering == SPQR_ORDERING_GIVEN) || (ordering == SPQR_ORDERING_NATURAL)) ; worksize = ((fill_reducing_ordering) ? 3:2) * n ; Work = (Long *) cholmod_l_malloc (worksize, sizeof (Long), cc) ; Degree = Work ; // size n Qrows = Work + n ; // size n Winv = Qrows ; // Winv and Qrows not needed at the same time W = Qrows + n ; // size n if fill-reducing ordering, else size 0 if (cc->status < CHOLMOD_OK) { // out of memory; free everything and return cholmod_l_free (worksize, sizeof (Long), Work, cc) ; cholmod_l_free (n+bncols, sizeof (Long), Q1fill, cc) ; return (FALSE) ; } // ------------------------------------------------------------------------- // initialze queue with empty columns, and columns with just one entry // ------------------------------------------------------------------------- n1cols = 0 ; n1rows = 0 ; for (j = 0 ; j < n ; j++) { p = Ap [j] ; d = Ap [j+1] - p ; if (d == 0) { // j is a dead column singleton PR (("initial dead %ld\n", j)) ; Q1fill [n1cols] = j ; Qrows [n1cols] = EMPTY ; n1cols++ ; Degree [j] = EMPTY ; } else if (d == 1 && spqr_abs (Ax [p], cc) > tol) { // j is a column singleton, live or dead PR (("initial live %ld %ld\n", j, Ai [p])) ; Q1fill [n1cols] = j ; Qrows [n1cols] = Ai [p] ; // this might be a duplicate n1cols++ ; Degree [j] = EMPTY ; } else { // j has degree > 1, it is not (yet) a singleton Degree [j] = d ; } } // Degree [j] = EMPTY if j is in the singleton queue, or the Degree [j] > 1 // is the degree of column j otherwise // ------------------------------------------------------------------------- // create AT = spones (A') // ------------------------------------------------------------------------- AT = cholmod_l_transpose (A, 0, cc) ; // [ if (cc->status < CHOLMOD_OK) { // out of memory; free everything and return cholmod_l_free (worksize, sizeof (Long), Work, cc) ; cholmod_l_free (n+bncols, sizeof (Long), Q1fill, cc) ; return (FALSE) ; } ATp = (Long *) AT->p ; ATj = (Long *) AT->i ; // ------------------------------------------------------------------------- // remove column singletons via breadth-first-search // ------------------------------------------------------------------------- for (k = 0 ; k < n1cols ; k++) { // --------------------------------------------------------------------- // get a new singleton from the queue // --------------------------------------------------------------------- col = Q1fill [k] ; row = Qrows [k] ; PR (("\n---- singleton col %ld row %ld\n", col, row)) ; ASSERT (Degree [col] == EMPTY) ; if (row == EMPTY || ATp [row] < 0) { // ----------------------------------------------------------------- // col is a dead column singleton; remove duplicate row index // ----------------------------------------------------------------- Qrows [k] = EMPTY ; row = EMPTY ; PR (("dead: %ld\n", col)) ; } else { // ----------------------------------------------------------------- // col is a live col singleton; remove its row from matrix // ----------------------------------------------------------------- n1rows++ ; p = ATp [row] ; ATp [row] = FLIP (p) ; // flag the singleton row pend = UNFLIP (ATp [row+1]) ; PR (("live: %ld row %ld\n", col, row)) ; for ( ; p < pend ; p++) { // look for new column singletons after row is removed j = ATj [p] ; d = Degree [j] ; if (d == EMPTY) { // j is already in the singleton queue continue ; } ASSERT (d >= 1) ; ASSERT2 (spqrDebug_listcount (j, Q1fill, n1cols, 0) == 0) ; d-- ; Degree [j] = d ; if (d == 0) { // a new dead col singleton PR (("newly dead %ld\n", j)) ; Q1fill [n1cols] = j ; Qrows [n1cols] = EMPTY ; n1cols++ ; Degree [j] = EMPTY ; } else if (d == 1) { // a new live col singleton; find its single live row for (p2 = Ap [j] ; p2 < Ap [j+1] ; p2++) { i = Ai [p2] ; if (ATp [i] >= 0 && spqr_abs (Ax [p2], cc) > tol) { // i might appear in Qrows [k+1:n1cols-1] PR (("newly live %ld\n", j)) ; ASSERT2 (spqrDebug_listcount (i,Qrows,k+1,1) == 0) ; Q1fill [n1cols] = j ; Qrows [n1cols] = i ; n1cols++ ; Degree [j] = EMPTY ; break ; } } } } } // Q1fill [0:k] and Qrows [0:k] have no duplicates ASSERT2 (spqrDebug_listcount (col, Q1fill, n1cols, 0) == 1) ; ASSERT2 (IMPLIES (row >= 0, spqrDebug_listcount (row, Qrows, k+1, 1) == 1)) ; } // ------------------------------------------------------------------------- // Degree flags the column singletons, ATp flags their rows // ------------------------------------------------------------------------- #ifndef NDEBUG k = 0 ; for (j = 0 ; j < n ; j++) { PR (("j %ld Degree[j] %ld\n", j, Degree [j])) ; if (Degree [j] > 0) k++ ; // j is not a column singleton } PR (("k %ld n %ld n1cols %ld\n", k, n, n1cols)) ; ASSERT (k == n - n1cols) ; for (k = 0 ; k < n1cols ; k++) { col = Q1fill [k] ; ASSERT (Degree [col] <= 0) ; } k = 0 ; for (i = 0 ; i < m ; i++) { if (ATp [i] >= 0) k++ ; // i is not a row of a col singleton } ASSERT (k == m - n1rows) ; for (k = 0 ; k < n1cols ; k++) { row = Qrows [k] ; ASSERT (IMPLIES (row != EMPTY, ATp [row] < 0)) ; } #endif // ------------------------------------------------------------------------- // find the row ordering // ------------------------------------------------------------------------- if (n1cols == 0) { // --------------------------------------------------------------------- // no singletons in the matrix; no R1 matrix, no P1inv permutation // --------------------------------------------------------------------- ASSERT (n1rows == 0) ; R1p = NULL ; P1inv = NULL ; } else { // --------------------------------------------------------------------- // construct the row singleton permutation // --------------------------------------------------------------------- // allocate result arrays R1p and P1inv R1p = (Long *) cholmod_l_malloc (n1rows+1, sizeof (Long), cc) ; P1inv = (Long *) cholmod_l_malloc (m, sizeof (Long), cc) ; if (cc->status < CHOLMOD_OK) { // out of memory; free everything and return cholmod_l_free_sparse (&AT, cc) ; cholmod_l_free (worksize, sizeof (Long), Work, cc) ; cholmod_l_free (n+bncols, sizeof (Long), Q1fill, cc) ; cholmod_l_free (n1rows+1, sizeof (Long), R1p, cc) ; cholmod_l_free (m, sizeof (Long), P1inv, cc) ; return (FALSE) ; } #ifndef NDEBUG for (i = 0 ; i < m ; i++) P1inv [i] = EMPTY ; #endif kk = 0 ; for (k = 0 ; k < n1cols ; k++) { i = Qrows [k] ; PR (("singleton col %ld row %ld\n", Q1fill [k], i)) ; if (i != EMPTY) { // row i is the kk-th singleton row ASSERT (ATp [i] < 0) ; ASSERT (P1inv [i] == EMPTY) ; P1inv [i] = kk ; // also find # of entries in row kk of R1 R1p [kk] = UNFLIP (ATp [i+1]) - UNFLIP (ATp [i]) ; kk++ ; } } ASSERT (kk == n1rows) ; for (i = 0 ; i < m ; i++) { if (ATp [i] >= 0) { // row i is not a singleton row ASSERT (P1inv [i] == EMPTY) ; P1inv [i] = kk ; kk++ ; } } ASSERT (kk == m) ; } // Qrows is no longer needed. // ------------------------------------------------------------------------- // complete the column ordering // ------------------------------------------------------------------------- if (!fill_reducing_ordering) { // --------------------------------------------------------------------- // natural ordering // --------------------------------------------------------------------- if (n1cols == 0) { // no singletons, so natural ordering is 0:n-1 for now for (k = 0 ; k < n ; k++) { Q1fill [k] = k ; } } else { // singleton columns appear first, then non column singletons k = n1cols ; for (j = 0 ; j < n ; j++) { if (Degree [j] > 0) { // column j is not a column singleton Q1fill [k++] = j ; } } ASSERT (k == n) ; } } else { // --------------------------------------------------------------------- // fill-reducing ordering of pruned submatrix // --------------------------------------------------------------------- if (n1cols == 0) { // ----------------------------------------------------------------- // no singletons found; do fill-reducing on entire matrix // ----------------------------------------------------------------- n2cols = n ; n2rows = m ; } else { // ----------------------------------------------------------------- // create the pruned matrix for fill-reducing by removing singletons // ----------------------------------------------------------------- // find the mapping of original columns to pruned columns n2cols = 0 ; for (j = 0 ; j < n ; j++) { if (Degree [j] > 0) { // column j is not a column singleton W [j] = n2cols++ ; PR (("W [%ld] = %ld\n", j, W [j])) ; } else { // column j is a column singleton W [j] = EMPTY ; PR (("W [%ld] = %ld (j is col singleton)\n", j, W [j])) ; } } ASSERT (n2cols == n - n1cols) ; // W is now a mapping of the original columns to the columns in the // pruned matrix. W [col] == EMPTY if col is a column singleton. // Otherwise col2 = W [j] is a column of the pruned matrix. // ----------------------------------------------------------------- // delete row and column singletons from A' // ----------------------------------------------------------------- // compact A' by removing row and column singletons nz2 = 0 ; n2rows = 0 ; for (i = 0 ; i < m ; i++) { p = ATp [i] ; if (p >= 0) { // row i is not a row of a column singleton ATp [n2rows++] = nz2 ; pend = UNFLIP (ATp [i+1]) ; for (p = ATp [i] ; p < pend ; p++) { j = ATj [p] ; ASSERT (W [j] >= 0 && W [j] < n-n1cols) ; ATj [nz2++] = W [j] ; } } } ATp [n2rows] = nz2 ; ASSERT (n2rows == m - n1rows) ; } // --------------------------------------------------------------------- // fill-reducing ordering of the transpose of the pruned A' matrix // --------------------------------------------------------------------- PR (("n1cols %ld n1rows %ld n2cols %ld n2rows %ld\n", n1cols, n1rows, n2cols, n2rows)) ; ASSERT ((Long) AT->nrow == n) ; ASSERT ((Long) AT->ncol == m) ; AT->nrow = n2cols ; AT->ncol = n2rows ; // save the current CHOLMOD settings Long save [6] ; save [0] = cc->supernodal ; save [1] = cc->nmethods ; save [2] = cc->postorder ; save [3] = cc->method [0].ordering ; save [4] = cc->method [1].ordering ; save [5] = cc->method [2].ordering ; // follow the ordering with a postordering of the column etree cc->postorder = TRUE ; // 8:best: best of COLAMD(A), AMD(A'A), and METIS (if available) if (ordering == SPQR_ORDERING_BEST) { ordering = SPQR_ORDERING_CHOLMOD ; cc->nmethods = 2 ; cc->method [0].ordering = CHOLMOD_COLAMD ; cc->method [1].ordering = CHOLMOD_AMD ; #ifndef NPARTITION cc->nmethods = 3 ; cc->method [2].ordering = CHOLMOD_METIS ; #endif } // 9:bestamd: best of COLAMD(A) and AMD(A'A) if (ordering == SPQR_ORDERING_BESTAMD) { // if METIS is not installed, this option is the same as 8:best ordering = SPQR_ORDERING_CHOLMOD ; cc->nmethods = 2 ; cc->method [0].ordering = CHOLMOD_COLAMD ; cc->method [1].ordering = CHOLMOD_AMD ; } #ifdef NPARTITION if (ordering == SPQR_ORDERING_METIS) { // METIS not installed; use default ordering ordering = SPQR_ORDERING_DEFAULT ; } #endif if (ordering == SPQR_ORDERING_DEFAULT) { // Version 1.2.0: just use COLAMD ordering = SPQR_ORDERING_COLAMD ; #if 0 // Version 1.1.2 and earlier: if (n2rows <= 2*n2cols) { // just use COLAMD; do not try AMD or METIS ordering = SPQR_ORDERING_COLAMD ; } else { #ifndef NPARTITION // use CHOLMOD's default ordering: try AMD and then METIS // if AMD gives high fill-in, and take the best ordering found ordering = SPQR_ORDERING_CHOLMOD ; cc->nmethods = 0 ; #else // METIS is not installed, so just use AMD ordering = SPQR_ORDERING_AMD ; #endif } #endif } if (ordering == SPQR_ORDERING_AMD) { // use CHOLMOD's interface to AMD to order A'*A cholmod_l_amd (AT, NULL, 0, (Long *) (Q1fill + n1cols), cc) ; } #ifndef NPARTITION else if (ordering == SPQR_ORDERING_METIS) { // use CHOLMOD's interface to METIS to order A'*A (if installed) cholmod_l_metis (AT, NULL, 0, TRUE, (Long *) (Q1fill + n1cols), cc) ; } #endif else if (ordering == SPQR_ORDERING_CHOLMOD) { // use CHOLMOD's internal ordering (defined by cc) to order AT PR (("Using CHOLMOD, nmethods %d\n", cc->nmethods)) ; cc->supernodal = CHOLMOD_SIMPLICIAL ; cc->postorder = TRUE ; cholmod_factor *Sc ; Sc = cholmod_l_analyze_p2 (FALSE, AT, NULL, NULL, 0, cc) ; if (Sc != NULL) { // copy perm from Sc->Perm [0:n2cols-1] to Q1fill (n1cols:n) Long *Sc_perm = (Long *) Sc->Perm ; for (k = 0 ; k < n2cols ; k++) { Q1fill [k + n1cols] = Sc_perm [k] ; } // CHOLMOD selected an ordering; determine the ordering used switch (Sc->ordering) { case CHOLMOD_AMD: ordering = SPQR_ORDERING_AMD ;break; case CHOLMOD_COLAMD: ordering = SPQR_ORDERING_COLAMD ;break; case CHOLMOD_METIS: ordering = SPQR_ORDERING_METIS ;break; } } cholmod_l_free_factor (&Sc, cc) ; PR (("CHOLMOD used method %d : ordering: %d\n", cc->selected, cc->method [cc->selected].ordering)) ; } else // SPQR_ORDERING_DEFAULT or SPQR_ORDERING_COLAMD { // use CHOLMOD's interface to COLAMD to order AT ordering = SPQR_ORDERING_COLAMD ; cholmod_l_colamd (AT, NULL, 0, TRUE, (Long *) (Q1fill + n1cols), cc) ; } cc->SPQR_istat [7] = ordering ; // restore the CHOLMOD settings cc->supernodal = save [0] ; cc->nmethods = save [1] ; cc->postorder = save [2] ; cc->method [0].ordering = save [3] ; cc->method [1].ordering = save [4] ; cc->method [2].ordering = save [5] ; AT->nrow = n ; AT->ncol = m ; } // ------------------------------------------------------------------------- // free AT // ------------------------------------------------------------------------- cholmod_l_free_sparse (&AT, cc) ; // ] // ------------------------------------------------------------------------- // check if the method succeeded // ------------------------------------------------------------------------- if (cc->status < CHOLMOD_OK) { // out of memory; free everything and return cholmod_l_free (worksize, sizeof (Long), Work, cc) ; cholmod_l_free (n+bncols, sizeof (Long), Q1fill, cc) ; cholmod_l_free (n1rows+1, sizeof (Long), R1p, cc) ; cholmod_l_free (m, sizeof (Long), P1inv, cc) ; return (FALSE) ; } // ------------------------------------------------------------------------- // map the fill-reducing ordering ordering back to A // ------------------------------------------------------------------------- if (n1cols > 0 && fill_reducing_ordering) { // Winv is workspace of size n2cols <= n #ifndef NDEBUG for (j = 0 ; j < n2cols ; j++) Winv [j] = EMPTY ; #endif for (j = 0 ; j < n ; j++) { // j is a column of A. col2 = W [j] is either EMPTY, or it is // the corresponding column of the pruned matrix col2 = W [j] ; if (col2 != EMPTY) { ASSERT (col2 >= 0 && col2 < n2cols) ; Winv [col2] = j ; } } for (k = n1cols ; k < n ; k++) { // col2 is a column of the pruned matrix col2 = Q1fill [k] ; // j is the corresonding column of the A matrix j = Winv [col2] ; ASSERT (j >= 0 && j < n) ; Q1fill [k] = j ; } } // ------------------------------------------------------------------------- // identity permutation of the columns of B // ------------------------------------------------------------------------- for (k = n ; k < n+bncols ; k++) { // tack on the identity permutation for columns of B Q1fill [k] = k ; } // ------------------------------------------------------------------------- // find column pointers for Y = [A2 B2]; columns of A2 // ------------------------------------------------------------------------- if (n1cols == 0 && bncols == 0) { // A will be factorized instead of Y Y = NULL ; } else { // Y has no entries yet; nnz(Y) will be determined later Y = cholmod_l_allocate_sparse (m-n1rows, n-n1cols+bncols, 0, FALSE, TRUE, 0, xtype, cc) ; if (cc->status < CHOLMOD_OK) { // out of memory; free everything and return cholmod_l_free (worksize, sizeof (Long), Work, cc) ; cholmod_l_free (n+bncols, sizeof (Long), Q1fill, cc) ; cholmod_l_free (n1rows+1, sizeof (Long), R1p, cc) ; cholmod_l_free (m, sizeof (Long), P1inv, cc) ; return (FALSE) ; } Yp = (Long *) Y->p ; ynz = 0 ; PR (("1c wrapup: n1cols %ld n %ld\n", n1cols, n)) ; for (k = n1cols ; k < n ; k++) { j = Q1fill [k] ; d = Degree [j] ; ASSERT (d >= 1 && d <= m) ; Yp [k-n1cols] = ynz ; ynz += d ; } Yp [n-n1cols] = ynz ; } // ------------------------------------------------------------------------- // free workspace and return results // ------------------------------------------------------------------------- cholmod_l_free (worksize, sizeof (Long), Work, cc) ; *p_Q1fill = Q1fill ; *p_R1p = R1p ; *p_P1inv = P1inv ; *p_Y = Y ; *p_n1cols = n1cols ; *p_n1rows = n1rows ; return (TRUE) ; }
template <typename Entry> int spqr_1fixed ( // inputs, not modified double tol, // only accept singletons above tol Long bncols, // number of columns of B cholmod_sparse *A, // m-by-n sparse matrix // output arrays, neither allocated nor defined on input. Long **p_R1p, // size n1rows+1, R1p [k] = # of nonzeros in kth // row of R1. NULL if n1cols == 0. Long **p_P1inv, // size m, singleton row inverse permutation. // If row i of A is the kth singleton row, then // P1inv [i] = k. NULL if n1cols is zero. cholmod_sparse **p_Y, // on output, only the first n-n1cols+1 entries of // Y->p are defined (if Y is not NULL), where // Y = [A B] or Y = [A2 B2]. If B is empty and // there are no column singletons, Y is NULL Long *p_n1cols, // number of column singletons found Long *p_n1rows, // number of corresponding rows found // workspace and parameters cholmod_common *cc ) { cholmod_sparse *Y ; Long *P1inv, *R1p, *Yp, *Qrows, *Ap, *Ai ; char *Mark ; Entry *Ax ; Long i, j, k, p, d, row, n1rows, n1cols, ynz, iold, inew, kk, m, n, xtype ; // ------------------------------------------------------------------------- // get inputs // ------------------------------------------------------------------------- xtype = spqr_type <Entry> ( ) ; m = A->nrow ; n = A->ncol ; Ap = (Long *) A->p ; Ai = (Long *) A->i ; Ax = (Entry *) A->x ; // set outputs to NULL in case of early return *p_R1p = NULL ; *p_P1inv = NULL ; *p_Y = NULL ; *p_n1cols = EMPTY ; *p_n1rows = EMPTY ; // ------------------------------------------------------------------------- // allocate workspace // ------------------------------------------------------------------------- Mark = (char *) cholmod_l_calloc (m, sizeof (char), cc) ; Qrows = (Long *) cholmod_l_malloc (n, sizeof (Long), cc) ; if (cc->status < CHOLMOD_OK) { // out of memory cholmod_l_free (m, sizeof (char), Mark, cc) ; cholmod_l_free (n, sizeof (Long), Qrows, cc) ; return (FALSE) ; } // ------------------------------------------------------------------------- // find singletons; no column permutations allowed // ------------------------------------------------------------------------- n1cols = 0 ; // number of column singletons found n1rows = 0 ; // number of corresponding singleton rows for (j = 0 ; j < n ; j++) { // count the number of unmarked rows in column j Entry aij = 0 ; d = 0 ; row = EMPTY ; for (p = Ap [j] ; d < 2 && p < Ap [j+1] ; p++) { i = Ai [p] ; if (!Mark [i]) { // row i is not taken by a prior column singleton. If this // is the only unflagged row and the value is large enough, // it will become the row for this column singleton. aij = Ax [p] ; row = i ; d++ ; } } if (d == 0) { // j is a dead column singleton Qrows [n1cols++] = EMPTY ; } else if (d == 1 && spqr_abs (aij, cc) > tol) { // j is a live column singleton Qrows [n1cols++] = row ; // flag row i as taken Mark [row] = TRUE ; n1rows++ ; } else { // j is not a singleton; quit searching break ; } } // ------------------------------------------------------------------------- // construct P1inv permutation, row counts R1p, and col pointers Yp // ------------------------------------------------------------------------- if (n1cols == 0 && bncols == 0) { // --------------------------------------------------------------------- // no singletons, and B empty; Y=A will be done via pointer alias // --------------------------------------------------------------------- Y = NULL ; Yp = NULL ; P1inv = NULL ; R1p = NULL ; } else if (n1cols == 0) { // --------------------------------------------------------------------- // no singletons in the matrix; no R1 matrix, no P1inv permutation // --------------------------------------------------------------------- // Y has no entries yet; nnz(Y) will be determined later Y = cholmod_l_allocate_sparse (m, n+bncols, 0, FALSE, TRUE, 0, xtype, cc) ; if (cc->status < CHOLMOD_OK) { // out of memory cholmod_l_free (m, sizeof (char), Mark, cc) ; cholmod_l_free (n, sizeof (Long), Qrows, cc) ; return (FALSE) ; } Yp = (Long *) Y->p ; ASSERT (n1rows == 0) ; P1inv = NULL ; R1p = NULL ; // --------------------------------------------------------------------- // copy the column pointers of A for the first part of Y = [A B] // --------------------------------------------------------------------- ynz = Ap [n] ; for (k = 0 ; k <= n ; k++) { Yp [k] = Ap [k] ; } } else { // --------------------------------------------------------------------- // construct the row singleton permutation // --------------------------------------------------------------------- // Y has no entries yet; nnz(Y) will be determined later Y = cholmod_l_allocate_sparse (m-n1rows, n-n1cols+bncols, 0, TRUE, TRUE, 0, xtype, cc) ; P1inv = (Long *) cholmod_l_malloc (m, sizeof (Long), cc) ; R1p = (Long *) cholmod_l_calloc (n1rows+1, sizeof (Long), cc) ; if (cc->status < CHOLMOD_OK) { // out of memory cholmod_l_free_sparse (&Y, cc) ; cholmod_l_free (m, sizeof (Long), P1inv, cc) ; cholmod_l_free (n1rows+1, sizeof (Long), R1p, cc) ; cholmod_l_free (m, sizeof (char), Mark, cc) ; cholmod_l_free (n, sizeof (Long), Qrows, cc) ; return (FALSE) ; } Yp = (Long *) Y->p ; #ifndef NDEBUG for (i = 0 ; i < m ; i++) P1inv [i] = EMPTY ; #endif kk = 0 ; for (k = 0 ; k < n1cols ; k++) { i = Qrows [k] ; if (i != EMPTY) { // row i is the kk-th singleton row ASSERT (Mark [i]) ; ASSERT (P1inv [i] == EMPTY) ; P1inv [i] = kk ; kk++ ; } } for (i = 0 ; i < m ; i++) { if (!Mark [i]) { // row i is not a singleton row ASSERT (P1inv [i] == EMPTY) ; P1inv [i] = kk ; kk++ ; } } ASSERT (kk == m) ; // --------------------------------------------------------------------- // find row counts for R11 // --------------------------------------------------------------------- for (k = 0 ; k < n1cols ; k++) { for (p = Ap [k] ; p < Ap [k+1] ; p++) { iold = Ai [p] ; inew = P1inv [iold] ; ASSERT (inew < n1rows) ; R1p [inew]++ ; // a singleton row; in R1 } } // --------------------------------------------------------------------- // find row counts for R12 and column pointers for A2 part of Y // --------------------------------------------------------------------- ynz = 0 ; for ( ; k < n ; k++) { Yp [k-n1cols] = ynz ; for (p = Ap [k] ; p < Ap [k+1] ; p++) { iold = Ai [p] ; inew = P1inv [iold] ; if (inew < n1rows) { R1p [inew]++ ; // a singleton row; in R1 } else { ynz++ ; // not a singleton row; in A2 } } } Yp [n-n1cols] = ynz ; #ifndef NDEBUG PR (("n1cols: %ld\n", n1cols)) ; for (i = 0 ; i < n1rows ; i++) { PR (("R1p [%ld] is %ld\n", i, R1p [i])) ; ASSERT (R1p [i] > 0) ; } #endif } // ------------------------------------------------------------------------- // free workspace and return results // ------------------------------------------------------------------------- cholmod_l_free (n, sizeof (Long), Qrows, cc) ; cholmod_l_free (m, sizeof (char), Mark, cc) ; *p_R1p = R1p ; *p_P1inv = P1inv ; *p_Y = Y ; *p_n1cols = n1cols ; *p_n1rows = n1rows ; return (TRUE) ; }