|
Matrix $Rev: 2718 $ at $LastChangedDate: 2011-10-06 11:45:17 +0200 (Thu, 06 Oct 2011) $
|
00001 #include "dgCMatrix.h" 00002 00003 /* for Csparse_transpose() : */ 00004 #include "Csparse.h" 00005 #include "chm_common.h" 00006 /* -> Mutils.h / SPQR ... */ 00007 00008 /* FIXME -- we "forget" about dimnames almost everywhere : */ 00009 00010 /* for dgCMatrix _and_ lgCMatrix and others (but *not* ngC...) : */ 00011 SEXP xCMatrix_validate(SEXP x) 00012 { 00013 /* Almost everything now in Csparse_validate ( ./Csparse.c ) 00014 * *but* the checking of the 'x' slot : */ 00015 if (length(GET_SLOT(x, Matrix_iSym)) != 00016 length(GET_SLOT(x, Matrix_xSym))) 00017 return mkString(_("lengths of slots 'i' and 'x' must match")); 00018 00019 return ScalarLogical(1); 00020 } 00021 00022 /* for dgRMatrix _and_ lgRMatrix and others (but *not* ngC...) : */ 00023 SEXP xRMatrix_validate(SEXP x) 00024 { 00025 /* Almost everything now in Rsparse_validate ( ./Csparse.c ) 00026 * *but* the checking of the 'x' slot : */ 00027 if (length(GET_SLOT(x, Matrix_jSym)) != 00028 length(GET_SLOT(x, Matrix_xSym))) 00029 return mkString(_("lengths of slots 'j' and 'x' must match")); 00030 00031 return ScalarLogical(1); 00032 } 00033 00034 /* This and the following R_to_CMatrix() lead to memory-not-mapped seg.faults 00035 * only with {32bit + R-devel + enable-R-shlib} -- no idea why */ 00036 SEXP compressed_to_TMatrix(SEXP x, SEXP colP) 00037 { 00038 int col = asLogical(colP); /* 1 if "C"olumn compressed; 0 if "R"ow */ 00039 /* however, for Csparse, we now effectively use the cholmod-based 00040 * Csparse_to_Tsparse() in ./Csparse.c ; maybe should simply write 00041 * an as_cholmod_Rsparse() function and then do "as there" ...*/ 00042 SEXP indSym = col ? Matrix_iSym : Matrix_jSym, 00043 ans, indP = GET_SLOT(x, indSym), 00044 pP = GET_SLOT(x, Matrix_pSym); 00045 int npt = length(pP) - 1; 00046 char *ncl = strdup(class_P(x)); 00047 static const char *valid[] = { MATRIX_VALID_Csparse, MATRIX_VALID_Rsparse, ""}; 00048 int ctype = Matrix_check_class_etc(x, valid); 00049 00050 if (ctype < 0) 00051 error(_("invalid class(x) '%s' in compressed_to_TMatrix(x)"), ncl); 00052 00053 /* replace 'C' or 'R' with 'T' :*/ 00054 ncl[2] = 'T'; 00055 ans = PROTECT(NEW_OBJECT(MAKE_CLASS(ncl))); 00056 00057 slot_dup(ans, x, Matrix_DimSym); 00058 if((ctype / 3) % 4 != 2) /* not n..Matrix */ 00059 slot_dup(ans, x, Matrix_xSym); 00060 if(ctype % 3) { /* s(ymmetric) or t(riangular) : */ 00061 slot_dup(ans, x, Matrix_uploSym); 00062 if(ctype % 3 == 2) /* t(riangular) : */ 00063 slot_dup(ans, x, Matrix_diagSym); 00064 } 00065 SET_DimNames(ans, x); 00066 SET_SLOT(ans, indSym, duplicate(indP)); 00067 expand_cmprPt(npt, INTEGER(pP), 00068 INTEGER(ALLOC_SLOT(ans, col ? Matrix_jSym : Matrix_iSym, 00069 INTSXP, length(indP)))); 00070 free(ncl); 00071 UNPROTECT(1); 00072 return ans; 00073 } 00074 00075 SEXP R_to_CMatrix(SEXP x) 00076 { 00077 SEXP ans, tri = PROTECT(allocVector(LGLSXP, 1)); 00078 char *ncl = strdup(class_P(x)); 00079 static const char *valid[] = { MATRIX_VALID_Rsparse, ""}; 00080 int ctype = Matrix_check_class_etc(x, valid); 00081 int *x_dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), *a_dims; 00082 PROTECT_INDEX ipx; 00083 00084 if (ctype < 0) 00085 error(_("invalid class(x) '%s' in R_to_CMatrix(x)"), ncl); 00086 00087 /* replace 'R' with 'C' : */ 00088 ncl[2] = 'C'; 00089 PROTECT_WITH_INDEX(ans = NEW_OBJECT(MAKE_CLASS(ncl)), &ipx); 00090 00091 a_dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); 00092 /* reversed dim() since we will transpose: */ 00093 a_dims[0] = x_dims[1]; 00094 a_dims[1] = x_dims[0]; 00095 00096 /* triangular: */ LOGICAL(tri)[0] = 0; 00097 if((ctype / 3) != 2) /* not n..Matrix */ 00098 slot_dup(ans, x, Matrix_xSym); 00099 if(ctype % 3) { /* s(ymmetric) or t(riangular) : */ 00100 SET_SLOT(ans, Matrix_uploSym, 00101 mkString((*uplo_P(x) == 'U') ? "L" : "U")); 00102 if(ctype % 3 == 2) { /* t(riangular) : */ 00103 LOGICAL(tri)[0] = 1; 00104 slot_dup(ans, x, Matrix_diagSym); 00105 } 00106 } 00107 SET_SLOT(ans, Matrix_iSym, duplicate(GET_SLOT(x, Matrix_jSym))); 00108 slot_dup(ans, x, Matrix_pSym); 00109 REPROTECT(ans = Csparse_transpose(ans, tri), ipx); 00110 SET_DimNames(ans, x); 00111 free(ncl); 00112 UNPROTECT(2); 00113 return ans; 00114 } 00115 00116 SEXP compressed_non_0_ij(SEXP x, SEXP colP) 00117 { 00118 int col = asLogical(colP); /* 1 if "C"olumn compressed; 0 if "R"ow */ 00119 SEXP ans, indSym = col ? Matrix_iSym : Matrix_jSym; 00120 SEXP indP = GET_SLOT(x, indSym), 00121 pP = GET_SLOT(x, Matrix_pSym); 00122 int i, *ij; 00123 int ncol = INTEGER(GET_SLOT(x, Matrix_DimSym))[1], 00124 n_el = INTEGER(pP)[ncol]; /* is only == length(indP), if the 00125 i-slot is not over-allocated */ 00126 00127 ij = INTEGER(ans = PROTECT(allocMatrix(INTSXP, n_el, 2))); 00128 /* expand the compressed margin to 'i' or 'j' : */ 00129 expand_cmprPt(ncol, INTEGER(pP), &ij[col ? n_el : 0]); 00130 /* and copy the other one: */ 00131 if (col) 00132 for(i = 0; i < n_el; i++) 00133 ij[i] = INTEGER(indP)[i]; 00134 else /* row compressed */ 00135 for(i = 0; i < n_el; i++) 00136 ij[i + n_el] = INTEGER(indP)[i]; 00137 00138 UNPROTECT(1); 00139 return ans; 00140 } 00141 00142 #if 0 /* unused */ 00143 SEXP dgCMatrix_lusol(SEXP x, SEXP y) 00144 { 00145 SEXP ycp = PROTECT((TYPEOF(y) == REALSXP) ? 00146 duplicate(y) : coerceVector(y, REALSXP)); 00147 CSP xc = AS_CSP__(x); 00148 R_CheckStack(); 00149 00150 if (xc->m != xc->n || xc->m <= 0) 00151 error(_("dgCMatrix_lusol requires a square, non-empty matrix")); 00152 if (LENGTH(ycp) != xc->m) 00153 error(_("Dimensions of system to be solved are inconsistent")); 00154 if (!cs_lusol(/*order*/ 1, xc, REAL(ycp), /*tol*/ 1e-7)) 00155 error(_("cs_lusol failed")); 00156 00157 UNPROTECT(1); 00158 return ycp; 00159 } 00160 #endif 00161 00162 SEXP dgCMatrix_qrsol(SEXP x, SEXP y, SEXP ord) 00163 { 00164 /* FIXME: extend this to work in multivariate case, i.e. y a matrix with > 1 column ! */ 00165 SEXP ycp = PROTECT((TYPEOF(y) == REALSXP) ? 00166 duplicate(y) : coerceVector(y, REALSXP)); 00167 CSP xc = AS_CSP(x); /* <--> x may be dgC* or dtC* */ 00168 int order = INTEGER(ord)[0]; 00169 #ifdef _not_yet_do_FIXME__ 00170 const char *nms[] = {"L", "coef", "Xty", "resid", ""}; 00171 SEXP ans = PROTECT(Rf_mkNamed(VECSXP, nms)); 00172 #endif 00173 R_CheckStack(); 00174 00175 if (order < 0 || order > 3) 00176 error(_("dgCMatrix_qrsol(., order) needs order in {0,..,3}")); 00177 /* --> cs_amd() --- order 0: natural, 1: Chol, 2: LU, 3: QR */ 00178 if (LENGTH(ycp) != xc->m) 00179 error(_("Dimensions of system to be solved are inconsistent")); 00180 /* FIXME? Note that qr_sol() would allow *under-determined systems; 00181 * In general, we'd need LENGTH(ycp) = max(n,m) 00182 * FIXME also: multivariate y (see above) 00183 */ 00184 if (xc->m < xc->n || xc->n <= 0) 00185 error(_("dgCMatrix_qrsol(<%d x %d>-matrix) requires a 'tall' rectangular matrix"), 00186 xc->m, xc->n); 00187 00188 /* cs_qrsol(): Tim Davis (2006) .. "8.2 Using a QR factorization", p.136f , calling 00189 * ------- cs_sqr(order, ..), see p.76 */ 00190 /* MM: FIXME: write our *OWN* version of - the first case (m >= n) - of cs_qrsol() 00191 * --------- which will (1) work with a *multivariate* y 00192 * (2) compute coefficients properly, not overwriting RHS 00193 */ 00194 if (!cs_qrsol(order, xc, REAL(ycp))) 00195 /* return value really is 0 or 1 - no more info there */ 00196 error(_("cs_qrsol() failed inside dgCMatrix_qrsol()")); 00197 00198 /* Solution is only in the first part of ycp -- cut its length back to n : */ 00199 { 00200 SEXP nms = getAttrib(ycp, R_NamesSymbol); 00201 SETLENGTH(ycp, xc->n); 00202 if(nms != R_NilValue) { 00203 SETLENGTH(nms, xc->n); 00204 setAttrib(ycp, R_NamesSymbol, nms); 00205 } 00206 } 00207 UNPROTECT(1); 00208 return ycp; 00209 } 00210 00211 /* Modified version of Tim Davis's cs_qr_mex.c file for MATLAB */ 00212 SEXP dgCMatrix_QR(SEXP Ap, SEXP order) 00213 { 00214 SEXP ans = PROTECT(NEW_OBJECT(MAKE_CLASS("sparseQR"))); 00215 CSP A = AS_CSP__(Ap), D; 00216 css *S; 00217 csn *N; 00218 int m = A->m, n = A->n, ord = asLogical(order) ? 3 : 0, *p; 00219 int *dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); 00220 R_CheckStack(); 00221 00222 if (m < n) error(_("A must have #{rows} >= #{columns}")) ; 00223 dims[0] = m; dims[1] = n; 00224 S = cs_sqr(ord, A, 1); /* symbolic QR ordering & analysis*/ 00225 if (!S) error(_("cs_sqr failed")); 00226 N = cs_qr(A, S); /* numeric QR factorization */ 00227 if (!N) error(_("cs_qr failed")) ; 00228 cs_dropzeros(N->L); /* drop zeros from V and sort */ 00229 D = cs_transpose(N->L, 1); 00230 cs_spfree(N->L); 00231 N->L = cs_transpose(D, 1); 00232 cs_spfree(D); 00233 cs_dropzeros(N->U); /* drop zeros from R and sort */ 00234 D = cs_transpose(N->U, 1); 00235 cs_spfree(N->U) ; 00236 N->U = cs_transpose(D, 1); 00237 cs_spfree(D); 00238 m = N->L->m; /* m may be larger now */ 00239 p = cs_pinv(S->pinv, m); /* p = pinv' */ 00240 SET_SLOT(ans, install("V"), 00241 Matrix_cs_to_SEXP(N->L, "dgCMatrix", 0)); 00242 Memcpy(REAL(ALLOC_SLOT(ans, install("beta"), 00243 REALSXP, n)), N->B, n); 00244 Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, 00245 INTSXP, m)), p, m); 00246 SET_SLOT(ans, install("R"), 00247 Matrix_cs_to_SEXP(N->U, "dgCMatrix", 0)); 00248 if (ord) 00249 Memcpy(INTEGER(ALLOC_SLOT(ans, install("q"), 00250 INTSXP, n)), S->q, n); 00251 else 00252 ALLOC_SLOT(ans, install("q"), INTSXP, 0); 00253 cs_nfree(N); 00254 cs_sfree(S); 00255 cs_free(p); 00256 UNPROTECT(1); 00257 return ans; 00258 } 00259 00260 #ifdef Matrix_with_SPQR 00261 00281 SEXP dgCMatrix_SPQR(SEXP Ap, SEXP ordering, SEXP econ, SEXP tol) 00282 { 00283 /* SEXP ans = PROTECT(allocVector(VECSXP, 4)); */ 00284 SEXP ans = PROTECT(NEW_OBJECT(MAKE_CLASS("SPQR"))); 00285 00286 CHM_SP A = AS_CHM_SP(Ap), Q, R; 00287 UF_long *E, rank;/* not always = int FIXME (Windows_64 ?) */ 00288 00289 if ((rank = SuiteSparseQR_C_QR(asInteger(ordering), 00290 asReal(tol),/* originally had SPQR_DEFAULT_TOL */ 00291 (UF_long)asInteger(econ),/* originally had 0 */ 00292 A, &Q, &R, &E, &cl)) == -1) 00293 error(_("SuiteSparseQR_C_QR returned an error code")); 00294 00295 SET_SLOT(ans, Matrix_DimSym, duplicate(GET_SLOT(Ap, Matrix_DimSym))); 00296 /* SET_VECTOR_ELT(ans, 0, */ 00297 /* chm_sparse_to_SEXP(Q, 0, 0, 0, "", R_NilValue)); */ 00298 SET_SLOT(ans, install("Q"), 00299 chm_sparse_to_SEXP(Q, 0, 0, 0, "", R_NilValue)); 00300 00301 /* Also gives a dgCMatrix (not a dtC* *triangular*) : 00302 * may make sense if to be used in the "spqr_solve" routines .. ?? */ 00303 /* SET_VECTOR_ELT(ans, 1, */ 00304 /* chm_sparse_to_SEXP(R, 0, 0, 0, "", R_NilValue)); */ 00305 SET_SLOT(ans, install("R"), 00306 chm_sparse_to_SEXP(R, 0, 0, 0, "", R_NilValue)); 00307 cholmod_free_sparse(&Al, &cl); 00308 cholmod_free_sparse(&R, &cl); 00309 cholmod_free_sparse(&Q, &cl); 00310 if (E) { 00311 int *Er; 00312 SET_VECTOR_ELT(ans, 2, allocVector(INTSXP, A->ncol)); 00313 Er = INTEGER(VECTOR_ELT(ans, 2)); 00314 for (int i = 0; i < A->ncol; i++) Er[i] = (int) E[i]; 00315 Free(E); 00316 } else SET_VECTOR_ELT(ans, 2, allocVector(INTSXP, 0)); 00317 SET_VECTOR_ELT(ans, 3, ScalarInteger((int)rank)); 00318 UNPROTECT(1); 00319 return ans; 00320 } 00321 #endif 00322 /* Matrix_with_SPQR */ 00323 00324 /* Modified version of Tim Davis's cs_lu_mex.c file for MATLAB */ 00325 void install_lu(SEXP Ap, int order, double tol, Rboolean err_sing) 00326 { 00327 // (order, tol) == (1, 1) by default, when called from R. 00328 SEXP ans; 00329 css *S; 00330 csn *N; 00331 int n, *p, *dims; 00332 CSP A = AS_CSP__(Ap), D; 00333 R_CheckStack(); 00334 00335 n = A->n; 00336 if (A->m != n) 00337 error(_("LU decomposition applies only to square matrices")); 00338 if (order) { /* not using natural order */ 00339 order = (tol == 1) ? 2 /* amd(S'*S) w/dense rows or I */ 00340 : 1; /* amd (A+A'), or natural */ 00341 } 00342 S = cs_sqr(order, A, /*qr = */ 0); /* symbolic ordering */ 00343 N = cs_lu(A, S, tol); /* numeric factorization */ 00344 if (!N) { 00345 if(err_sing) 00346 error(_("cs_lu(A) failed: near-singular A (or out of memory)")); 00347 else { 00348 /* No warning: The useR should be careful : 00349 * Put NA into "LU" factor */ 00350 set_factors(Ap, ScalarLogical(NA_LOGICAL), "LU"); 00351 return; 00352 } 00353 } 00354 cs_dropzeros(N->L); /* drop zeros from L and sort it */ 00355 D = cs_transpose(N->L, 1); 00356 cs_spfree(N->L); 00357 N->L = cs_transpose(D, 1); 00358 cs_spfree(D); 00359 cs_dropzeros(N->U); /* drop zeros from U and sort it */ 00360 D = cs_transpose(N->U, 1); 00361 cs_spfree(N->U); 00362 N->U = cs_transpose(D, 1); 00363 cs_spfree(D); 00364 p = cs_pinv(N->pinv, n); /* p=pinv' */ 00365 ans = PROTECT(NEW_OBJECT(MAKE_CLASS("sparseLU"))); 00366 dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); 00367 dims[0] = n; dims[1] = n; 00368 SET_SLOT(ans, install("L"), 00369 Matrix_cs_to_SEXP(N->L, "dtCMatrix", 0)); 00370 SET_SLOT(ans, install("U"), 00371 Matrix_cs_to_SEXP(N->U, "dtCMatrix", 0)); 00372 Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, /* "p" */ 00373 INTSXP, n)), p, n); 00374 if (order) 00375 Memcpy(INTEGER(ALLOC_SLOT(ans, install("q"), 00376 INTSXP, n)), S->q, n); 00377 cs_nfree(N); 00378 cs_sfree(S); 00379 cs_free(p); 00380 UNPROTECT(1); 00381 set_factors(Ap, ans, "LU"); 00382 } 00383 00384 SEXP dgCMatrix_LU(SEXP Ap, SEXP orderp, SEXP tolp, SEXP error_on_sing) 00385 { 00386 SEXP ans; 00387 Rboolean err_sing = asLogical(error_on_sing); 00388 /* FIXME: dgCMatrix_LU should check ans for consistency in 00389 * permutation type with the requested value - Should have two 00390 * classes or two different names in the factors list for LU with 00391 * permuted columns or not. 00392 * OTOH, currently (order, tol) === (1, 1) always. 00393 * It is true that length(LU@q) does flag the order argument. 00394 */ 00395 if (!isNull(ans = get_factors(Ap, "LU"))) 00396 return ans; 00397 install_lu(Ap, asInteger(orderp), asReal(tolp), err_sing); 00398 return get_factors(Ap, "LU"); 00399 } 00400 00401 SEXP dgCMatrix_matrix_solve(SEXP Ap, SEXP b) 00402 { 00403 SEXP ans = PROTECT(dup_mMatrix_as_dgeMatrix(b)), 00404 lu, qslot; 00405 CSP L, U; 00406 int *bdims = INTEGER(GET_SLOT(ans, Matrix_DimSym)), *p, *q; 00407 int j, n = bdims[0], nrhs = bdims[1]; 00408 double *ax = REAL(GET_SLOT(ans, Matrix_xSym)), 00409 *x = Alloca(n, double); 00410 R_CheckStack(); 00411 00412 if (isNull(lu = get_factors(Ap, "LU"))) { 00413 install_lu(Ap, /* order = */ 1, /* tol = */ 1.0, /* err_sing = */ TRUE); 00414 lu = get_factors(Ap, "LU"); 00415 } 00416 qslot = GET_SLOT(lu, install("q")); 00417 L = AS_CSP__(GET_SLOT(lu, install("L"))); 00418 U = AS_CSP__(GET_SLOT(lu, install("U"))); 00419 R_CheckStack(); 00420 00421 p = INTEGER(GET_SLOT(lu, Matrix_pSym)); 00422 q = LENGTH(qslot) ? INTEGER(qslot) : (int *) NULL; 00423 00424 if (U->n != n || nrhs < 1 || n < 1) 00425 error(_("Dimensions of system to be solved are inconsistent")); 00426 for (j = 0; j < nrhs; j++) { 00427 cs_pvec(p, ax + j * n, x, n); /* x = b(p) */ 00428 cs_lsolve(L, x); /* x = L\x */ 00429 cs_usolve(U, x); /* x = U\x */ 00430 if (q) /* b(q) = x */ 00431 cs_ipvec(q, x, ax + j * n, n); 00432 else 00433 Memcpy(ax + j * n, x, n); 00434 } 00435 UNPROTECT(1); 00436 return ans; 00437 } 00438 00439 SEXP dgCMatrix_cholsol(SEXP x, SEXP y) 00440 { 00441 /* Solve Sparse Least Squares X %*% beta ~= y with dense RHS y, 00442 * where X = t(x) i.e. we pass x = t(X) as argument, 00443 * via "Cholesky(X'X)" .. well not really: 00444 * cholmod_factorize("x", ..) finds L in X'X = L'L directly */ 00445 CHM_SP cx = AS_CHM_SP(x); 00446 /* FIXME: extend this to work in multivariate case, i.e. y a matrix with > 1 column ! */ 00447 CHM_DN cy = AS_CHM_DN(coerceVector(y, REALSXP)), rhs, cAns, resid; 00448 CHM_FR L; 00449 int n = cx->ncol;/* #{obs.} {x = t(X) !} */ 00450 double one[] = {1,0}, zero[] = {0,0}, neg1[] = {-1,0}; 00451 const char *nms[] = {"L", "coef", "Xty", "resid", ""}; 00452 SEXP ans = PROTECT(Rf_mkNamed(VECSXP, nms)); 00453 R_CheckStack(); 00454 00455 if (n < cx->nrow || n <= 0) 00456 error(_("dgCMatrix_cholsol requires a 'short, wide' rectangular matrix")); 00457 if (cy->nrow != n) 00458 error(_("Dimensions of system to be solved are inconsistent")); 00459 rhs = cholmod_allocate_dense(cx->nrow, 1, cx->nrow, CHOLMOD_REAL, &c); 00460 /* cholmod_sdmult(A, transp, alpha, beta, X, Y, &c): 00461 * Y := alpha*(A*X) + beta*Y or alpha*(A'*X) + beta*Y ; 00462 * here: rhs := 1 * x %*% y + 0 = x %*% y = X'y */ 00463 if (!(cholmod_sdmult(cx, 0 /* trans */, one, zero, cy, rhs, &c))) 00464 error(_("cholmod_sdmult error (rhs)")); 00465 L = cholmod_analyze(cx, &c); 00466 if (!cholmod_factorize(cx, L, &c)) 00467 error(_("cholmod_factorize failed: status %d, minor %d from ncol %d"), 00468 c.status, L->minor, L->n); 00469 /* FIXME: Do this in stages so an "effects" vector can be calculated */ 00470 if (!(cAns = cholmod_solve(CHOLMOD_A, L, rhs, &c))) 00471 error(_("cholmod_solve (CHOLMOD_A) failed: status %d, minor %d from ncol %d"), 00472 c.status, L->minor, L->n); 00473 /* L : */ 00474 SET_VECTOR_ELT(ans, 0, chm_factor_to_SEXP(L, 0)); 00475 /* coef : */ 00476 SET_VECTOR_ELT(ans, 1, allocVector(REALSXP, cx->nrow)); 00477 Memcpy(REAL(VECTOR_ELT(ans, 1)), (double*)(cAns->x), cx->nrow); 00478 /* X'y : */ 00479 /* FIXME: Change this when the "effects" vector is available */ 00480 SET_VECTOR_ELT(ans, 2, allocVector(REALSXP, cx->nrow)); 00481 Memcpy(REAL(VECTOR_ELT(ans, 2)), (double*)(rhs->x), cx->nrow); 00482 /* resid := y */ 00483 resid = cholmod_copy_dense(cy, &c); 00484 /* cholmod_sdmult(A, transp, alp, bet, X, Y, &c): 00485 * Y := alp*(A*X) + bet*Y or alp*(A'*X) + beta*Y ; 00486 * here: resid := -1 * x' %*% coef + 1 * y = y - X %*% coef */ 00487 if (!(cholmod_sdmult(cx, 1/* trans */, neg1, one, cAns, resid, &c))) 00488 error(_("cholmod_sdmult error (resid)")); 00489 /* FIXME: for multivariate case, i.e. resid *matrix* with > 1 column ! */ 00490 SET_VECTOR_ELT(ans, 3, allocVector(REALSXP, n)); 00491 Memcpy(REAL(VECTOR_ELT(ans, 3)), (double*)(resid->x), n); 00492 00493 cholmod_free_factor(&L, &c); 00494 cholmod_free_dense(&rhs, &c); 00495 cholmod_free_dense(&cAns, &c); 00496 UNPROTECT(1); 00497 return ans; 00498 } 00499 00500 00501 /* Define all of 00502 * dgCMatrix_colSums(....) 00503 * igCMatrix_colSums(....) 00504 * lgCMatrix_colSums_d(....) 00505 * lgCMatrix_colSums_i(....) 00506 * ngCMatrix_colSums_d(....) 00507 * ngCMatrix_colSums_i(....) 00508 */ 00509 #define _dgC_ 00510 #include "t_gCMatrix_colSums.c" 00511 00512 #define _igC_ 00513 #include "t_gCMatrix_colSums.c" 00514 00515 #define _lgC_ 00516 #include "t_gCMatrix_colSums.c" 00517 00518 #define _ngC_ 00519 #include "t_gCMatrix_colSums.c" 00520 00521 #define _lgC_mn 00522 #include "t_gCMatrix_colSums.c" 00523 00524 #define _ngC_mn 00525 #include "t_gCMatrix_colSums.c" 00526 00527 00528 SEXP lgCMatrix_colSums(SEXP x, SEXP NArm, SEXP spRes, SEXP trans, SEXP means) 00529 { 00530 if(asLogical(means)) /* ==> result will be "double" / "dsparseVector" */ 00531 return lgCMatrix_colSums_d(x, NArm, spRes, trans, means); 00532 else 00533 return lgCMatrix_colSums_i(x, NArm, spRes, trans, means); 00534 } 00535 00536 SEXP ngCMatrix_colSums(SEXP x, SEXP NArm, SEXP spRes, SEXP trans, SEXP means) 00537 { 00538 if(asLogical(means)) /* ==> result will be "double" / "dsparseVector" */ 00539 return ngCMatrix_colSums_d(x, NArm, spRes, trans, means); 00540 else 00541 return ngCMatrix_colSums_i(x, NArm, spRes, trans, means); 00542 }