csn *cs_chol(const cs *A, const css *S) { double d, lki; double *Lx, *x, *Cx; int top, i, p, k, n, *Li, *Lp, *cp, *pinv, *s, *c, *parent, *Cp, *Ci; cs *L, *C, *E; csn *N; if (!CS_CSC (A) || !S || !S->cp || !S->parent) return (NULL); n = A->n; N = (csn *) cs_calloc(1, sizeof(csn)); /* allocate result */ c = (int *) cs_malloc(2 * n, sizeof(int)); /* get int workspace */ x = (double *) cs_malloc(n, sizeof(double)); /* get double workspace */ cp = S->cp; pinv = S->pinv; parent = S->parent; C = pinv ? cs_symperm(A, pinv, 1) : ((cs *) A); E = pinv ? C : NULL; /* E is alias for A, or a copy E=A(p,p) */ if (!N || !c || !x || !C) return (cs_ndone(N, E, c, x, 0)); s = c + n; Cp = C->p; Ci = C->i; Cx = C->x; N->L = L = cs_spalloc(n, n, cp[n], 1, 0); /* allocate result */ if (!L) return (cs_ndone(N, E, c, x, 0)); Lp = L->p; Li = L->i; Lx = L->x; for (k = 0; k < n; k++) Lp[k] = c[k] = cp[k]; for (k = 0; k < n; k++) /* compute L(k,:) for L*L' = C */ { /* --- Nonzero pattern of L(k,:) ------------------------------------ */ top = cs_ereach(C, k, parent, s, c); /* find pattern of L(k,:) */ x[k] = 0; /* x (0:k) is now zero */ for (p = Cp[k]; p < Cp[k + 1]; p++) /* x = full(triu(C(:,k))) */ { if (Ci[p] <= k) x[Ci[p]] = Cx[p]; } d = x[k]; /* d = C(k,k) */ x[k] = 0; /* clear x for k+1st iteration */ /* --- Triangular solve --------------------------------------------- */ for (; top < n; top++) /* solve L(0:k-1,0:k-1) * x = C(:,k) */ { i = s[top]; /* s [top..n-1] is pattern of L(k,:) */ lki = x[i] / Lx[Lp[i]]; /* L(k,i) = x (i) / L(i,i) */ x[i] = 0; /* clear x for k+1st iteration */ for (p = Lp[i] + 1; p < c[i]; p++) { x[Li[p]] -= Lx[p] * lki; } d -= lki * lki; /* d = d - L(k,i)*L(k,i) */ p = c[i]++; Li[p] = k; /* store L(k,i) in column i */ Lx[p] = lki; } /* --- Compute L(k,k) ----------------------------------------------- */ if (d <= 0) return (cs_ndone(N, E, c, x, 0)); /* not pos def */ p = c[k]++; Li[p] = k; /* store L(k,k) = sqrt (d) in column k */ Lx[p] = sqrt(d); } Lp[n] = cp[n]; /* finalize L */ return (cs_ndone(N, E, c, x, 1)); /* success: free E,s,x; return N */ }
/* [L,U,pinv]=lu(A, [q lnz unz]). lnz and unz can be guess */ csn *cs_lu (const cs *A, const css *S, double tol) { cs *L, *U ; csn *N ; double pivot, *Lx, *Ux, *x, a, t ; int *Lp, *Li, *Up, *Ui, *pinv, *xi, *q, n, ipiv, k, top, p, i, col, lnz,unz; if (!CS_CSC (A) || !S) return (NULL) ; /* check inputs */ n = A->n ; q = S->q ; lnz = S->lnz ; unz = S->unz ; x = cs_malloc (n, sizeof (double)) ; /* get double workspace */ xi = cs_malloc (2*n, sizeof (int)) ; /* get int workspace */ N = cs_calloc (1, sizeof (csn)) ; /* allocate result */ if (!x || !xi || !N) return (cs_ndone (N, NULL, xi, x, 0)) ; N->L = L = cs_spalloc (n, n, lnz, 1, 0) ; /* allocate result L */ N->U = U = cs_spalloc (n, n, unz, 1, 0) ; /* allocate result U */ N->pinv = pinv = cs_malloc (n, sizeof (int)) ; /* allocate result pinv */ if (!L || !U || !pinv) return (cs_ndone (N, NULL, xi, x, 0)) ; Lp = L->p ; Up = U->p ; for (i = 0 ; i < n ; i++) x [i] = 0 ; /* clear workspace */ for (i = 0 ; i < n ; i++) pinv [i] = -1 ; /* no rows pivotal yet */ for (k = 0 ; k <= n ; k++) Lp [k] = 0 ; /* no cols of L yet */ lnz = unz = 0 ; for (k = 0 ; k < n ; k++) /* compute L(:,k) and U(:,k) */ { /* --- Triangular solve --------------------------------------------- */ Lp [k] = lnz ; /* L(:,k) starts here */ Up [k] = unz ; /* U(:,k) starts here */ if ((lnz + n > L->nzmax && !cs_sprealloc (L, 2*L->nzmax + n)) || (unz + n > U->nzmax && !cs_sprealloc (U, 2*U->nzmax + n))) { return (cs_ndone (N, NULL, xi, x, 0)) ; } Li = L->i ; Lx = L->x ; Ui = U->i ; Ux = U->x ; col = q ? (q [k]) : k ; top = cs_spsolve (L, A, col, xi, x, pinv, 1) ; /* x = L\A(:,col) */ /* --- Find pivot --------------------------------------------------- */ ipiv = -1 ; a = -1 ; for (p = top ; p < n ; p++) { i = xi [p] ; /* x(i) is nonzero */ if (pinv [i] < 0) /* row i is not yet pivotal */ { if ((t = fabs (x [i])) > a) { a = t ; /* largest pivot candidate so far */ ipiv = i ; } } else /* x(i) is the entry U(pinv[i],k) */ { Ui [unz] = pinv [i] ; Ux [unz++] = x [i] ; } } if (ipiv == -1 || a <= 0) return (cs_ndone (N, NULL, xi, x, 0)) ; if (pinv [col] < 0 && fabs (x [col]) >= a*tol) ipiv = col ; /* --- Divide by pivot ---------------------------------------------- */ pivot = x [ipiv] ; /* the chosen pivot */ Ui [unz] = k ; /* last entry in U(:,k) is U(k,k) */ Ux [unz++] = pivot ; pinv [ipiv] = k ; /* ipiv is the kth pivot row */ Li [lnz] = ipiv ; /* first entry in L(:,k) is L(k,k) = 1 */ Lx [lnz++] = 1 ; for (p = top ; p < n ; p++) /* L(k+1:n,k) = x / pivot */ { i = xi [p] ; if (pinv [i] < 0) /* x(i) is an entry in L(:,k) */ { Li [lnz] = i ; /* save unpermuted row in L */ Lx [lnz++] = x [i] / pivot ; /* scale pivot column */ } x [i] = 0 ; /* x [0..n-1] = 0 for next k */ } } /* --- Finalize L and U ------------------------------------------------- */ Lp [n] = lnz ; Up [n] = unz ; Li = L->i ; /* fix row indices of L for final pinv */ for (p = 0 ; p < lnz ; p++) Li [p] = pinv [Li [p]] ; cs_sprealloc (L, 0) ; /* remove extra space from L and U */ cs_sprealloc (U, 0) ; return (cs_ndone (N, NULL, xi, x, 1)) ; /* success */ }
/* sparse QR factorization [V,beta,pinv,R] = qr (A) */ csn *cs_qr (const cs *A, const css *S) { CS_ENTRY *Rx, *Vx, *Ax, *x ; double *Beta ; CS_INT i, k, p, m, n, vnz, p1, top, m2, len, col, rnz, *s, *leftmost, *Ap, *Ai, *parent, *Rp, *Ri, *Vp, *Vi, *w, *pinv, *q ; cs *R, *V ; csn *N ; if (!CS_CSC (A) || !S) return (NULL) ; m = A->m ; n = A->n ; Ap = A->p ; Ai = A->i ; Ax = A->x ; q = S->q ; parent = S->parent ; pinv = S->pinv ; m2 = S->m2 ; vnz = S->lnz ; rnz = S->unz ; leftmost = S->leftmost ; w = cs_malloc (m2+n, sizeof (CS_INT)) ; /* get CS_INT workspace */ x = cs_malloc (m2, sizeof (CS_ENTRY)) ; /* get CS_ENTRY workspace */ N = cs_calloc (1, sizeof (csn)) ; /* allocate result */ if (!w || !x || !N) return (cs_ndone (N, NULL, w, x, 0)) ; s = w + m2 ; /* s is size n */ for (k = 0 ; k < m2 ; k++) x [k] = 0 ; /* clear workspace x */ N->L = V = cs_spalloc (m2, n, vnz, 1, 0) ; /* allocate result V */ N->U = R = cs_spalloc (m2, n, rnz, 1, 0) ; /* allocate result R */ N->B = Beta = cs_malloc (n, sizeof (double)) ; /* allocate result Beta */ if (!R || !V || !Beta) return (cs_ndone (N, NULL, w, x, 0)) ; Rp = R->p ; Ri = R->i ; Rx = R->x ; Vp = V->p ; Vi = V->i ; Vx = V->x ; for (i = 0 ; i < m2 ; i++) w [i] = -1 ; /* clear w, to mark nodes */ rnz = 0 ; vnz = 0 ; for (k = 0 ; k < n ; k++) /* compute V and R */ { Rp [k] = rnz ; /* R(:,k) starts here */ Vp [k] = p1 = vnz ; /* V(:,k) starts here */ w [k] = k ; /* add V(k,k) to pattern of V */ Vi [vnz++] = k ; top = n ; col = q ? q [k] : k ; for (p = Ap [col] ; p < Ap [col+1] ; p++) /* find R(:,k) pattern */ { i = leftmost [Ai [p]] ; /* i = min(find(A(i,q))) */ for (len = 0 ; w [i] != k ; i = parent [i]) /* traverse up to k */ { s [len++] = i ; w [i] = k ; } while (len > 0) s [--top] = s [--len] ; /* push path on stack */ i = pinv [Ai [p]] ; /* i = permuted row of A(:,col) */ x [i] = Ax [p] ; /* x (i) = A(:,col) */ if (i > k && w [i] < k) /* pattern of V(:,k) = x (k+1:m) */ { Vi [vnz++] = i ; /* add i to pattern of V(:,k) */ w [i] = k ; } } for (p = top ; p < n ; p++) /* for each i in pattern of R(:,k) */ { i = s [p] ; /* R(i,k) is nonzero */ cs_happly (V, i, Beta [i], x) ; /* apply (V(i),Beta(i)) to x */ Ri [rnz] = i ; /* R(i,k) = x(i) */ Rx [rnz++] = x [i] ; x [i] = 0 ; if (parent [i] == k) vnz = cs_scatter (V, i, 0, w, NULL, k, V, vnz); } for (p = p1 ; p < vnz ; p++) /* gather V(:,k) = x */ { Vx [p] = x [Vi [p]] ; x [Vi [p]] = 0 ; } Ri [rnz] = k ; /* R(k,k) = norm (x) */ Rx [rnz++] = cs_house (Vx+p1, Beta+k, vnz-p1) ; /* [v,beta]=house(x) */ } Rp [n] = rnz ; /* finalize R */ Vp [n] = vnz ; /* finalize V */ return (cs_ndone (N, NULL, w, x, 1)) ; /* success */ }