|
Matrix $Rev: 2718 $ at $LastChangedDate: 2011-10-06 11:45:17 +0200 (Thu, 06 Oct 2011) $
|
00001 /* Sparse matrices in compressed column-oriented form */ 00002 00003 #include "Csparse.h" 00004 #include "Tsparse.h" 00005 #include "chm_common.h" 00006 00008 Rboolean isValid_Csparse(SEXP x) 00009 { 00010 /* NB: we do *NOT* check a potential 'x' slot here, at all */ 00011 SEXP pslot = GET_SLOT(x, Matrix_pSym), 00012 islot = GET_SLOT(x, Matrix_iSym); 00013 int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), j, 00014 nrow = dims[0], 00015 ncol = dims[1], 00016 *xp = INTEGER(pslot), 00017 *xi = INTEGER(islot); 00018 00019 if (length(pslot) != dims[1] + 1) 00020 return FALSE; 00021 if (xp[0] != 0) 00022 return FALSE; 00023 if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/ 00024 return FALSE; 00025 for (j = 0; j < xp[ncol]; j++) { 00026 if (xi[j] < 0 || xi[j] >= nrow) 00027 return FALSE; 00028 } 00029 for (j = 0; j < ncol; j++) { 00030 if (xp[j] > xp[j + 1]) 00031 return FALSE; 00032 } 00033 return TRUE; 00034 } 00035 00036 SEXP Csparse_validate(SEXP x) { 00037 return Csparse_validate_(x, FALSE); 00038 } 00039 00040 SEXP Csparse_validate2(SEXP x, SEXP maybe_modify) { 00041 return Csparse_validate_(x, asLogical(maybe_modify)); 00042 } 00043 00044 SEXP Csparse_validate_(SEXP x, Rboolean maybe_modify) 00045 { 00046 /* NB: we do *NOT* check a potential 'x' slot here, at all */ 00047 SEXP pslot = GET_SLOT(x, Matrix_pSym), 00048 islot = GET_SLOT(x, Matrix_iSym); 00049 Rboolean sorted, strictly; 00050 int j, k, 00051 *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), 00052 nrow = dims[0], 00053 ncol = dims[1], 00054 *xp = INTEGER(pslot), 00055 *xi = INTEGER(islot); 00056 00057 if (length(pslot) != dims[1] + 1) 00058 return mkString(_("slot p must have length = ncol(.) + 1")); 00059 if (xp[0] != 0) 00060 return mkString(_("first element of slot p must be zero")); 00061 if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/ 00062 return 00063 mkString(_("last element of slot p must match length of slots i and x")); 00064 for (j = 0; j < xp[ncol]; j++) { 00065 if (xi[j] < 0 || xi[j] >= nrow) 00066 return mkString(_("all row indices must be between 0 and nrow-1")); 00067 } 00068 sorted = TRUE; strictly = TRUE; 00069 for (j = 0; j < ncol; j++) { 00070 if (xp[j] > xp[j + 1]) 00071 return mkString(_("slot p must be non-decreasing")); 00072 if(sorted) /* only act if >= 2 entries in column j : */ 00073 for (k = xp[j] + 1; k < xp[j + 1]; k++) { 00074 if (xi[k] < xi[k - 1]) 00075 sorted = FALSE; 00076 else if (xi[k] == xi[k - 1]) 00077 strictly = FALSE; 00078 } 00079 } 00080 if (!sorted) { 00081 if(maybe_modify) { 00082 CHM_SP chx = (CHM_SP) alloca(sizeof(cholmod_sparse)); 00083 R_CheckStack(); 00084 as_cholmod_sparse(chx, x, FALSE, TRUE);/*-> cholmod_l_sort() ! */ 00085 /* as chx = AS_CHM_SP__(x) but ^^^^ sorting x in_place !!! */ 00086 00087 /* Now re-check that row indices are *strictly* increasing 00088 * (and not just increasing) within each column : */ 00089 for (j = 0; j < ncol; j++) { 00090 for (k = xp[j] + 1; k < xp[j + 1]; k++) 00091 if (xi[k] == xi[k - 1]) 00092 return mkString(_("slot i is not *strictly* increasing inside a column (even after cholmod_l_sort)")); 00093 } 00094 } else { /* no modifying sorting : */ 00095 return mkString(_("row indices are not sorted within columns")); 00096 } 00097 } else if(!strictly) { /* sorted, but not strictly */ 00098 return mkString(_("slot i is not *strictly* increasing inside a column")); 00099 } 00100 return ScalarLogical(1); 00101 } 00102 00103 SEXP Rsparse_validate(SEXP x) 00104 { 00105 /* NB: we do *NOT* check a potential 'x' slot here, at all */ 00106 SEXP pslot = GET_SLOT(x, Matrix_pSym), 00107 jslot = GET_SLOT(x, Matrix_jSym); 00108 Rboolean sorted, strictly; 00109 int i, k, 00110 *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), 00111 nrow = dims[0], 00112 ncol = dims[1], 00113 *xp = INTEGER(pslot), 00114 *xj = INTEGER(jslot); 00115 00116 if (length(pslot) != dims[0] + 1) 00117 return mkString(_("slot p must have length = nrow(.) + 1")); 00118 if (xp[0] != 0) 00119 return mkString(_("first element of slot p must be zero")); 00120 if (length(jslot) < xp[nrow]) /* allow larger slots from over-allocation!*/ 00121 return 00122 mkString(_("last element of slot p must match length of slots j and x")); 00123 for (i = 0; i < length(jslot); i++) { 00124 if (xj[i] < 0 || xj[i] >= ncol) 00125 return mkString(_("all column indices must be between 0 and ncol-1")); 00126 } 00127 sorted = TRUE; strictly = TRUE; 00128 for (i = 0; i < nrow; i++) { 00129 if (xp[i] > xp[i+1]) 00130 return mkString(_("slot p must be non-decreasing")); 00131 if(sorted) 00132 for (k = xp[i] + 1; k < xp[i + 1]; k++) { 00133 if (xj[k] < xj[k - 1]) 00134 sorted = FALSE; 00135 else if (xj[k] == xj[k - 1]) 00136 strictly = FALSE; 00137 } 00138 } 00139 if (!sorted) 00140 /* cannot easily use cholmod_sort(.) ... -> "error out" :*/ 00141 return mkString(_("slot j is not increasing inside a column")); 00142 else if(!strictly) /* sorted, but not strictly */ 00143 return mkString(_("slot j is not *strictly* increasing inside a column")); 00144 00145 return ScalarLogical(1); 00146 } 00147 00148 00149 /* Called from ../R/Csparse.R : */ 00150 /* Can only return [dln]geMatrix (no symm/triang); 00151 * FIXME: replace by non-CHOLMOD code ! */ 00152 SEXP Csparse_to_dense(SEXP x) 00153 { 00154 CHM_SP chxs = AS_CHM_SP__(x); 00155 /* This loses the symmetry property, since cholmod_dense has none, 00156 * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices 00157 * to numeric (CHOLMOD_REAL) ones : */ 00158 CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c); 00159 int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x); 00160 R_CheckStack(); 00161 00162 return chm_dense_to_SEXP(chxd, 1, Rkind, GET_SLOT(x, Matrix_DimNamesSym)); 00163 } 00164 00165 // FIXME: do not go via CHM (should not be too hard, to just *drop* the x-slot, right? 00166 SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri) 00167 { 00168 CHM_SP chxs = AS_CHM_SP__(x); 00169 CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c); 00170 int tr = asLogical(tri); 00171 R_CheckStack(); 00172 00173 return chm_sparse_to_SEXP(chxcp, 1/*do_free*/, 00174 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, 00175 0, tr ? diag_P(x) : "", 00176 GET_SLOT(x, Matrix_DimNamesSym)); 00177 } 00178 00179 // n.CMatrix --> [dli].CMatrix (not going through CHM!) 00180 SEXP nz_pattern_to_Csparse(SEXP x, SEXP res_kind) 00181 { 00182 return nz2Csparse(x, asInteger(res_kind)); 00183 } 00184 // n.CMatrix --> [dli].CMatrix (not going through CHM!) 00185 SEXP nz2Csparse(SEXP x, enum x_slot_kind r_kind) 00186 { 00187 const char *cl_x = class_P(x); 00188 if(cl_x[0] != 'n') error(_("not a 'n.CMatrix'")); 00189 if(cl_x[2] != 'C') error(_("not a CsparseMatrix")); 00190 int nnz = LENGTH(GET_SLOT(x, Matrix_iSym)); 00191 SEXP ans; 00192 char *ncl = strdup(cl_x); 00193 double *dx_x; int *ix_x; 00194 ncl[0] = (r_kind == x_double ? 'd' : 00195 (r_kind == x_logical ? 'l' : 00196 /* else (for now): r_kind == x_integer : */ 'i')); 00197 PROTECT(ans = NEW_OBJECT(MAKE_CLASS(ncl))); 00198 // create a correct 'x' slot: 00199 switch(r_kind) { 00200 int i; 00201 case x_double: // 'd' 00202 dx_x = REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)); 00203 for (i=0; i < nnz; i++) dx_x[i] = 1.; 00204 break; 00205 case x_logical: // 'l' 00206 ix_x = LOGICAL(ALLOC_SLOT(ans, Matrix_xSym, LGLSXP, nnz)); 00207 for (i=0; i < nnz; i++) ix_x[i] = TRUE; 00208 break; 00209 case x_integer: // 'i' 00210 ix_x = INTEGER(ALLOC_SLOT(ans, Matrix_xSym, INTSXP, nnz)); 00211 for (i=0; i < nnz; i++) ix_x[i] = 1; 00212 break; 00213 00214 default: 00215 error(_("nz2Csparse(): invalid/non-implemented r_kind = %d"), 00216 r_kind); 00217 } 00218 00219 // now copy all other slots : 00220 slot_dup(ans, x, Matrix_iSym); 00221 slot_dup(ans, x, Matrix_pSym); 00222 slot_dup(ans, x, Matrix_DimSym); 00223 slot_dup(ans, x, Matrix_DimNamesSym); 00224 if(ncl[1] != 'g') { // symmetric or triangular ... 00225 slot_dup_if_has(ans, x, Matrix_uploSym); 00226 slot_dup_if_has(ans, x, Matrix_diagSym); 00227 } 00228 UNPROTECT(1); 00229 return ans; 00230 } 00231 00232 SEXP Csparse_to_matrix(SEXP x) 00233 { 00234 return chm_dense_to_matrix(cholmod_sparse_to_dense(AS_CHM_SP__(x), &c), 00235 1 /*do_free*/, GET_SLOT(x, Matrix_DimNamesSym)); 00236 } 00237 00238 SEXP Csparse_to_Tsparse(SEXP x, SEXP tri) 00239 { 00240 CHM_SP chxs = AS_CHM_SP__(x); 00241 CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c); 00242 int tr = asLogical(tri); 00243 int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00244 R_CheckStack(); 00245 00246 return chm_triplet_to_SEXP(chxt, 1, 00247 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, 00248 Rkind, tr ? diag_P(x) : "", 00249 GET_SLOT(x, Matrix_DimNamesSym)); 00250 } 00251 00252 /* this used to be called sCMatrix_to_gCMatrix(..) [in ./dsCMatrix.c ]: */ 00253 SEXP Csparse_symmetric_to_general(SEXP x) 00254 { 00255 CHM_SP chx = AS_CHM_SP__(x), chgx; 00256 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00257 R_CheckStack(); 00258 00259 if (!(chx->stype)) 00260 error(_("Nonsymmetric matrix in Csparse_symmetric_to_general")); 00261 chgx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c); 00262 /* xtype: pattern, "real", complex or .. */ 00263 return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "", 00264 GET_SLOT(x, Matrix_DimNamesSym)); 00265 } 00266 00267 SEXP Csparse_general_to_symmetric(SEXP x, SEXP uplo) 00268 { 00269 CHM_SP chx = AS_CHM_SP__(x), chgx; 00270 int uploT = (*CHAR(STRING_ELT(uplo,0)) == 'U') ? 1 : -1; 00271 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00272 R_CheckStack(); 00273 00274 chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c); 00275 /* xtype: pattern, "real", complex or .. */ 00276 return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "", 00277 GET_SLOT(x, Matrix_DimNamesSym)); 00278 } 00279 00280 SEXP Csparse_transpose(SEXP x, SEXP tri) 00281 { 00282 /* TODO: lgCMatrix & igC* currently go via double prec. cholmod - 00283 * since cholmod (& cs) lacks sparse 'int' matrices */ 00284 CHM_SP chx = AS_CHM_SP__(x); 00285 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00286 CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c); 00287 SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp; 00288 int tr = asLogical(tri); 00289 R_CheckStack(); 00290 00291 tmp = VECTOR_ELT(dn, 0); /* swap the dimnames */ 00292 SET_VECTOR_ELT(dn, 0, VECTOR_ELT(dn, 1)); 00293 SET_VECTOR_ELT(dn, 1, tmp); 00294 UNPROTECT(1); 00295 return chm_sparse_to_SEXP(chxt, 1, /* SWAP 'uplo' for triangular */ 00296 tr ? ((*uplo_P(x) == 'U') ? -1 : 1) : 0, 00297 Rkind, tr ? diag_P(x) : "", dn); 00298 } 00299 00300 SEXP Csparse_Csparse_prod(SEXP a, SEXP b) 00301 { 00302 CHM_SP 00303 cha = AS_CHM_SP(a), 00304 chb = AS_CHM_SP(b), 00305 chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0, 00306 /* values:= is_numeric (T/F) */ cha->xtype > 0, 00307 /*out sorted:*/ 1, &c); 00308 const char *cl_a = class_P(a), *cl_b = class_P(b); 00309 char diag[] = {'\0', '\0'}; 00310 int uploT = 0; 00311 SEXP dn = PROTECT(allocVector(VECSXP, 2)); 00312 R_CheckStack(); 00313 00314 #ifdef DEBUG_Matrix_verbose 00315 Rprintf("DBG Csparse_C*_prod(%s, %s)\n", cl_a, cl_b); 00316 #endif 00317 00318 /* Preserve triangularity and even unit-triangularity if appropriate. 00319 * Note that in that case, the multiplication itself should happen 00320 * faster. But there's no support for that in CHOLMOD */ 00321 00322 /* UGLY hack -- rather should have (fast!) C-level version of 00323 * is(a, "triangularMatrix") etc */ 00324 if (cl_a[1] == 't' && cl_b[1] == 't') 00325 /* FIXME: fails for "Cholesky","BunchKaufmann"..*/ 00326 if(*uplo_P(a) == *uplo_P(b)) { /* both upper, or both lower tri. */ 00327 uploT = (*uplo_P(a) == 'U') ? 1 : -1; 00328 if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */ 00329 /* "remove the diagonal entries": */ 00330 chm_diagN2U(chc, uploT, /* do_realloc */ FALSE); 00331 diag[0]= 'U'; 00332 } 00333 else diag[0]= 'N'; 00334 } 00335 SET_VECTOR_ELT(dn, 0, /* establish dimnames */ 00336 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0))); 00337 SET_VECTOR_ELT(dn, 1, 00338 duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), 1))); 00339 UNPROTECT(1); 00340 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn); 00341 } 00342 00343 SEXP Csparse_Csparse_crossprod(SEXP a, SEXP b, SEXP trans) 00344 { 00345 int tr = asLogical(trans); 00346 CHM_SP 00347 cha = AS_CHM_SP(a), 00348 chb = AS_CHM_SP(b), 00349 chTr, chc; 00350 const char *cl_a = class_P(a), *cl_b = class_P(b); 00351 char diag[] = {'\0', '\0'}; 00352 int uploT = 0; 00353 SEXP dn = PROTECT(allocVector(VECSXP, 2)); 00354 R_CheckStack(); 00355 00356 chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c); 00357 chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb, 00358 /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c); 00359 cholmod_free_sparse(&chTr, &c); 00360 00361 /* Preserve triangularity and unit-triangularity if appropriate; 00362 * see Csparse_Csparse_prod() for comments */ 00363 if (cl_a[1] == 't' && cl_b[1] == 't') 00364 if(*uplo_P(a) != *uplo_P(b)) { /* one 'U', the other 'L' */ 00365 uploT = (*uplo_P(b) == 'U') ? 1 : -1; 00366 if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */ 00367 chm_diagN2U(chc, uploT, /* do_realloc */ FALSE); 00368 diag[0]= 'U'; 00369 } 00370 else diag[0]= 'N'; 00371 } 00372 SET_VECTOR_ELT(dn, 0, /* establish dimnames */ 00373 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), (tr) ? 0 : 1))); 00374 SET_VECTOR_ELT(dn, 1, 00375 duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), (tr) ? 0 : 1))); 00376 UNPROTECT(1); 00377 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn); 00378 } 00379 00380 SEXP Csparse_dense_prod(SEXP a, SEXP b) 00381 { 00382 CHM_SP cha = AS_CHM_SP(a); 00383 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b)); 00384 CHM_DN chb = AS_CHM_DN(b_M); 00385 CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow, 00386 chb->xtype, &c); 00387 SEXP dn = PROTECT(allocVector(VECSXP, 2)); 00388 double one[] = {1,0}, zero[] = {0,0}; 00389 int nprot = 2; 00390 R_CheckStack(); 00391 /* Tim Davis, please FIXME: currently (2010-11) *fails* when a is a pattern matrix:*/ 00392 if(cha->xtype == CHOLMOD_PATTERN) { 00393 /* warning(_("Csparse_dense_prod(): cholmod_sdmult() not yet implemented for pattern./ ngCMatrix" */ 00394 /* " --> slightly inefficient coercion")); */ 00395 00396 // This *fails* to produce a CHOLMOD_REAL .. 00397 // CHM_SP chd = cholmod_l_copy(cha, cha->stype, CHOLMOD_REAL, &c); 00398 // --> use our Matrix-classes 00399 SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++; 00400 cha = AS_CHM_SP(da); 00401 } 00402 cholmod_sdmult(cha, 0, one, zero, chb, chc, &c); 00403 SET_VECTOR_ELT(dn, 0, /* establish dimnames */ 00404 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0))); 00405 SET_VECTOR_ELT(dn, 1, 00406 duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1))); 00407 UNPROTECT(nprot); 00408 return chm_dense_to_SEXP(chc, 1, 0, dn); 00409 } 00410 00411 SEXP Csparse_dense_crossprod(SEXP a, SEXP b) 00412 { 00413 CHM_SP cha = AS_CHM_SP(a); 00414 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b)); 00415 CHM_DN chb = AS_CHM_DN(b_M); 00416 CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol, 00417 chb->xtype, &c); 00418 SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2; 00419 double one[] = {1,0}, zero[] = {0,0}; 00420 R_CheckStack(); 00421 // -- see Csparse_dense_prod() above : 00422 if(cha->xtype == CHOLMOD_PATTERN) { 00423 SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++; 00424 cha = AS_CHM_SP(da); 00425 } 00426 cholmod_sdmult(cha, 1, one, zero, chb, chc, &c); 00427 SET_VECTOR_ELT(dn, 0, /* establish dimnames */ 00428 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1))); 00429 SET_VECTOR_ELT(dn, 1, 00430 duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1))); 00431 UNPROTECT(nprot); 00432 return chm_dense_to_SEXP(chc, 1, 0, dn); 00433 } 00434 00435 /* Computes x'x or x x' -- *also* for Tsparse (triplet = TRUE) 00436 see Csparse_Csparse_crossprod above for x'y and x y' */ 00437 SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet) 00438 { 00439 int trip = asLogical(triplet), 00440 tr = asLogical(trans); /* gets reversed because _aat is tcrossprod */ 00441 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY 00442 CHM_TR cht = trip ? AS_CHM_TR(x) : (CHM_TR) NULL; 00443 #else /* workaround needed:*/ 00444 SEXP xx = PROTECT(Tsparse_diagU2N(x)); 00445 CHM_TR cht = trip ? AS_CHM_TR__(xx) : (CHM_TR) NULL; 00446 #endif 00447 CHM_SP chcp, chxt, 00448 chx = (trip ? 00449 cholmod_triplet_to_sparse(cht, cht->nnz, &c) : 00450 AS_CHM_SP(x)); 00451 SEXP dn = PROTECT(allocVector(VECSXP, 2)); 00452 R_CheckStack(); 00453 00454 if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c); 00455 chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c); 00456 if(!chcp) { 00457 UNPROTECT(1); 00458 error(_("Csparse_crossprod(): error return from cholmod_aat()")); 00459 } 00460 cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c); 00461 chcp->stype = 1; 00462 if (trip) cholmod_free_sparse(&chx, &c); 00463 if (!tr) cholmod_free_sparse(&chxt, &c); 00464 SET_VECTOR_ELT(dn, 0, /* establish dimnames */ 00465 duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym), 00466 (tr) ? 0 : 1))); 00467 SET_VECTOR_ELT(dn, 1, duplicate(VECTOR_ELT(dn, 0))); 00468 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY 00469 UNPROTECT(1); 00470 #else 00471 UNPROTECT(2); 00472 #endif 00473 return chm_sparse_to_SEXP(chcp, 1, 0, 0, "", dn); 00474 } 00475 00476 /* Csparse_drop(x, tol): drop entries with absolute value < tol, i.e, 00477 * at least all "explicit" zeros */ 00478 SEXP Csparse_drop(SEXP x, SEXP tol) 00479 { 00480 const char *cl = class_P(x); 00481 /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ 00482 int tr = (cl[1] == 't'); 00483 CHM_SP chx = AS_CHM_SP__(x); 00484 CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c); 00485 double dtol = asReal(tol); 00486 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00487 R_CheckStack(); 00488 00489 if(!cholmod_drop(dtol, ans, &c)) 00490 error(_("cholmod_drop() failed")); 00491 return chm_sparse_to_SEXP(ans, 1, 00492 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, 00493 Rkind, tr ? diag_P(x) : "", 00494 GET_SLOT(x, Matrix_DimNamesSym)); 00495 } 00496 00497 SEXP Csparse_horzcat(SEXP x, SEXP y) 00498 { 00499 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y); 00500 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0, 00501 Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0, 00502 Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0; 00503 R_CheckStack(); 00504 00505 /* TODO: currently drops dimnames - and we fix at R level */ 00506 return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c), 00507 1, 0, Rkind, "", R_NilValue); 00508 } 00509 00510 SEXP Csparse_vertcat(SEXP x, SEXP y) 00511 { 00512 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y); 00513 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0, 00514 Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0, 00515 Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0; 00516 R_CheckStack(); 00517 00518 /* TODO: currently drops dimnames - and we fix at R level */ 00519 return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c), 00520 1, 0, Rkind, "", R_NilValue); 00521 } 00522 00523 SEXP Csparse_band(SEXP x, SEXP k1, SEXP k2) 00524 { 00525 CHM_SP chx = AS_CHM_SP__(x); 00526 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00527 CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c); 00528 R_CheckStack(); 00529 00530 return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "", 00531 GET_SLOT(x, Matrix_DimNamesSym)); 00532 } 00533 00534 SEXP Csparse_diagU2N(SEXP x) 00535 { 00536 const char *cl = class_P(x); 00537 /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ 00538 if (cl[1] != 't' || *diag_P(x) != 'U') { 00539 /* "trivially fast" when not triangular (<==> no 'diag' slot), 00540 or not *unit* triangular */ 00541 return (x); 00542 } 00543 else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */ 00544 CHM_SP chx = AS_CHM_SP__(x); 00545 CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c); 00546 double one[] = {1, 0}; 00547 CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c); 00548 int uploT = (*uplo_P(x) == 'U') ? 1 : -1; 00549 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00550 00551 R_CheckStack(); 00552 cholmod_free_sparse(&eye, &c); 00553 return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N", 00554 GET_SLOT(x, Matrix_DimNamesSym)); 00555 } 00556 } 00557 00558 SEXP Csparse_diagN2U(SEXP x) 00559 { 00560 const char *cl = class_P(x); 00561 /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ 00562 if (cl[1] != 't' || *diag_P(x) != 'N') { 00563 /* "trivially fast" when not triangular (<==> no 'diag' slot), 00564 or already *unit* triangular */ 00565 return (x); 00566 } 00567 else { /* triangular with diag='N'): now drop the diagonal */ 00568 /* duplicate, since chx will be modified: */ 00569 CHM_SP chx = AS_CHM_SP__(duplicate(x)); 00570 int uploT = (*uplo_P(x) == 'U') ? 1 : -1, 00571 Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00572 R_CheckStack(); 00573 00574 chm_diagN2U(chx, uploT, /* do_realloc */ FALSE); 00575 00576 return chm_sparse_to_SEXP(chx, /*dofree*/ 0/* or 1 ?? */, 00577 uploT, Rkind, "U", 00578 GET_SLOT(x, Matrix_DimNamesSym)); 00579 } 00580 } 00581 00591 SEXP Csparse_submatrix(SEXP x, SEXP i, SEXP j) 00592 { 00593 CHM_SP chx = AS_CHM_SP(x); /* << does diagU2N() when needed */ 00594 int rsize = (isNull(i)) ? -1 : LENGTH(i), 00595 csize = (isNull(j)) ? -1 : LENGTH(j); 00596 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; 00597 R_CheckStack(); 00598 00599 if (rsize >= 0 && !isInteger(i)) 00600 error(_("Index i must be NULL or integer")); 00601 if (csize >= 0 && !isInteger(j)) 00602 error(_("Index j must be NULL or integer")); 00603 00604 if (chx->stype) /* symmetricMatrix */ 00605 /* for now, cholmod_submatrix() only accepts "generalMatrix" */ 00606 chx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c); 00607 00608 return chm_sparse_to_SEXP(cholmod_submatrix(chx, 00609 (rsize < 0) ? NULL : INTEGER(i), rsize, 00610 (csize < 0) ? NULL : INTEGER(j), csize, 00611 TRUE, TRUE, &c), 00612 1, 0, Rkind, "", 00613 /* FIXME: drops dimnames */ R_NilValue); 00614 } 00615 00616 #define _d_Csp_ 00617 #include "t_Csparse_subassign.c" 00618 00619 #define _l_Csp_ 00620 #include "t_Csparse_subassign.c" 00621 00622 #define _i_Csp_ 00623 #include "t_Csparse_subassign.c" 00624 00625 #define _n_Csp_ 00626 #include "t_Csparse_subassign.c" 00627 00628 #define _z_Csp_ 00629 #include "t_Csparse_subassign.c" 00630 00631 00632 00633 SEXP Csparse_MatrixMarket(SEXP x, SEXP fname) 00634 { 00635 FILE *f = fopen(CHAR(asChar(fname)), "w"); 00636 00637 if (!f) 00638 error(_("failure to open file \"%s\" for writing"), 00639 CHAR(asChar(fname))); 00640 if (!cholmod_write_sparse(f, AS_CHM_SP(x), 00641 (CHM_SP)NULL, (char*) NULL, &c)) 00642 error(_("cholmod_write_sparse returned error code")); 00643 fclose(f); 00644 return R_NilValue; 00645 } 00646 00647 00660 SEXP diag_tC_ptr(int n, int *x_p, double *x_x, int *perm, SEXP resultKind) 00661 /* ^^^^^^ FIXME[Generalize] to int / ... */ 00662 { 00663 const char* res_ch = CHAR(STRING_ELT(resultKind,0)); 00664 enum diag_kind { diag, diag_backpermuted, trace, prod, sum_log 00665 } res_kind = ((!strcmp(res_ch, "trace")) ? trace : 00666 ((!strcmp(res_ch, "sumLog")) ? sum_log : 00667 ((!strcmp(res_ch, "prod")) ? prod : 00668 ((!strcmp(res_ch, "diag")) ? diag : 00669 ((!strcmp(res_ch, "diagBack")) ? diag_backpermuted : 00670 -1))))); 00671 int i, n_x, i_from = 0; 00672 SEXP ans = PROTECT(allocVector(REALSXP, 00673 /* ^^^^ FIXME[Generalize] */ 00674 (res_kind == diag || 00675 res_kind == diag_backpermuted) ? n : 1)); 00676 double *v = REAL(ans); 00677 /* ^^^^^^ ^^^^ FIXME[Generalize] */ 00678 00679 #define for_DIAG(v_ASSIGN) \ 00680 for(i = 0; i < n; i++, i_from += n_x) { \ 00681 /* looking at i-th column */ \ 00682 n_x = x_p[i+1] - x_p[i];/* #{entries} in this column */ \ 00683 v_ASSIGN; \ 00684 } 00685 00686 /* NOTA BENE: we assume -- uplo = "L" i.e. lower triangular matrix 00687 * for uplo = "U" (makes sense with a "dtCMatrix" !), 00688 * should use x_x[i_from + (nx - 1)] instead of x_x[i_from], 00689 * where nx = (x_p[i+1] - x_p[i]) 00690 */ 00691 00692 switch(res_kind) { 00693 case trace: 00694 v[0] = 0.; 00695 for_DIAG(v[0] += x_x[i_from]); 00696 break; 00697 00698 case sum_log: 00699 v[0] = 0.; 00700 for_DIAG(v[0] += log(x_x[i_from])); 00701 break; 00702 00703 case prod: 00704 v[0] = 1.; 00705 for_DIAG(v[0] *= x_x[i_from]); 00706 break; 00707 00708 case diag: 00709 for_DIAG(v[i] = x_x[i_from]); 00710 break; 00711 00712 case diag_backpermuted: 00713 for_DIAG(v[i] = x_x[i_from]); 00714 00715 warning(_("resultKind = 'diagBack' (back-permuted) is experimental")); 00716 /* now back_permute : */ 00717 for(i = 0; i < n; i++) { 00718 double tmp = v[i]; v[i] = v[perm[i]]; v[perm[i]] = tmp; 00719 /*^^^^ FIXME[Generalize] */ 00720 } 00721 break; 00722 00723 default: /* -1 from above */ 00724 error(_("diag_tC(): invalid 'resultKind'")); 00725 /* Wall: */ ans = R_NilValue; v = REAL(ans); 00726 } 00727 00728 UNPROTECT(1); 00729 return ans; 00730 } 00731 00744 SEXP diag_tC(SEXP pslot, SEXP xslot, SEXP perm_slot, SEXP resultKind) 00745 { 00746 int n = length(pslot) - 1, /* n = ncol(.) = nrow(.) */ 00747 *x_p = INTEGER(pslot), 00748 *perm = INTEGER(perm_slot); 00749 double *x_x = REAL(xslot); 00750 /* ^^^^^^ ^^^^ FIXME[Generalize] to INTEGER(.) / LOGICAL(.) / ... xslot !*/ 00751 00752 return diag_tC_ptr(n, x_p, x_x, perm, resultKind); 00753 } 00754 00773 SEXP create_Csparse(char* cls, int* i, int* j, int* p, int np, 00774 void* x, int nnz, int* dims, SEXP dimnames, 00775 int index1) 00776 { 00777 SEXP ans; 00778 int *ij = (int*)NULL, *tri, *trj, 00779 mi, mj, mp, nrow = -1, ncol = -1; 00780 int xtype = -1; /* -Wall */ 00781 CHM_TR T; 00782 CHM_SP A; 00783 00784 if (np < 0 || nnz < 0) 00785 error(_("negative vector lengths not allowed: np = %d, nnz = %d"), 00786 np, nnz); 00787 if (1 != ((mi = (i == (int*)NULL)) + 00788 (mj = (j == (int*)NULL)) + 00789 (mp = (p == (int*)NULL)))) 00790 error(_("exactly 1 of 'i', 'j' or 'p' must be NULL")); 00791 if (mp) { 00792 if (np) error(_("np = %d, must be zero when p is NULL"), np); 00793 } else { 00794 if (np) { /* Expand p to form i or j */ 00795 if (!(p[0])) error(_("p[0] = %d, should be zero"), p[0]); 00796 for (int ii = 0; ii < np; ii++) 00797 if (p[ii] > p[ii + 1]) 00798 error(_("p must be non-decreasing")); 00799 if (p[np] != nnz) 00800 error("p[np] = %d != nnz = %d", p[np], nnz); 00801 ij = Calloc(nnz, int); 00802 if (mi) { 00803 i = ij; 00804 nrow = np; 00805 } else { 00806 j = ij; 00807 ncol = np; 00808 } 00809 /* Expand p to 0-based indices */ 00810 for (int ii = 0; ii < np; ii++) 00811 for (int jj = p[ii]; jj < p[ii + 1]; jj++) ij[jj] = ii; 00812 } else { 00813 if (nnz) 00814 error(_("Inconsistent dimensions: np = 0 and nnz = %d"), 00815 nnz); 00816 } 00817 } 00818 /* calculate nrow and ncol */ 00819 if (nrow < 0) { 00820 for (int ii = 0; ii < nnz; ii++) { 00821 int i1 = i[ii] + (index1 ? 0 : 1); /* 1-based index */ 00822 if (i1 < 1) error(_("invalid row index at position %d"), ii); 00823 if (i1 > nrow) nrow = i1; 00824 } 00825 } 00826 if (ncol < 0) { 00827 for (int jj = 0; jj < nnz; jj++) { 00828 int j1 = j[jj] + (index1 ? 0 : 1); 00829 if (j1 < 1) error(_("invalid column index at position %d"), jj); 00830 if (j1 > ncol) ncol = j1; 00831 } 00832 } 00833 if (dims != (int*)NULL) { 00834 if (dims[0] > nrow) nrow = dims[0]; 00835 if (dims[1] > ncol) ncol = dims[1]; 00836 } 00837 /* check the class name */ 00838 if (strlen(cls) != 8) 00839 error(_("strlen of cls argument = %d, should be 8"), strlen(cls)); 00840 if (!strcmp(cls + 2, "CMatrix")) 00841 error(_("cls = \"%s\" does not end in \"CMatrix\""), cls); 00842 switch(cls[0]) { 00843 case 'd': 00844 case 'l': 00845 xtype = CHOLMOD_REAL; 00846 break; 00847 case 'n': 00848 xtype = CHOLMOD_PATTERN; 00849 break; 00850 default: 00851 error(_("cls = \"%s\" must begin with 'd', 'l' or 'n'"), cls); 00852 } 00853 if (cls[1] != 'g') 00854 error(_("Only 'g'eneral sparse matrix types allowed")); 00855 /* allocate and populate the triplet */ 00856 T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0, 00857 xtype, &c); 00858 T->x = x; 00859 tri = (int*)T->i; 00860 trj = (int*)T->j; 00861 for (int ii = 0; ii < nnz; ii++) { 00862 tri[ii] = i[ii] - ((!mi && index1) ? 1 : 0); 00863 trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0); 00864 } 00865 /* create the cholmod_sparse structure */ 00866 A = cholmod_triplet_to_sparse(T, nnz, &c); 00867 cholmod_free_triplet(&T, &c); 00868 /* copy the information to the SEXP */ 00869 ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls))); 00870 /* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */ 00871 /* allocate and copy common slots */ 00872 nnz = cholmod_nnz(A, &c); 00873 dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); 00874 dims[0] = A->nrow; dims[1] = A->ncol; 00875 Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1); 00876 Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), (int*)A->i, nnz); 00877 switch(cls[1]) { 00878 case 'd': 00879 Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), (double*)A->x, nnz); 00880 break; 00881 case 'l': 00882 error(_("code not yet written for cls = \"lgCMatrix\"")); 00883 } 00884 /* FIXME: dimnames are *NOT* put there yet (if non-NULL) */ 00885 cholmod_free_sparse(&A, &c); 00886 UNPROTECT(1); 00887 return ans; 00888 }