SEXP csc_tcrossprod(SEXP x) { cholmod_sparse *cha = cholmod_aat(as_cholmod_sparse(x), (int *) NULL, 0, 1, &c); cha->stype = -1; /* set the symmetry */ cholmod_sort(cha, &c); /* drop redundant entries */ return chm_sparse_to_SEXP(cha, -1); }
/** * Populate ans with the pointers from x and modify its scalar * elements accordingly. Note that later changes to the contents of * ans will change the contents of the SEXP. * * In most cases this function is called through the macros * AS_CHM_SP() or AS_CHM_SP__(). It is unusual to call it directly. * * @param ans a CHM_SP pointer * @param x pointer to an object that inherits from CsparseMatrix * @param check_Udiag boolean - should a check for (and consequent * expansion of) a unit diagonal be performed. * @param sort_in_place boolean - if the i and x slots are to be sorted * should they be sorted in place? If the i and x slots are pointers * to an input SEXP they should not be modified. * * @return ans containing pointers to the slots of x, *unless* * check_Udiag and x is unitriangular. */ CHM_SP as_cholmod_sparse(CHM_SP ans, SEXP x, Rboolean check_Udiag, Rboolean sort_in_place) { static const char *valid[] = { MATRIX_VALID_Csparse, ""}; int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), ctype = R_check_class_etc(x, valid); SEXP islot = GET_SLOT(x, Matrix_iSym); if (ctype < 0) error(_("invalid class of object to as_cholmod_sparse")); if (!isValid_Csparse(x)) error(_("invalid object passed to as_cholmod_sparse")); memset(ans, 0, sizeof(cholmod_sparse)); /* zero the struct */ ans->itype = CHOLMOD_INT; /* characteristics of the system */ ans->dtype = CHOLMOD_DOUBLE; ans->packed = TRUE; /* slots always present */ ans->i = INTEGER(islot); ans->p = INTEGER(GET_SLOT(x, Matrix_pSym)); /* dimensions and nzmax */ ans->nrow = dims[0]; ans->ncol = dims[1]; /* Allow for over-allocation of the i and x slots. Needed for * sparse X form in lme4. Right now it looks too difficult to * check for the length of the x slot, because of the xpt * utility, but the lengths of x and i should agree. */ ans->nzmax = LENGTH(islot); /* values depending on ctype */ ans->x = xpt(ctype, x); ans->stype = stype(ctype, x); ans->xtype = xtype(ctype); /* are the columns sorted (increasing row numbers) ?*/ ans->sorted = check_sorted_chm(ans); if (!(ans->sorted)) { /* sort columns */ if(sort_in_place) { if (!cholmod_sort(ans, &c)) error(_("in_place cholmod_sort returned an error code")); ans->sorted = 1; } else { CHM_SP tmp = cholmod_copy_sparse(ans, &c); if (!cholmod_sort(tmp, &c)) error(_("cholmod_sort returned an error code")); #ifdef DEBUG_Matrix /* This "triggers" exactly for return values of dtCMatrix_sparse_solve():*/ /* Don't want to translate this: want it report */ Rprintf("Note: as_cholmod_sparse() needed cholmod_sort()ing\n"); #endif chm2Ralloc(ans, tmp); cholmod_free_sparse(&tmp, &c); } } if (check_Udiag && ctype % 3 == 2 // triangular && (*diag_P(x) == 'U')) { /* diagU2N(.) "in place" : */ double one[] = {1, 0}; CHM_SP eye = cholmod_speye(ans->nrow, ans->ncol, ans->xtype, &c); CHM_SP tmp = cholmod_add(ans, eye, one, one, TRUE, TRUE, &c); #ifdef DEBUG_Matrix_verbose /* happens quite often, e.g. in ../tests/indexing.R : */ Rprintf("Note: as_cholmod_sparse(<ctype=%d>) - diagU2N\n", ctype); #endif chm2Ralloc(ans, tmp); cholmod_free_sparse(&tmp, &c); cholmod_free_sparse(&eye, &c); } /* else : * NOTE: if(*diag_P(x) == 'U'), the diagonal is lost (!); * ---- that may be ok, e.g. if we are just converting from/to Tsparse, * but is *not* at all ok, e.g. when used before matrix products */ return ans; }
/** * Copy the contents of a to an appropriate CsparseMatrix object and, * optionally, free a or free both a and its the pointers to its contents. * * @param a (cholmod_sparse) matrix to be converted * @param dofree 0 - don't free a; > 0 cholmod_free a; < 0 Free a * @param uploT 0 - not triangular; > 0 upper triangular; < 0 lower * @param Rkind - vector type to store for a->xtype == CHOLMOD_REAL, * 0 - REAL; 1 - LOGICAL [unused for other a->xtype] * @param diag character string suitable for the diag slot of a * triangular matrix (not accessed if uploT == 0). * @param dn either R_NilValue or an SEXP suitable for the Dimnames slot. * * @return SEXP containing a copy of a */ SEXP chm_sparse_to_SEXP(CHM_SP a, int dofree, int uploT, int Rkind, const char* diag, SEXP dn) { SEXP ans; char *cls = "";/* -Wall */ Rboolean longi = (a->itype) == CHOLMOD_LONG; int *dims, nnz, *ansp, *ansi; // if (longi) : SuiteSparse_long *ail = (SuiteSparse_long*)(a->i), *apl = (SuiteSparse_long*)(a->p); // else ((a->itype) == CHOLMOD_INT) : int *aii = (int*)(a->i), *api = (int*)(a->p); PROTECT(dn); /* dn is usually UNPROTECTed before the call */ /* ensure a is sorted and packed */ if (!a->sorted || !a->packed) longi ? cholmod_l_sort(a, &cl) : cholmod_sort(a, &c); /* determine the class of the result */ #define DOFREE_MAYBE \ if (dofree > 0) \ longi ? cholmod_l_free_sparse(&a, &cl) : cholmod_free_sparse(&a, &c); \ else if (dofree < 0) Free(a) switch(a->xtype) { case CHOLMOD_PATTERN: cls = uploT ? "ntCMatrix": ((a->stype) ? "nsCMatrix" : "ngCMatrix"); break; case CHOLMOD_REAL: switch(Rkind) { case 0: cls = uploT ? "dtCMatrix": ((a->stype) ? "dsCMatrix" : "dgCMatrix"); break; case 1: cls = uploT ? "ltCMatrix": ((a->stype) ? "lsCMatrix" : "lgCMatrix"); break; default: DOFREE_MAYBE; error(_("chm_sparse_to_SEXP(<real>, *): invalid 'Rkind' (real kind code)")); } break; case CHOLMOD_COMPLEX: cls = uploT ? "ztCMatrix": ((a->stype) ? "zsCMatrix" : "zgCMatrix"); break; default: DOFREE_MAYBE; error(_("unknown xtype in cholmod_sparse object")); } ans = PROTECT(NEW_OBJECT_OF_CLASS(cls)); /* allocate and copy common slots */ nnz = longi ? cholmod_l_nnz(a, &cl) : cholmod_nnz(a, &c); dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); dims[0] = a->nrow; dims[1] = a->ncol; ansp = INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, a->ncol + 1)); ansi = INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)); for (int j = 0; j <= a->ncol; j++) ansp[j] = longi ? (int)(apl[j]) : api[j]; for (int p = 0; p < nnz; p++) ansi[p] = longi ? (int)(ail[p]) : aii[p]; /* copy data slot if present */ if (a->xtype == CHOLMOD_REAL) { int i, *m_x; double *a_x = (double *) a->x; switch(Rkind) { case 0: Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), a_x, nnz); break; case 1: m_x = LOGICAL(ALLOC_SLOT(ans, Matrix_xSym, LGLSXP, nnz)); for (i=0; i < nnz; i++) m_x[i] = ISNAN(a_x[i]) ? NA_LOGICAL : (a_x[i] != 0); break; } } else if (a->xtype == CHOLMOD_COMPLEX) { DOFREE_MAYBE; error(_("complex sparse matrix code not yet written")); /* Memcpy(COMPLEX(ALLOC_SLOT(ans, Matrix_xSym, CPLXSXP, nnz)), */ /* (complex *) a->x, nnz); */ } if (uploT) { /* slots for triangularMatrix */ if (a->stype) error(_("Symmetric and triangular both set")); SET_SLOT(ans, Matrix_uploSym, mkString((uploT > 0) ? "U" : "L")); SET_SLOT(ans, Matrix_diagSym, mkString(diag)); } if (a->stype) /* slot for symmetricMatrix */ SET_SLOT(ans, Matrix_uploSym, mkString((a->stype > 0) ? "U" : "L")); DOFREE_MAYBE; if (dn != R_NilValue) SET_SLOT(ans, Matrix_DimNamesSym, duplicate(dn)); UNPROTECT(2); return ans; }
/** * Drop the (unit) diagonal entries from a cholmod_sparse matrix * * @param chx cholmod_sparse matrix. * Note that the matrix "slots" are modified _in place_ * @param uploT integer code (= +/- 1) indicating if chx is * upper (+1) or lower (-1) triangular * @param do_realloc Rboolean indicating, if a cholmod_sprealloc() should * finalize the procedure; not needed, e.g. when the * result is converted to a SEXP immediately afterwards. */ void chm_diagN2U(CHM_SP chx, int uploT, Rboolean do_realloc) { int i, n = chx->nrow, nnz = (int)cholmod_nnz(chx, &c), n_nnz = nnz - n, /* the new nnz : we will have removed n entries */ i_to = 0, i_from = 0; if(chx->ncol != n) error(_("chm_diagN2U(<non-square matrix>): nrow=%d, ncol=%d"), n, chx->ncol); if (!chx->sorted || !chx->packed) cholmod_sort(chx, &c); /* dimensions and nzmax */ #define _i(I) ( (int*) chx->i)[I] #define _x(I) ((double*) chx->x)[I] #define _p(I) ( (int*) chx->p)[I] /* work by copying from i_from to i_to ==> MUST i_to <= i_from */ if(uploT == 1) { /* "U" : upper triangular */ for(i = 0; i < n; i++) { /* looking at i-th column */ int j, n_i = _p(i+1) - _p(i); /* = #{entries} in this column */ /* 1) copy all but the last _above-diagonal_ column-entries: */ for(j = 1; j < n_i; j++, i_to++, i_from++) { _i(i_to) = _i(i_from); _x(i_to) = _x(i_from); } /* 2) drop the last column-entry == diagonal entry */ i_from++; } } else if(uploT == -1) { /* "L" : lower triangular */ for(i = 0; i < n; i++) { /* looking at i-th column */ int j, n_i = _p(i+1) - _p(i); /* = #{entries} in this column */ /* 1) drop the first column-entry == diagonal entry */ i_from++; /* 2) copy the other _below-diagonal_ column-entries: */ for(j = 1; j < n_i; j++, i_to++, i_from++) { _i(i_to) = _i(i_from); _x(i_to) = _x(i_from); } } } else { error(_("chm_diagN2U(x, uploT = %d): uploT should be +- 1"), uploT); } /* the column pointers are modified the same in both cases :*/ for(i=1; i <= n; i++) _p(i) -= i; #undef _i #undef _x #undef _p if(do_realloc) /* shorten (i- and x-slots from nnz to n_nnz */ cholmod_reallocate_sparse(n_nnz, chx, &c); return; }