static Int analyze_worker /* returns KLU_OK or < 0 if error */ ( /* inputs, not modified */ Int n, /* A is n-by-n */ Int Ap [ ], /* size n+1, column pointers */ Int Ai [ ], /* size nz, row indices */ Int nblocks, /* # of blocks */ Int Pbtf [ ], /* BTF row permutation */ Int Qbtf [ ], /* BTF col permutation */ Int R [ ], /* size n+1, but only Rbtf [0..nblocks] is used */ Int ordering, /* what ordering to use (0, 1, or 3 for this routine) */ /* output only, not defined on input */ Int P [ ], /* size n */ Int Q [ ], /* size n */ double Lnz [ ], /* size n, but only Lnz [0..nblocks-1] is used */ /* workspace, not defined on input or output */ Int Pblk [ ], /* size maxblock */ Int Cp [ ], /* size maxblock+1 */ Int Ci [ ], /* size MAX (nz+1, Cilen) */ Int Cilen, /* nz+1, or COLAMD_recommend(nz,n,n) for COLAMD */ Int Pinv [ ], /* size maxblock */ /* input/output */ KLU_symbolic *Symbolic, KLU_common *Common ) { double amd_Info [AMD_INFO], lnz, lnz1, flops, flops1 ; Int k1, k2, nk, k, block, oldcol, pend, newcol, result, pc, p, newrow, maxnz, nzoff, cstats [COLAMD_STATS], ok, err = KLU_INVALID ; /* ---------------------------------------------------------------------- */ /* initializations */ /* ---------------------------------------------------------------------- */ /* compute the inverse of Pbtf */ #ifndef NDEBUG for (k = 0 ; k < n ; k++) { P [k] = EMPTY ; Q [k] = EMPTY ; Pinv [k] = EMPTY ; } #endif for (k = 0 ; k < n ; k++) { ASSERT (Pbtf [k] >= 0 && Pbtf [k] < n) ; Pinv [Pbtf [k]] = k ; } #ifndef NDEBUG for (k = 0 ; k < n ; k++) ASSERT (Pinv [k] != EMPTY) ; #endif nzoff = 0 ; lnz = 0 ; maxnz = 0 ; flops = 0 ; Symbolic->symmetry = EMPTY ; /* only computed by AMD */ /* ---------------------------------------------------------------------- */ /* order each block */ /* ---------------------------------------------------------------------- */ for (block = 0 ; block < nblocks ; block++) { /* ------------------------------------------------------------------ */ /* the block is from rows/columns k1 to k2-1 */ /* ------------------------------------------------------------------ */ k1 = R [block] ; k2 = R [block+1] ; nk = k2 - k1 ; PRINTF (("BLOCK %d, k1 %d k2-1 %d nk %d\n", block, k1, k2-1, nk)) ; /* ------------------------------------------------------------------ */ /* construct the kth block, C */ /* ------------------------------------------------------------------ */ Lnz [block] = EMPTY ; pc = 0 ; for (k = k1 ; k < k2 ; k++) { newcol = k-k1 ; Cp [newcol] = pc ; oldcol = Qbtf [k] ; pend = Ap [oldcol+1] ; for (p = Ap [oldcol] ; p < pend ; p++) { newrow = Pinv [Ai [p]] ; if (newrow < k1) { nzoff++ ; } else { /* (newrow,newcol) is an entry in the block */ ASSERT (newrow < k2) ; newrow -= k1 ; Ci [pc++] = newrow ; } } } Cp [nk] = pc ; maxnz = MAX (maxnz, pc) ; ASSERT (KLU_valid (nk, Cp, Ci, NULL)) ; /* ------------------------------------------------------------------ */ /* order the block C */ /* ------------------------------------------------------------------ */ if (nk <= 3) { /* -------------------------------------------------------------- */ /* use natural ordering for tiny blocks (3-by-3 or less) */ /* -------------------------------------------------------------- */ for (k = 0 ; k < nk ; k++) { Pblk [k] = k ; } lnz1 = nk * (nk + 1) / 2 ; flops1 = nk * (nk - 1) / 2 + (nk-1)*nk*(2*nk-1) / 6 ; ok = TRUE ; } else if (ordering == 0) { /* -------------------------------------------------------------- */ /* order the block with AMD (C+C') */ /* -------------------------------------------------------------- */ result = AMD_order (nk, Cp, Ci, Pblk, NULL, amd_Info) ; ok = (result >= AMD_OK) ; if (result == AMD_OUT_OF_MEMORY) { err = KLU_OUT_OF_MEMORY ; } /* account for memory usage in AMD */ Common->mempeak = MAX (Common->mempeak, Common->memusage + amd_Info [AMD_MEMORY]) ; /* get the ordering statistics from AMD */ lnz1 = (Int) (amd_Info [AMD_LNZ]) + nk ; flops1 = 2 * amd_Info [AMD_NMULTSUBS_LU] + amd_Info [AMD_NDIV] ; if (pc == maxnz) { /* get the symmetry of the biggest block */ Symbolic->symmetry = amd_Info [AMD_SYMMETRY] ; } } else if (ordering == 1) { /* -------------------------------------------------------------- */ /* order the block with COLAMD (C) */ /* -------------------------------------------------------------- */ /* order (and destroy) Ci, returning column permutation in Cp. * COLAMD "cannot" fail since the matrix has already been checked, * and Ci allocated. */ ok = COLAMD (nk, nk, Cilen, Ci, Cp, NULL, cstats) ; lnz1 = EMPTY ; flops1 = EMPTY ; /* copy the permutation from Cp to Pblk */ for (k = 0 ; k < nk ; k++) { Pblk [k] = Cp [k] ; } } else { /* -------------------------------------------------------------- */ /* pass the block to the user-provided ordering function */ /* -------------------------------------------------------------- */ lnz1 = (Common->user_order) (nk, Cp, Ci, Pblk, Common) ; flops1 = EMPTY ; ok = (lnz1 != 0) ; } if (!ok) { return (err) ; /* ordering method failed */ } /* ------------------------------------------------------------------ */ /* keep track of nnz(L) and flops statistics */ /* ------------------------------------------------------------------ */ Lnz [block] = lnz1 ; lnz = (lnz == EMPTY || lnz1 == EMPTY) ? EMPTY : (lnz + lnz1) ; flops = (flops == EMPTY || flops1 == EMPTY) ? EMPTY : (flops + flops1) ; /* ------------------------------------------------------------------ */ /* combine the preordering with the BTF ordering */ /* ------------------------------------------------------------------ */ PRINTF (("Pblk, 1-based:\n")) ; for (k = 0 ; k < nk ; k++) { ASSERT (k + k1 < n) ; ASSERT (Pblk [k] + k1 < n) ; Q [k + k1] = Qbtf [Pblk [k] + k1] ; } for (k = 0 ; k < nk ; k++) { ASSERT (k + k1 < n) ; ASSERT (Pblk [k] + k1 < n) ; P [k + k1] = Pbtf [Pblk [k] + k1] ; } } PRINTF (("nzoff %d Ap[n] %d\n", nzoff, Ap [n])) ; ASSERT (nzoff >= 0 && nzoff <= Ap [n]) ; /* return estimates of # of nonzeros in L including diagonal */ Symbolic->lnz = lnz ; /* EMPTY if COLAMD used */ Symbolic->unz = lnz ; Symbolic->nzoff = nzoff ; Symbolic->est_flops = flops ; /* EMPTY if COLAMD or user-ordering used */ return (KLU_OK) ; }
Int KLU_solve ( /* inputs, not modified */ KLU_symbolic *Symbolic, KLU_numeric *Numeric, Int d, /* leading dimension of B */ Int nrhs, /* number of right-hand-sides */ /* right-hand-side on input, overwritten with solution to Ax=b on output */ double B [ ], /* size n*nrhs, in column-oriented form, with * leading dimension d. */ /* --------------- */ KLU_common *Common ) { Entry x [4], offik, s ; double rs, *Rs ; Entry *Offx, *X, *Bz, *Udiag ; Int *Q, *R, *Pnum, *Offp, *Offi, *Lip, *Uip, *Llen, *Ulen ; Unit **LUbx ; Int k1, k2, nk, k, block, pend, n, p, nblocks, chunk, nr, i ; /* ---------------------------------------------------------------------- */ /* check inputs */ /* ---------------------------------------------------------------------- */ if (Common == NULL) { return (FALSE) ; } if (Numeric == NULL || Symbolic == NULL || d < Symbolic->n || nrhs < 0 || B == NULL) { Common->status = KLU_INVALID ; return (FALSE) ; } Common->status = KLU_OK ; /* ---------------------------------------------------------------------- */ /* get the contents of the Symbolic object */ /* ---------------------------------------------------------------------- */ Bz = (Entry *) B ; n = Symbolic->n ; nblocks = Symbolic->nblocks ; Q = Symbolic->Q ; R = Symbolic->R ; /* ---------------------------------------------------------------------- */ /* get the contents of the Numeric object */ /* ---------------------------------------------------------------------- */ ASSERT (nblocks == Numeric->nblocks) ; Pnum = Numeric->Pnum ; Offp = Numeric->Offp ; Offi = Numeric->Offi ; Offx = (Entry *) Numeric->Offx ; Lip = Numeric->Lip ; Llen = Numeric->Llen ; Uip = Numeric->Uip ; Ulen = Numeric->Ulen ; LUbx = (Unit **) Numeric->LUbx ; Udiag = Numeric->Udiag ; Rs = Numeric->Rs ; X = (Entry *) Numeric->Xwork ; ASSERT (KLU_valid (n, Offp, Offi, Offx)) ; /* ---------------------------------------------------------------------- */ /* solve in chunks of 4 columns at a time */ /* ---------------------------------------------------------------------- */ for (chunk = 0 ; chunk < nrhs ; chunk += 4) { /* ------------------------------------------------------------------ */ /* get the size of the current chunk */ /* ------------------------------------------------------------------ */ nr = MIN (nrhs - chunk, 4) ; /* ------------------------------------------------------------------ */ /* scale and permute the right hand side, X = P*(R\B) */ /* ------------------------------------------------------------------ */ if (Rs == NULL) { /* no scaling */ switch (nr) { case 1: for (k = 0 ; k < n ; k++) { X [k] = Bz [Pnum [k]] ; } break ; case 2: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; X [2*k ] = Bz [i ] ; X [2*k + 1] = Bz [i + d ] ; } break ; case 3: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; X [3*k ] = Bz [i ] ; X [3*k + 1] = Bz [i + d ] ; X [3*k + 2] = Bz [i + d*2] ; } break ; case 4: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; X [4*k ] = Bz [i ] ; X [4*k + 1] = Bz [i + d ] ; X [4*k + 2] = Bz [i + d*2] ; X [4*k + 3] = Bz [i + d*3] ; } break ; } } else { switch (nr) { case 1: for (k = 0 ; k < n ; k++) { SCALE_DIV_ASSIGN (X [k], Bz [Pnum [k]], Rs [k]) ; } break ; case 2: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; rs = Rs [k] ; SCALE_DIV_ASSIGN (X [2*k], Bz [i], rs) ; SCALE_DIV_ASSIGN (X [2*k + 1], Bz [i + d], rs) ; } break ; case 3: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; rs = Rs [k] ; SCALE_DIV_ASSIGN (X [3*k], Bz [i], rs) ; SCALE_DIV_ASSIGN (X [3*k + 1], Bz [i + d], rs) ; SCALE_DIV_ASSIGN (X [3*k + 2], Bz [i + d*2], rs) ; } break ; case 4: for (k = 0 ; k < n ; k++) { i = Pnum [k] ; rs = Rs [k] ; SCALE_DIV_ASSIGN (X [4*k], Bz [i], rs) ; SCALE_DIV_ASSIGN (X [4*k + 1], Bz [i + d], rs) ; SCALE_DIV_ASSIGN (X [4*k + 2], Bz [i + d*2], rs) ; SCALE_DIV_ASSIGN (X [4*k + 3], Bz [i + d*3], rs) ; } break ; } } /* ------------------------------------------------------------------ */ /* solve X = (L*U + Off)\X */ /* ------------------------------------------------------------------ */ for (block = nblocks-1 ; block >= 0 ; block--) { /* -------------------------------------------------------------- */ /* the block of size nk is from rows/columns k1 to k2-1 */ /* -------------------------------------------------------------- */ k1 = R [block] ; k2 = R [block+1] ; nk = k2 - k1 ; PRINTF (("solve %d, k1 %d k2-1 %d nk %d\n", block, k1,k2-1,nk)) ; /* solve the block system */ if (nk == 1) { s = Udiag [k1] ; switch (nr) { case 1: DIV (X [k1], X [k1], s) ; break ; case 2: DIV (X [2*k1], X [2*k1], s) ; DIV (X [2*k1 + 1], X [2*k1 + 1], s) ; break ; case 3: DIV (X [3*k1], X [3*k1], s) ; DIV (X [3*k1 + 1], X [3*k1 + 1], s) ; DIV (X [3*k1 + 2], X [3*k1 + 2], s) ; break ; case 4: DIV (X [4*k1], X [4*k1], s) ; DIV (X [4*k1 + 1], X [4*k1 + 1], s) ; DIV (X [4*k1 + 2], X [4*k1 + 2], s) ; DIV (X [4*k1 + 3], X [4*k1 + 3], s) ; break ; } } else { KLU_lsolve (nk, Lip + k1, Llen + k1, LUbx [block], nr, X + nr*k1) ; KLU_usolve (nk, Uip + k1, Ulen + k1, LUbx [block], Udiag + k1, nr, X + nr*k1) ; } /* -------------------------------------------------------------- */ /* block back-substitution for the off-diagonal-block entries */ /* -------------------------------------------------------------- */ if (block > 0) { switch (nr) { case 1: for (k = k1 ; k < k2 ; k++) { pend = Offp [k+1] ; x [0] = X [k] ; for (p = Offp [k] ; p < pend ; p++) { MULT_SUB (X [Offi [p]], Offx [p], x [0]) ; } } break ; case 2: for (k = k1 ; k < k2 ; k++) { pend = Offp [k+1] ; x [0] = X [2*k ] ; x [1] = X [2*k + 1] ; for (p = Offp [k] ; p < pend ; p++) { i = Offi [p] ; offik = Offx [p] ; MULT_SUB (X [2*i], offik, x [0]) ; MULT_SUB (X [2*i + 1], offik, x [1]) ; } } break ; case 3: for (k = k1 ; k < k2 ; k++) { pend = Offp [k+1] ; x [0] = X [3*k ] ; x [1] = X [3*k + 1] ; x [2] = X [3*k + 2] ; for (p = Offp [k] ; p < pend ; p++) { i = Offi [p] ; offik = Offx [p] ; MULT_SUB (X [3*i], offik, x [0]) ; MULT_SUB (X [3*i + 1], offik, x [1]) ; MULT_SUB (X [3*i + 2], offik, x [2]) ; } } break ; case 4: for (k = k1 ; k < k2 ; k++) { pend = Offp [k+1] ; x [0] = X [4*k ] ; x [1] = X [4*k + 1] ; x [2] = X [4*k + 2] ; x [3] = X [4*k + 3] ; for (p = Offp [k] ; p < pend ; p++) { i = Offi [p] ; offik = Offx [p] ; MULT_SUB (X [4*i], offik, x [0]) ; MULT_SUB (X [4*i + 1], offik, x [1]) ; MULT_SUB (X [4*i + 2], offik, x [2]) ; MULT_SUB (X [4*i + 3], offik, x [3]) ; } } break ; } } } /* ------------------------------------------------------------------ */ /* permute the result, Bz = Q*X */ /* ------------------------------------------------------------------ */ switch (nr) { case 1: for (k = 0 ; k < n ; k++) { Bz [Q [k]] = X [k] ; } break ; case 2: for (k = 0 ; k < n ; k++) { i = Q [k] ; Bz [i ] = X [2*k ] ; Bz [i + d ] = X [2*k + 1] ; } break ; case 3: for (k = 0 ; k < n ; k++) { i = Q [k] ; Bz [i ] = X [3*k ] ; Bz [i + d ] = X [3*k + 1] ; Bz [i + d*2] = X [3*k + 2] ; } break ; case 4: for (k = 0 ; k < n ; k++) { i = Q [k] ; Bz [i ] = X [4*k ] ; Bz [i + d ] = X [4*k + 1] ; Bz [i + d*2] = X [4*k + 2] ; Bz [i + d*3] = X [4*k + 3] ; } break ; } /* ------------------------------------------------------------------ */ /* go to the next chunk of B */ /* ------------------------------------------------------------------ */ Bz += d*4 ; } return (TRUE) ; }
Int KLU_refactor /* returns TRUE if successful, FALSE otherwise */ ( /* inputs, not modified */ Int Ap [ ], /* size n+1, column pointers */ Int Ai [ ], /* size nz, row indices */ double Ax [ ], KLU_symbolic<Entry, Int> *Symbolic, /* input/output */ KLU_numeric<Entry, Int> *Numeric, KLU_common<Entry, Int> *Common ) { Entry ukk, ujk, s ; Entry *Offx, *Lx, *Ux, *X, *Az, *Udiag ; double *Rs ; Int *P, *Q, *R, *Pnum, *Offp, *Offi, *Ui, *Li, *Pinv, *Lip, *Uip, *Llen, *Ulen ; Unit **LUbx ; Unit *LU ; Int k1, k2, nk, k, block, oldcol, pend, oldrow, n, p, newrow, scale, nblocks, poff, i, j, up, ulen, llen, maxblock, nzoff ; /* ---------------------------------------------------------------------- */ /* check inputs */ /* ---------------------------------------------------------------------- */ if (Common == NULL) { return (FALSE) ; } Common->status = KLU_OK ; if (Numeric == NULL) { /* invalid Numeric object */ Common->status = KLU_INVALID ; return (FALSE) ; } Common->numerical_rank = EMPTY ; Common->singular_col = EMPTY ; Az = (Entry *) Ax ; /* ---------------------------------------------------------------------- */ /* get the contents of the Symbolic object */ /* ---------------------------------------------------------------------- */ n = Symbolic->n ; P = Symbolic->P ; Q = Symbolic->Q ; R = Symbolic->R ; nblocks = Symbolic->nblocks ; maxblock = Symbolic->maxblock ; /* ---------------------------------------------------------------------- */ /* get the contents of the Numeric object */ /* ---------------------------------------------------------------------- */ Pnum = Numeric->Pnum ; Offp = Numeric->Offp ; Offi = Numeric->Offi ; Offx = (Entry *) Numeric->Offx ; LUbx = (Unit **) Numeric->LUbx ; scale = Common->scale ; if (scale > 0) { /* factorization was not scaled, but refactorization is scaled */ if (Numeric->Rs == NULL) { Numeric->Rs = (double *)KLU_malloc (n, sizeof (double), Common) ; if (Common->status < KLU_OK) { Common->status = KLU_OUT_OF_MEMORY ; return (FALSE) ; } } } else { /* no scaling for refactorization; ensure Numeric->Rs is freed. This * does nothing if Numeric->Rs is already NULL. */ Numeric->Rs = (double *) KLU_free (Numeric->Rs, n, sizeof (double), Common) ; } Rs = Numeric->Rs ; Pinv = Numeric->Pinv ; X = (Entry *) Numeric->Xwork ; Common->nrealloc = 0 ; Udiag = (Entry *) Numeric->Udiag ; nzoff = Symbolic->nzoff ; /* ---------------------------------------------------------------------- */ /* check the input matrix compute the row scale factors, Rs */ /* ---------------------------------------------------------------------- */ /* do no scale, or check the input matrix, if scale < 0 */ if (scale >= 0) { /* check for out-of-range indices, but do not check for duplicates */ if (!KLU_scale (scale, n, Ap, Ai, Ax, Rs, NULL, Common)) { return (FALSE) ; } } /* ---------------------------------------------------------------------- */ /* clear workspace X */ /* ---------------------------------------------------------------------- */ for (k = 0 ; k < maxblock ; k++) { /* X [k] = 0 */ CLEAR (X [k]) ; } poff = 0 ; /* ---------------------------------------------------------------------- */ /* factor each block */ /* ---------------------------------------------------------------------- */ if (scale <= 0) { /* ------------------------------------------------------------------ */ /* no scaling */ /* ------------------------------------------------------------------ */ for (block = 0 ; block < nblocks ; block++) { /* -------------------------------------------------------------- */ /* the block is from rows/columns k1 to k2-1 */ /* -------------------------------------------------------------- */ k1 = R [block] ; k2 = R [block+1] ; nk = k2 - k1 ; if (nk == 1) { /* ---------------------------------------------------------- */ /* singleton case */ /* ---------------------------------------------------------- */ oldcol = Q [k1] ; pend = Ap [oldcol+1] ; CLEAR (s) ; for (p = Ap [oldcol] ; p < pend ; p++) { newrow = Pinv [Ai [p]] - k1 ; if (newrow < 0 && poff < nzoff) { /* entry in off-diagonal block */ Offx [poff] = Az [p] ; poff++ ; } else { /* singleton */ s = Az [p] ; } } Udiag [k1] = s ; } else { /* ---------------------------------------------------------- */ /* construct and factor the kth block */ /* ---------------------------------------------------------- */ Lip = Numeric->Lip + k1 ; Llen = Numeric->Llen + k1 ; Uip = Numeric->Uip + k1 ; Ulen = Numeric->Ulen + k1 ; LU = LUbx [block] ; for (k = 0 ; k < nk ; k++) { /* ------------------------------------------------------ */ /* scatter kth column of the block into workspace X */ /* ------------------------------------------------------ */ oldcol = Q [k+k1] ; pend = Ap [oldcol+1] ; for (p = Ap [oldcol] ; p < pend ; p++) { newrow = Pinv [Ai [p]] - k1 ; if (newrow < 0 && poff < nzoff) { /* entry in off-diagonal block */ Offx [poff] = Az [p] ; poff++ ; } else { /* (newrow,k) is an entry in the block */ X [newrow] = Az [p] ; } } /* ------------------------------------------------------ */ /* compute kth column of U, and update kth column of A */ /* ------------------------------------------------------ */ GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, ulen) ; for (up = 0 ; up < ulen ; up++) { j = Ui [up] ; ujk = X [j] ; /* X [j] = 0 */ CLEAR (X [j]) ; Ux [up] = ujk ; GET_POINTER (LU, Lip, Llen, Li, Lx, j, llen) ; for (p = 0 ; p < llen ; p++) { /* X [Li [p]] -= Lx [p] * ujk */ MULT_SUB (X [Li [p]], Lx [p], ujk) ; } } /* get the diagonal entry of U */ ukk = X [k] ; /* X [k] = 0 */ CLEAR (X [k]) ; /* singular case */ if (IS_ZERO (ukk)) { /* matrix is numerically singular */ Common->status = KLU_SINGULAR ; if (Common->numerical_rank == EMPTY) { Common->numerical_rank = k+k1 ; Common->singular_col = Q [k+k1] ; } if (Common->halt_if_singular) { /* do not continue the factorization */ return (FALSE) ; } } Udiag [k+k1] = ukk ; /* gather and divide by pivot to get kth column of L */ GET_POINTER (LU, Lip, Llen, Li, Lx, k, llen) ; for (p = 0 ; p < llen ; p++) { i = Li [p] ; DIV (Lx [p], X [i], ukk) ; CLEAR (X [i]) ; } } } } } else { /* ------------------------------------------------------------------ */ /* scaling */ /* ------------------------------------------------------------------ */ for (block = 0 ; block < nblocks ; block++) { /* -------------------------------------------------------------- */ /* the block is from rows/columns k1 to k2-1 */ /* -------------------------------------------------------------- */ k1 = R [block] ; k2 = R [block+1] ; nk = k2 - k1 ; if (nk == 1) { /* ---------------------------------------------------------- */ /* singleton case */ /* ---------------------------------------------------------- */ oldcol = Q [k1] ; pend = Ap [oldcol+1] ; CLEAR (s) ; for (p = Ap [oldcol] ; p < pend ; p++) { oldrow = Ai [p] ; newrow = Pinv [oldrow] - k1 ; if (newrow < 0 && poff < nzoff) { /* entry in off-diagonal block */ /* Offx [poff] = Az [p] / Rs [oldrow] */ SCALE_DIV_ASSIGN (Offx [poff], Az [p], Rs [oldrow]) ; poff++ ; } else { /* singleton */ /* s = Az [p] / Rs [oldrow] */ SCALE_DIV_ASSIGN (s, Az [p], Rs [oldrow]) ; } } Udiag [k1] = s ; } else { /* ---------------------------------------------------------- */ /* construct and factor the kth block */ /* ---------------------------------------------------------- */ Lip = Numeric->Lip + k1 ; Llen = Numeric->Llen + k1 ; Uip = Numeric->Uip + k1 ; Ulen = Numeric->Ulen + k1 ; LU = LUbx [block] ; for (k = 0 ; k < nk ; k++) { /* ------------------------------------------------------ */ /* scatter kth column of the block into workspace X */ /* ------------------------------------------------------ */ oldcol = Q [k+k1] ; pend = Ap [oldcol+1] ; for (p = Ap [oldcol] ; p < pend ; p++) { oldrow = Ai [p] ; newrow = Pinv [oldrow] - k1 ; if (newrow < 0 && poff < nzoff) { /* entry in off-diagonal part */ /* Offx [poff] = Az [p] / Rs [oldrow] */ SCALE_DIV_ASSIGN (Offx [poff], Az [p], Rs [oldrow]); poff++ ; } else { /* (newrow,k) is an entry in the block */ /* X [newrow] = Az [p] / Rs [oldrow] */ SCALE_DIV_ASSIGN (X [newrow], Az [p], Rs [oldrow]) ; } } /* ------------------------------------------------------ */ /* compute kth column of U, and update kth column of A */ /* ------------------------------------------------------ */ GET_POINTER (LU, Uip, Ulen, Ui, Ux, k, ulen) ; for (up = 0 ; up < ulen ; up++) { j = Ui [up] ; ujk = X [j] ; /* X [j] = 0 */ CLEAR (X [j]) ; Ux [up] = ujk ; GET_POINTER (LU, Lip, Llen, Li, Lx, j, llen) ; for (p = 0 ; p < llen ; p++) { /* X [Li [p]] -= Lx [p] * ujk */ MULT_SUB (X [Li [p]], Lx [p], ujk) ; } } /* get the diagonal entry of U */ ukk = X [k] ; /* X [k] = 0 */ CLEAR (X [k]) ; /* singular case */ if (IS_ZERO (ukk)) { /* matrix is numerically singular */ Common->status = KLU_SINGULAR ; if (Common->numerical_rank == EMPTY) { Common->numerical_rank = k+k1 ; Common->singular_col = Q [k+k1] ; } if (Common->halt_if_singular) { /* do not continue the factorization */ return (FALSE) ; } } Udiag [k+k1] = ukk ; /* gather and divide by pivot to get kth column of L */ GET_POINTER (LU, Lip, Llen, Li, Lx, k, llen) ; for (p = 0 ; p < llen ; p++) { i = Li [p] ; DIV (Lx [p], X [i], ukk) ; CLEAR (X [i]) ; } } } } } /* ---------------------------------------------------------------------- */ /* permute scale factors Rs according to pivotal row order */ /* ---------------------------------------------------------------------- */ if (scale > 0) { for (k = 0 ; k < n ; k++) { /* TODO : Check. REAL(X[k]) Can be just X[k] */ /* REAL (X [k]) = Rs [Pnum [k]] ; */ X [k] = Rs [Pnum [k]] ; } for (k = 0 ; k < n ; k++) { Rs [k] = REAL (X [k]) ; } } #ifndef NDEBUGKLU2 ASSERT (Offp [n] == poff) ; ASSERT (Symbolic->nzoff == poff) ; PRINTF (("\n------------------- Off diagonal entries, new:\n")) ; ASSERT (KLU_valid (n, Offp, Offi, Offx)) ; if (Common->status == KLU_OK) { PRINTF (("\n ########### KLU_BTF_REFACTOR done, nblocks %d\n",nblocks)); for (block = 0 ; block < nblocks ; block++) { k1 = R [block] ; k2 = R [block+1] ; nk = k2 - k1 ; PRINTF (( "\n================KLU_refactor output: k1 %d k2 %d nk %d\n", k1, k2, nk)) ; if (nk == 1) { PRINTF (("singleton ")) ; PRINT_ENTRY (Udiag [k1]) ; } else { Lip = Numeric->Lip + k1 ; Llen = Numeric->Llen + k1 ; LU = (Unit *) Numeric->LUbx [block] ; PRINTF (("\n---- L block %d\n", block)) ; ASSERT (KLU_valid_LU (nk, TRUE, Lip, Llen, LU)) ; Uip = Numeric->Uip + k1 ; Ulen = Numeric->Ulen + k1 ; PRINTF (("\n---- U block %d\n", block)) ; ASSERT (KLU_valid_LU (nk, FALSE, Uip, Ulen, LU)) ; } } } #endif return (TRUE) ; }