// called from package MatrixModels's R code: SEXP dgCMatrix_cholsol(SEXP x, SEXP y) { /* Solve Sparse Least Squares X %*% beta ~= y with dense RHS y, * where X = t(x) i.e. we pass x = t(X) as argument, * via "Cholesky(X'X)" .. well not really: * cholmod_factorize("x", ..) finds L in X'X = L'L directly */ CHM_SP cx = AS_CHM_SP(x); /* FIXME: extend this to work in multivariate case, i.e. y a matrix with > 1 column ! */ CHM_DN cy = AS_CHM_DN(coerceVector(y, REALSXP)), rhs, cAns, resid; CHM_FR L; int n = cx->ncol;/* #{obs.} {x = t(X) !} */ double one[] = {1,0}, zero[] = {0,0}, neg1[] = {-1,0}; const char *nms[] = {"L", "coef", "Xty", "resid", ""}; SEXP ans = PROTECT(Rf_mkNamed(VECSXP, nms)); R_CheckStack(); if (n < cx->nrow || n <= 0) error(_("dgCMatrix_cholsol requires a 'short, wide' rectangular matrix")); if (cy->nrow != n) error(_("Dimensions of system to be solved are inconsistent")); rhs = cholmod_allocate_dense(cx->nrow, 1, cx->nrow, CHOLMOD_REAL, &c); /* cholmod_sdmult(A, transp, alpha, beta, X, Y, &c): * Y := alpha*(A*X) + beta*Y or alpha*(A'*X) + beta*Y ; * here: rhs := 1 * x %*% y + 0 = x %*% y = X'y */ if (!(cholmod_sdmult(cx, 0 /* trans */, one, zero, cy, rhs, &c))) error(_("cholmod_sdmult error (rhs)")); L = cholmod_analyze(cx, &c); if (!cholmod_factorize(cx, L, &c)) error(_("cholmod_factorize failed: status %d, minor %d from ncol %d"), c.status, L->minor, L->n); /* FIXME: Do this in stages so an "effects" vector can be calculated */ if (!(cAns = cholmod_solve(CHOLMOD_A, L, rhs, &c))) error(_("cholmod_solve (CHOLMOD_A) failed: status %d, minor %d from ncol %d"), c.status, L->minor, L->n); /* L : */ SET_VECTOR_ELT(ans, 0, chm_factor_to_SEXP(L, 0)); /* coef : */ SET_VECTOR_ELT(ans, 1, allocVector(REALSXP, cx->nrow)); Memcpy(REAL(VECTOR_ELT(ans, 1)), (double*)(cAns->x), cx->nrow); /* X'y : */ /* FIXME: Change this when the "effects" vector is available */ SET_VECTOR_ELT(ans, 2, allocVector(REALSXP, cx->nrow)); Memcpy(REAL(VECTOR_ELT(ans, 2)), (double*)(rhs->x), cx->nrow); /* resid := y */ resid = cholmod_copy_dense(cy, &c); /* cholmod_sdmult(A, transp, alp, bet, X, Y, &c): * Y := alp*(A*X) + bet*Y or alp*(A'*X) + beta*Y ; * here: resid := -1 * x' %*% coef + 1 * y = y - X %*% coef */ if (!(cholmod_sdmult(cx, 1/* trans */, neg1, one, cAns, resid, &c))) error(_("cholmod_sdmult error (resid)")); /* FIXME: for multivariate case, i.e. resid *matrix* with > 1 column ! */ SET_VECTOR_ELT(ans, 3, allocVector(REALSXP, n)); Memcpy(REAL(VECTOR_ELT(ans, 3)), (double*)(resid->x), n); cholmod_free_factor(&L, &c); cholmod_free_dense(&rhs, &c); cholmod_free_dense(&cAns, &c); UNPROTECT(1); return ans; }
SEXP destructive_CHM_update(SEXP object, SEXP parent, SEXP mult) { CHM_FR L = AS_CHM_FR(object); CHM_SP A = AS_CHM_SP__(parent); R_CheckStack(); return chm_factor_to_SEXP(chm_factor_update(L, A, asReal(mult)), 0); }
SEXP CHMfactor_update(SEXP object, SEXP parent, SEXP mult) { CHM_FR L = AS_CHM_FR(object), Lcp; CHM_SP A = AS_CHM_SP__(parent); R_CheckStack(); Lcp = cholmod_copy_factor(L, &c); return chm_factor_to_SEXP(chm_factor_update(Lcp, A, asReal(mult)), 1); }
/** * Return a CHOLMOD copy of the cached Cholesky decomposition with the * required perm, LDL and super attributes. If Imult is nonzero, * update the numeric values before returning. * * If no cached copy is available then evaluate one, cache it (for * zero Imult), and return a copy. * * @param Ap dsCMatrix object * @param perm integer indicating if permutation is required (>0), * forbidden (0) or optional (<0) * @param LDL integer indicating if the LDL' form is required (>0), * forbidden (0) or optional (<0) * @param super integer indicating if the supernodal form is required (>0), * forbidden (0) or optional (<0) * @param Imult numeric multiplier of I in |A + Imult * I| */ static CHM_FR internal_chm_factor(SEXP Ap, int perm, int LDL, int super, double Imult) { SEXP facs = GET_SLOT(Ap, Matrix_factorSym); SEXP nms = getAttrib(facs, R_NamesSymbol); int sup, ll; CHM_FR L; CHM_SP A = AS_CHM_SP__(Ap); R_CheckStack(); if (LENGTH(facs)) { for (int i = 0; i < LENGTH(nms); i++) { /* look for a match in cache */ if (chk_nm(CHAR(STRING_ELT(nms, i)), perm, LDL, super)) { L = AS_CHM_FR(VECTOR_ELT(facs, i)); R_CheckStack(); /* copy the factor so later it can safely be cholmod_l_free'd */ L = cholmod_l_copy_factor(L, &c); if (Imult) cholmod_l_factorize_p(A, &Imult, (int*)NULL, 0, L, &c); return L; } } } /* No cached factor - create one */ sup = c.supernodal; /* save current settings */ ll = c.final_ll; c.final_ll = (LDL == 0) ? 1 : 0; c.supernodal = (super > 0) ? CHOLMOD_SUPERNODAL : CHOLMOD_SIMPLICIAL; if (perm) { /* obtain fill-reducing permutation */ L = cholmod_l_analyze(A, &c); } else { /* require identity permutation */ /* save current settings */ int nmethods = c.nmethods, ord0 = c.method[0].ordering, postorder = c.postorder; c.nmethods = 1; c.method[0].ordering = CHOLMOD_NATURAL; c.postorder = FALSE; L = cholmod_l_analyze(A, &c); /* and now restore */ c.nmethods = nmethods; c.method[0].ordering = ord0; c.postorder = postorder; } if (!cholmod_l_factorize_p(A, &Imult, (int*)NULL, 0 /*fsize*/, L, &c)) error(_("Cholesky factorization failed")); c.supernodal = sup; /* restore previous settings */ c.final_ll = ll; if (!Imult) { /* cache the factor */ char fnm[12] = "sPDCholesky"; if (super > 0) fnm[0] = 'S'; if (perm == 0) fnm[1] = 'p'; if (LDL == 0) fnm[2] = 'd'; set_factors(Ap, chm_factor_to_SEXP(L, 0), fnm); } return L; }
SEXP CHMfactor_updown(SEXP upd, SEXP C_, SEXP L_) { CHM_FR L = AS_CHM_FR(L_), Lcp; CHM_SP C = AS_CHM_SP__(C_); int update = asInteger(upd); R_CheckStack(); Lcp = cholmod_copy_factor(L, &c); int r = cholmod_updown(update, C, Lcp, &c); if(!r) error(_("cholmod_updown() returned %d"), r); return chm_factor_to_SEXP(Lcp, 1); }
SEXP dsCMatrix_Cholesky(SEXP Ap, SEXP perm, SEXP LDL, SEXP super, SEXP Imult) { int c_pr = c.print; c.print = 0;/* stop CHOLMOD printing; we cannot suppress it (in R), and have error handler already */ SEXP r = chm_factor_to_SEXP(internal_chm_factor(Ap, asLogical(perm), asLogical(LDL), asLogical(super), asReal(Imult)), 1 /* dofree */); c.print = c_pr; return r; }
SEXP dsCMatrix_Cholesky(SEXP Ap, SEXP permP, SEXP LDLp, SEXP superP) { char *fname = strdup("spdCholesky"); /* template for factorization name */ int perm = asLogical(permP), LDL = asLogical(LDLp), super = asLogical(superP); SEXP Chol; cholmod_sparse *A; cholmod_factor *L; int sup, ll; if (super) fname[0] = 'S'; if (perm) fname[1] = 'P'; if (LDL) fname[2] = 'D'; Chol = get_factors(Ap, "fname"); if (Chol != R_NilValue) return Chol; A = as_cholmod_sparse(Ap); sup = c.supernodal; ll = c.final_ll; if (!A->stype) error("Non-symmetric matrix passed to dsCMatrix_chol"); c.final_ll = !LDL; /* leave as LL' or form LDL' */ c.supernodal = super ? CHOLMOD_SUPERNODAL : CHOLMOD_SIMPLICIAL; if (perm) { L = cholmod_analyze(A, &c); /* get fill-reducing permutation */ } else { /* require identity permutation */ int nmethods = c.nmethods, ord0 = c.method[0].ordering, postorder = c.postorder; c.nmethods = 1; c.method[0].ordering = CHOLMOD_NATURAL; c.postorder = FALSE; L = cholmod_analyze(A, &c); c.nmethods = nmethods; c.method[0].ordering = ord0; c.postorder = postorder; } c.supernodal = sup; /* restore previous setting */ c.final_ll = ll; if (!cholmod_factorize(A, L, &c)) error(_("Cholesky factorization failed")); Free(A); Chol = set_factors(Ap, chm_factor_to_SEXP(L, 1), fname); free(fname); /* note, this must be free, not Free */ return Chol; }