SCM

SCM Repository

[matrix] View of /pkg/Matrix/src/Csparse.c
ViewVC logotype

View of /pkg/Matrix/src/Csparse.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2681 - (download) (as text) (annotate)
Mon Aug 1 20:43:16 2011 UTC (7 years, 8 months ago) by mmaechler
File size: 38129 byte(s)
Csparse_subassign() now also for dtCMatrix
			/* Sparse matrices in compressed column-oriented form */

#include <stdint.h> // C99 for int64_t
#include "Csparse.h"
#include "Tsparse.h"
#include "chm_common.h"

/** "Cheap" C version of  Csparse_validate() - *not* sorting : */
Rboolean isValid_Csparse(SEXP x)
{
    /* NB: we do *NOT* check a potential 'x' slot here, at all */
    SEXP pslot = GET_SLOT(x, Matrix_pSym),
	islot = GET_SLOT(x, Matrix_iSym);
    int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), j,
	nrow = dims[0],
	ncol = dims[1],
	*xp = INTEGER(pslot),
	*xi = INTEGER(islot);

    if (length(pslot) != dims[1] + 1)
	return FALSE;
    if (xp[0] != 0)
	return FALSE;
    if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
	return FALSE;
    for (j = 0; j < xp[ncol]; j++) {
	if (xi[j] < 0 || xi[j] >= nrow)
	    return FALSE;
    }
    for (j = 0; j < ncol; j++) {
	if (xp[j] > xp[j + 1])
	    return FALSE;
    }
    return TRUE;
}

SEXP Csparse_validate(SEXP x) {
    return Csparse_validate_(x, FALSE);
}

SEXP Csparse_validate2(SEXP x, SEXP maybe_modify) {
    return Csparse_validate_(x, asLogical(maybe_modify));
}

SEXP Csparse_validate_(SEXP x, Rboolean maybe_modify)
{
    /* NB: we do *NOT* check a potential 'x' slot here, at all */
    SEXP pslot = GET_SLOT(x, Matrix_pSym),
	islot = GET_SLOT(x, Matrix_iSym);
    Rboolean sorted, strictly;
    int j, k,
	*dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
	nrow = dims[0],
	ncol = dims[1],
	*xp = INTEGER(pslot),
	*xi = INTEGER(islot);

    if (length(pslot) != dims[1] + 1)
	return mkString(_("slot p must have length = ncol(.) + 1"));
    if (xp[0] != 0)
	return mkString(_("first element of slot p must be zero"));
    if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
	return
	    mkString(_("last element of slot p must match length of slots i and x"));
    for (j = 0; j < xp[ncol]; j++) {
	if (xi[j] < 0 || xi[j] >= nrow)
	    return mkString(_("all row indices must be between 0 and nrow-1"));
    }
    sorted = TRUE; strictly = TRUE;
    for (j = 0; j < ncol; j++) {
	if (xp[j] > xp[j + 1])
	    return mkString(_("slot p must be non-decreasing"));
	if(sorted) /* only act if >= 2 entries in column j : */
	    for (k = xp[j] + 1; k < xp[j + 1]; k++) {
		if (xi[k] < xi[k - 1])
		    sorted = FALSE;
		else if (xi[k] == xi[k - 1])
		    strictly = FALSE;
	    }
    }
    if (!sorted) {
	if(maybe_modify) {
	    CHM_SP chx = (CHM_SP) alloca(sizeof(cholmod_sparse));
	    R_CheckStack();
	    as_cholmod_sparse(chx, x, FALSE, TRUE);/*-> cholmod_l_sort() ! */
	    /* as chx = AS_CHM_SP__(x)  but  ^^^^ sorting x in_place !!! */

	    /* Now re-check that row indices are *strictly* increasing
	     * (and not just increasing) within each column : */
	    for (j = 0; j < ncol; j++) {
		for (k = xp[j] + 1; k < xp[j + 1]; k++)
		    if (xi[k] == xi[k - 1])
			return mkString(_("slot i is not *strictly* increasing inside a column (even after cholmod_l_sort)"));
	    }
	} else { /* no modifying sorting : */
	    return mkString(_("row indices are not sorted within columns"));
	}
    } else if(!strictly) {  /* sorted, but not strictly */
	return mkString(_("slot i is not *strictly* increasing inside a column"));
    }
    return ScalarLogical(1);
}

SEXP Rsparse_validate(SEXP x)
{
    /* NB: we do *NOT* check a potential 'x' slot here, at all */
    SEXP pslot = GET_SLOT(x, Matrix_pSym),
	jslot = GET_SLOT(x, Matrix_jSym);
    Rboolean sorted, strictly;
    int i, k,
	*dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
	nrow = dims[0],
	ncol = dims[1],
	*xp = INTEGER(pslot),
	*xj = INTEGER(jslot);

    if (length(pslot) != dims[0] + 1)
	return mkString(_("slot p must have length = nrow(.) + 1"));
    if (xp[0] != 0)
	return mkString(_("first element of slot p must be zero"));
    if (length(jslot) < xp[nrow]) /* allow larger slots from over-allocation!*/
	return
	    mkString(_("last element of slot p must match length of slots j and x"));
    for (i = 0; i < length(jslot); i++) {
	if (xj[i] < 0 || xj[i] >= ncol)
	    return mkString(_("all column indices must be between 0 and ncol-1"));
    }
    sorted = TRUE; strictly = TRUE;
    for (i = 0; i < nrow; i++) {
	if (xp[i] > xp[i+1])
	    return mkString(_("slot p must be non-decreasing"));
	if(sorted)
	    for (k = xp[i] + 1; k < xp[i + 1]; k++) {
		if (xj[k] < xj[k - 1])
		    sorted = FALSE;
		else if (xj[k] == xj[k - 1])
		    strictly = FALSE;
	    }
    }
    if (!sorted)
	/* cannot easily use cholmod_sort(.) ... -> "error out" :*/
	return mkString(_("slot j is not increasing inside a column"));
    else if(!strictly) /* sorted, but not strictly */
	return mkString(_("slot j is not *strictly* increasing inside a column"));

    return ScalarLogical(1);
}


/* Called from ../R/Csparse.R : */
/* Can only return [dln]geMatrix (no symm/triang);
 * FIXME: replace by non-CHOLMOD code ! */
SEXP Csparse_to_dense(SEXP x)
{
    CHM_SP chxs = AS_CHM_SP__(x);
    /* This loses the symmetry property, since cholmod_dense has none,
     * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices
     * to numeric (CHOLMOD_REAL) ones : */
    CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c);
    int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x);
    R_CheckStack();

    return chm_dense_to_SEXP(chxd, 1, Rkind, GET_SLOT(x, Matrix_DimNamesSym));
}

// FIXME: do not go via CHM (should not be too hard, to just *drop* the x-slot, right?
SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri)
{
    CHM_SP chxs = AS_CHM_SP__(x);
    CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c);
    int tr = asLogical(tri);
    R_CheckStack();

    return chm_sparse_to_SEXP(chxcp, 1/*do_free*/,
			      tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
			      0, tr ? diag_P(x) : "",
			      GET_SLOT(x, Matrix_DimNamesSym));
}

// n.CMatrix --> [dli].CMatrix  (not going through CHM!)
SEXP nz_pattern_to_Csparse(SEXP x, SEXP res_kind)
{
    return nz2Csparse(x, asInteger(res_kind));
}
// n.CMatrix --> [dli].CMatrix  (not going through CHM!)
SEXP nz2Csparse(SEXP x, enum x_slot_kind r_kind)
{
    const char *cl_x = class_P(x);
    if(cl_x[0] != 'n') error(_("not a 'n.CMatrix'"));
    if(cl_x[2] != 'C') error(_("not a CsparseMatrix"));
    int nnz = LENGTH(GET_SLOT(x, Matrix_iSym));
    SEXP ans;
    char *ncl = strdup(cl_x);
    double *dx_x; int *ix_x;
    ncl[0] = (r_kind == x_double ? 'd' :
	      (r_kind == x_logical ? 'l' :
	       /* else (for now):  r_kind == x_integer : */ 'i'));
    PROTECT(ans = NEW_OBJECT(MAKE_CLASS(ncl)));
    // create a correct 'x' slot:
    switch(r_kind) {
	int i;
    case x_double: // 'd'
	dx_x = REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz));
	for (i=0; i < nnz; i++) dx_x[i] = 1.;
	break;
    case x_logical: // 'l'
	ix_x = LOGICAL(ALLOC_SLOT(ans, Matrix_xSym, LGLSXP, nnz));
	for (i=0; i < nnz; i++) ix_x[i] = TRUE;
	break;
    case x_integer: // 'i'
	ix_x = INTEGER(ALLOC_SLOT(ans, Matrix_xSym, INTSXP, nnz));
	for (i=0; i < nnz; i++) ix_x[i] = 1;
	break;

    default:
	error(_("nz2Csparse(): invalid/non-implemented r_kind = %d"),
	      r_kind);
    }

    // now copy all other slots :
    slot_dup(ans, x, Matrix_iSym);
    slot_dup(ans, x, Matrix_pSym);
    slot_dup(ans, x, Matrix_DimSym);
    slot_dup(ans, x, Matrix_DimNamesSym);
    if(ncl[1] != 'g') { // symmetric or triangular ...
	slot_dup_if_has(ans, x, Matrix_uploSym);
	slot_dup_if_has(ans, x, Matrix_diagSym);
    }
    UNPROTECT(1);
    return ans;
}

SEXP Csparse_to_matrix(SEXP x)
{
    return chm_dense_to_matrix(cholmod_sparse_to_dense(AS_CHM_SP__(x), &c),
			       1 /*do_free*/, GET_SLOT(x, Matrix_DimNamesSym));
}

SEXP Csparse_to_Tsparse(SEXP x, SEXP tri)
{
    CHM_SP chxs = AS_CHM_SP__(x);
    CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c);
    int tr = asLogical(tri);
    int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    R_CheckStack();

    return chm_triplet_to_SEXP(chxt, 1,
			       tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
			       Rkind, tr ? diag_P(x) : "",
			       GET_SLOT(x, Matrix_DimNamesSym));
}

/* this used to be called  sCMatrix_to_gCMatrix(..)   [in ./dsCMatrix.c ]: */
SEXP Csparse_symmetric_to_general(SEXP x)
{
    CHM_SP chx = AS_CHM_SP__(x), chgx;
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    R_CheckStack();

    if (!(chx->stype))
	error(_("Nonsymmetric matrix in Csparse_symmetric_to_general"));
    chgx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
    /* xtype: pattern, "real", complex or .. */
    return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
			      GET_SLOT(x, Matrix_DimNamesSym));
}

SEXP Csparse_general_to_symmetric(SEXP x, SEXP uplo)
{
    CHM_SP chx = AS_CHM_SP__(x), chgx;
    int uploT = (*CHAR(STRING_ELT(uplo,0)) == 'U') ? 1 : -1;
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    R_CheckStack();

    chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c);
    /* xtype: pattern, "real", complex or .. */
    return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
			      GET_SLOT(x, Matrix_DimNamesSym));
}

SEXP Csparse_transpose(SEXP x, SEXP tri)
{
    /* TODO: lgCMatrix & igC* currently go via double prec. cholmod -
     *       since cholmod (& cs) lacks sparse 'int' matrices */
    CHM_SP chx = AS_CHM_SP__(x);
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c);
    SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp;
    int tr = asLogical(tri);
    R_CheckStack();

    tmp = VECTOR_ELT(dn, 0);	/* swap the dimnames */
    SET_VECTOR_ELT(dn, 0, VECTOR_ELT(dn, 1));
    SET_VECTOR_ELT(dn, 1, tmp);
    UNPROTECT(1);
    return chm_sparse_to_SEXP(chxt, 1, /* SWAP 'uplo' for triangular */
			      tr ? ((*uplo_P(x) == 'U') ? -1 : 1) : 0,
			      Rkind, tr ? diag_P(x) : "", dn);
}

SEXP Csparse_Csparse_prod(SEXP a, SEXP b)
{
    CHM_SP
	cha = AS_CHM_SP(a),
	chb = AS_CHM_SP(b),
	chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0,
			       /* values:= is_numeric (T/F) */ cha->xtype > 0,
			       /*out sorted:*/ 1, &c);
    const char *cl_a = class_P(a), *cl_b = class_P(b);
    char diag[] = {'\0', '\0'};
    int uploT = 0;
    SEXP dn = PROTECT(allocVector(VECSXP, 2));
    R_CheckStack();

#ifdef DEBUG_Matrix_verbose
    Rprintf("DBG Csparse_C*_prod(%s, %s)\n", cl_a, cl_b);
#endif

    /* Preserve triangularity and even unit-triangularity if appropriate.
     * Note that in that case, the multiplication itself should happen
     * faster.  But there's no support for that in CHOLMOD */

    /* UGLY hack -- rather should have (fast!) C-level version of
     *       is(a, "triangularMatrix") etc */
    if (cl_a[1] == 't' && cl_b[1] == 't')
	/* FIXME: fails for "Cholesky","BunchKaufmann"..*/
	if(*uplo_P(a) == *uplo_P(b)) { /* both upper, or both lower tri. */
	    uploT = (*uplo_P(a) == 'U') ? 1 : -1;
	    if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
		/* "remove the diagonal entries": */
		chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
		diag[0]= 'U';
	    }
	    else diag[0]= 'N';
	}
    SET_VECTOR_ELT(dn, 0,	/* establish dimnames */
		   duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
    SET_VECTOR_ELT(dn, 1,
		   duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), 1)));
    UNPROTECT(1);
    return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
}

SEXP Csparse_Csparse_crossprod(SEXP a, SEXP b, SEXP trans)
{
    int tr = asLogical(trans);
    CHM_SP
	cha = AS_CHM_SP(a),
	chb = AS_CHM_SP(b),
	chTr, chc;
    const char *cl_a = class_P(a), *cl_b = class_P(b);
    char diag[] = {'\0', '\0'};
    int uploT = 0;
    SEXP dn = PROTECT(allocVector(VECSXP, 2));
    R_CheckStack();

    chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c);
    chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb,
			 /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c);
    cholmod_free_sparse(&chTr, &c);

    /* Preserve triangularity and unit-triangularity if appropriate;
     * see Csparse_Csparse_prod() for comments */
    if (cl_a[1] == 't' && cl_b[1] == 't')
	if(*uplo_P(a) != *uplo_P(b)) { /* one 'U', the other 'L' */
	    uploT = (*uplo_P(b) == 'U') ? 1 : -1;
	    if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
		chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
		diag[0]= 'U';
	    }
	    else diag[0]= 'N';
	}
    SET_VECTOR_ELT(dn, 0,	/* establish dimnames */
		   duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), (tr) ? 0 : 1)));
    SET_VECTOR_ELT(dn, 1,
		   duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), (tr) ? 0 : 1)));
    UNPROTECT(1);
    return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
}

SEXP Csparse_dense_prod(SEXP a, SEXP b)
{
    CHM_SP cha = AS_CHM_SP(a);
    SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
    CHM_DN chb = AS_CHM_DN(b_M);
    CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow,
					chb->xtype, &c);
    SEXP dn = PROTECT(allocVector(VECSXP, 2));
    double one[] = {1,0}, zero[] = {0,0};
    int nprot = 2;
    R_CheckStack();
    /* Tim Davis, please FIXME:  currently (2010-11) *fails* when  a  is a pattern matrix:*/
    if(cha->xtype == CHOLMOD_PATTERN) {
	/* warning(_("Csparse_dense_prod(): cholmod_sdmult() not yet implemented for pattern./ ngCMatrix" */
	/* 	  " --> slightly inefficient coercion")); */

	// This *fails* to produce a CHOLMOD_REAL ..
	// CHM_SP chd = cholmod_l_copy(cha, cha->stype, CHOLMOD_REAL, &c);
	// --> use our Matrix-classes
	SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
	cha = AS_CHM_SP(da);
    }
    cholmod_sdmult(cha, 0, one, zero, chb, chc, &c);
    SET_VECTOR_ELT(dn, 0,	/* establish dimnames */
		   duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
    SET_VECTOR_ELT(dn, 1,
		   duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
    UNPROTECT(nprot);
    return chm_dense_to_SEXP(chc, 1, 0, dn);
}

SEXP Csparse_dense_crossprod(SEXP a, SEXP b)
{
    CHM_SP cha = AS_CHM_SP(a);
    SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
    CHM_DN chb = AS_CHM_DN(b_M);
    CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol,
					chb->xtype, &c);
    SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2;
    double one[] = {1,0}, zero[] = {0,0};
    R_CheckStack();
    // -- see Csparse_dense_prod() above :
    if(cha->xtype == CHOLMOD_PATTERN) {
	SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
	cha = AS_CHM_SP(da);
    }
    cholmod_sdmult(cha, 1, one, zero, chb, chc, &c);
    SET_VECTOR_ELT(dn, 0,	/* establish dimnames */
		   duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1)));
    SET_VECTOR_ELT(dn, 1,
		   duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
    UNPROTECT(nprot);
    return chm_dense_to_SEXP(chc, 1, 0, dn);
}

/* Computes   x'x  or  x x' -- *also* for Tsparse (triplet = TRUE)
   see Csparse_Csparse_crossprod above for  x'y and x y' */
SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet)
{
    int trip = asLogical(triplet),
	tr   = asLogical(trans); /* gets reversed because _aat is tcrossprod */
#ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
    CHM_TR cht = trip ? AS_CHM_TR(x) : (CHM_TR) NULL;
#else /* workaround needed:*/
    SEXP xx = PROTECT(Tsparse_diagU2N(x));
    CHM_TR cht = trip ? AS_CHM_TR__(xx) : (CHM_TR) NULL;
#endif
    CHM_SP chcp, chxt,
	chx = (trip ?
	       cholmod_triplet_to_sparse(cht, cht->nnz, &c) :
	       AS_CHM_SP(x));
    SEXP dn = PROTECT(allocVector(VECSXP, 2));
    R_CheckStack();

    if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c);
    chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c);
    if(!chcp) {
	UNPROTECT(1);
	error(_("Csparse_crossprod(): error return from cholmod_aat()"));
    }
    cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c);
    chcp->stype = 1;
    if (trip) cholmod_free_sparse(&chx, &c);
    if (!tr) cholmod_free_sparse(&chxt, &c);
    SET_VECTOR_ELT(dn, 0,	/* establish dimnames */
		   duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym),
					(tr) ? 0 : 1)));
    SET_VECTOR_ELT(dn, 1, duplicate(VECTOR_ELT(dn, 0)));
#ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
    UNPROTECT(1);
#else
    UNPROTECT(2);
#endif
    return chm_sparse_to_SEXP(chcp, 1, 0, 0, "", dn);
}

/* Csparse_drop(x, tol):  drop entries with absolute value < tol, i.e,
*  at least all "explicit" zeros */
SEXP Csparse_drop(SEXP x, SEXP tol)
{
    const char *cl = class_P(x);
    /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
    int tr = (cl[1] == 't');
    CHM_SP chx = AS_CHM_SP__(x);
    CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c);
    double dtol = asReal(tol);
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    R_CheckStack();

    if(!cholmod_drop(dtol, ans, &c))
	error(_("cholmod_drop() failed"));
    return chm_sparse_to_SEXP(ans, 1,
			      tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
			      Rkind, tr ? diag_P(x) : "",
			      GET_SLOT(x, Matrix_DimNamesSym));
}

SEXP Csparse_horzcat(SEXP x, SEXP y)
{
    CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
    int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
	Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
	Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
    R_CheckStack();

    /* TODO: currently drops dimnames - and we fix at R level */
    return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c),
			      1, 0, Rkind, "", R_NilValue);
}

SEXP Csparse_vertcat(SEXP x, SEXP y)
{
    CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
    int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
	Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
	Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
    R_CheckStack();

    /* TODO: currently drops dimnames - and we fix at R level */
    return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c),
			      1, 0, Rkind, "", R_NilValue);
}

SEXP Csparse_band(SEXP x, SEXP k1, SEXP k2)
{
    CHM_SP chx = AS_CHM_SP__(x);
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c);
    R_CheckStack();

    return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "",
			      GET_SLOT(x, Matrix_DimNamesSym));
}

SEXP Csparse_diagU2N(SEXP x)
{
    const char *cl = class_P(x);
    /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
    if (cl[1] != 't' || *diag_P(x) != 'U') {
	/* "trivially fast" when not triangular (<==> no 'diag' slot),
	   or not *unit* triangular */
	return (x);
    }
    else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */
	CHM_SP chx = AS_CHM_SP__(x);
	CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c);
	double one[] = {1, 0};
	CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c);
	int uploT = (*uplo_P(x) == 'U') ? 1 : -1;
	int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;

	R_CheckStack();
	cholmod_free_sparse(&eye, &c);
	return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N",
				  GET_SLOT(x, Matrix_DimNamesSym));
    }
}

SEXP Csparse_diagN2U(SEXP x)
{
    const char *cl = class_P(x);
    /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
    if (cl[1] != 't' || *diag_P(x) != 'N') {
	/* "trivially fast" when not triangular (<==> no 'diag' slot),
	   or already *unit* triangular */
	return (x);
    }
    else { /* triangular with diag='N'): now drop the diagonal */
	/* duplicate, since chx will be modified: */
	CHM_SP chx = AS_CHM_SP__(duplicate(x));
	int uploT = (*uplo_P(x) == 'U') ? 1 : -1,
	    Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
	R_CheckStack();

	chm_diagN2U(chx, uploT, /* do_realloc */ FALSE);

	return chm_sparse_to_SEXP(chx, /*dofree*/ 0/* or 1 ?? */,
				  uploT, Rkind, "U",
				  GET_SLOT(x, Matrix_DimNamesSym));
    }
}

/**
 * "Indexing" aka subsetting : Compute  x[i,j], also for vectors i and j
 * Working via CHOLMOD_submatrix, see ./CHOLMOD/MatrixOps/cholmod_submatrix.c
 * @param x CsparseMatrix
 * @param i row     indices (0-origin), or NULL (R's)
 * @param j columns indices (0-origin), or NULL
 *
 * @return x[i,j]  still CsparseMatrix --- currently, this loses dimnames
 */
SEXP Csparse_submatrix(SEXP x, SEXP i, SEXP j)
{
    CHM_SP chx = AS_CHM_SP(x); /* << does diagU2N() when needed */
    int rsize = (isNull(i)) ? -1 : LENGTH(i),
	csize = (isNull(j)) ? -1 : LENGTH(j);
    int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
    R_CheckStack();

    if (rsize >= 0 && !isInteger(i))
	error(_("Index i must be NULL or integer"));
    if (csize >= 0 && !isInteger(j))
	error(_("Index j must be NULL or integer"));

    if (chx->stype) /* symmetricMatrix */
	/* for now, cholmod_submatrix() only accepts "generalMatrix" */
	chx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);

    return chm_sparse_to_SEXP(cholmod_submatrix(chx,
				(rsize < 0) ? NULL : INTEGER(i), rsize,
			        (csize < 0) ? NULL : INTEGER(j), csize,
						  TRUE, TRUE, &c),
			      1, 0, Rkind, "",
			      /* FIXME: drops dimnames */ R_NilValue);
}

/**
 * Subassignment:  x[i,j]  <- value
 *
 * @param x
 * @param i_ integer row    index 0-origin vector (as returned from R .ind.prep2())
 * @param j_ integer column index 0-origin vector
 * @param value currently must be a dsparseVector {which is recycled if needed}
 *
 * @return a Csparse matrix like x, but with the values replaced
 */
SEXP Csparse_subassign(SEXP x, SEXP i_, SEXP j_, SEXP value)
{
    // TODO: for other classes consider using a trick as  RallocedReal() in ./chm_common.c
    static const char
	*valid_cM [] = { // the only ones, for "the moment". FIXME: extend (!)
	"dgCMatrix",// 0
	"dtCMatrix",// 1
	""},
	// value: assume a  "dsparseVector" for now -- slots: (i, length, x)
	*valid_spv[] = {"dsparseVector",
			""};

    int ctype_x = Matrix_check_class_etc(x, valid_cM),
	ctype_v = Matrix_check_class_etc(value, valid_spv);
    if (ctype_x < 0)
	error(_("invalid class of 'x' in Csparse_subassign()"));
    if (ctype_v < 0)
	error(_("invalid class of 'value' in Csparse_subassign()"));

    SEXP
	islot   = GET_SLOT(x, Matrix_iSym),
	dimslot = GET_SLOT(x, Matrix_DimSym),
	i_cp = PROTECT(coerceVector(i_, INTSXP)),
	j_cp = PROTECT(coerceVector(j_, INTSXP));
	// for d.CMatrix and l.CMatrix  but not n.CMatrix:

    int *dims = INTEGER(dimslot),
	ncol = dims[1],	/* nrow = dims[0], */
	*i = INTEGER(i_cp), len_i = LENGTH(i_cp),
	*j = INTEGER(j_cp), len_j = LENGTH(j_cp),
	k,
	nnz_x = LENGTH(islot);
    int nnz = nnz_x;

#define MATRIX_SUBASSIGN_VERBOSE
// Temporary hack for debugging --- remove eventually -- FIXME
#ifdef MATRIX_SUBASSIGN_VERBOSE
    Rboolean verbose = i[0] < 0;
    if(verbose) i[0] = -i[0];
#endif

    SEXP val_i_slot;
    PROTECT(val_i_slot = coerceVector(GET_SLOT(value, Matrix_iSym), REALSXP));
    double *val_i = REAL(val_i_slot);
    int nnz_val =  LENGTH(GET_SLOT(value, Matrix_iSym));
    // for dsparseVector only:
    double *val_x =   REAL (GET_SLOT(value, Matrix_xSym));
    int64_t len_val = (int64_t) asReal(GET_SLOT(value, Matrix_lengthSym));
    /* llen_i = (int64_t) len_i; */

    double z_ans = 0.;

    SEXP ans;
    /* Instead of simple "duplicate": PROTECT(ans = duplicate(x)) , build up: */
    // Assuming that ans will have the same basic Matrix type as x :
    ans = PROTECT(NEW_OBJECT(MAKE_CLASS(valid_cM[ctype_x])));
    SET_SLOT(ans, Matrix_DimSym,      duplicate(dimslot));
    SET_SLOT(ans, Matrix_DimNamesSym, duplicate(GET_SLOT(x, Matrix_DimNamesSym)));
    SET_SLOT(ans, Matrix_pSym,        duplicate(GET_SLOT(x, Matrix_pSym)));
    SEXP r_pslot = GET_SLOT(ans, Matrix_pSym);
    // and assign the i- and x- slots at the end, as they are potentially modified
    // not just in content, but also in their *length*
    int	*rp = INTEGER(r_pslot),
	*ri = Calloc(nnz_x, int);       // to contain the final i - slot
    // for d.CMatrix only:
    double *rx = Calloc(nnz_x, double); // to contain the final x - slot
    Memcpy(ri, INTEGER(islot), nnz_x);
    Memcpy(rx, REAL(GET_SLOT(x, Matrix_xSym)), nnz_x);
    // NB:  nnz_x : will always be the "current allocated length" of (i, x) slots
    // --   nnz   : the current *used* length; always   nnz <= nnz_x

    int jj, j_val = 0; // in "running" conceptionally through all value[i+ jj*len_i]
    // values, we are "below"/"before" the (j_val)-th non-zero one.
    // e.g. if value = (0,0,...,0), have nnz_val == 0, j_val must remain == 0
    int64_t ii_val;// == "running" index (i + jj*len_i) % len_val for value[]
    for(jj = 0, ii_val=0; jj < len_j; jj++) {
	int j__ = j[jj];
	/* int64_t j_l = jj * llen_i; */
	R_CheckUserInterrupt();
	for(int ii = 0; ii < len_i; ii++, ii_val++) {
	    int i__ = i[ii], p1, p2;
	    if(nnz_val && ii_val >= len_val) { // "recycle" indexing into value[]
		ii_val -= len_val; // = (ii + jj*len_i) % len_val
		j_val = 0;
	    }
	    int64_t ii_v1;//= ii_val + 1;
	    double v, /* := value[(ii + j_l) % len_val]
			 = dsparseVector_sub((ii + j_l) % len_val,
			                     nnz_val, val_i, val_x, len_val)
		      */
		M_ij;
	    int ind;
	    Rboolean have_entry = FALSE;

	    // note that rp[]'s may have *changed* even when 'j' remained!
	    // "FIXME": do this only *when* rp[] has changed
	    p1 = rp[j__], p2 = rp[j__ + 1];

	    // v :=  value[(ii + j_l) % len_val] = value[ii_val]
	    v = z_ans;
	    if(j_val < nnz_val) { // maybe find v := non-zero value[ii_val]
		ii_v1 = ii_val + 1;
		if(ii_v1 < val_i[j_val]) { // typical case: are still in zero-stretch
		    v = z_ans; // v = 0
		} else if(ii_v1 == val_i[j_val]) { // have a match
		    v = val_x[j_val];
		    j_val++;// from now on, look at the next non-zero entry
		} else { //  ii_v1 > val_i[j_val]
		    REprintf("programming thinko in Csparse_subassign(*, i=%d,j=%d): ii_v=%d, v@i[j_val=%ld]=%g\n",
			     i__,j__, ii_v1, j_val, val_i[j_val]);
		    j_val++;// from now on, look at the next non-zero entry
		}
	    }
	    // --------------- M_ij := getM(i., j.) --------------------------------
	    M_ij = z_ans; // as in  ./t_sparseVector.c
	    for(ind = p1; ind < p2; ind++) {
		if(ri[ind] >= i__) {
		    if(ri[ind] == i__) {
			M_ij = rx[ind];
#ifdef MATRIX_SUBASSIGN_VERBOSE
			if(verbose) REprintf("have entry x[%d, %d] = %g\n",
					     i__, j__, M_ij);
#endif
			have_entry = TRUE;
		    } else { // ri[ind] > i__
#ifdef MATRIX_SUBASSIGN_VERBOSE
			if(verbose)
			    REprintf("@i > i__ = %d --> ind-- = %d\n", i__, ind);
#endif
		    }
		    break;
		}
	    }

	    //-- R:  if(getM(i., j.) != (v <- getV(ii, jj)))

	    if(M_ij != v) { // contents differ ==> value needs to be changed
#ifdef MATRIX_SUBASSIGN_VERBOSE
		if(verbose)
		    REprintf("setting x[%d, %d] <- %g", i__,j__, v);
#endif
		// (otherwise: nothing to do):
		// setM(i__, j__, v)
		// ----------------------------------------------------------

		// Case I --------------------------------------------
/* 		if(v == z_ans) { // remove x[i, j] = M_ij  which we know is *non*-zero */
//-------- Better (memory-management) *NOT* to remove, but rather at the very end
//         currently using  drop0() in R code
/* 		    // we know : have_entry = TRUE ; */
/* 		    //  ri[ind] == i__; M_ij = rx[ind]; */
/* #ifdef MATRIX_SUBASSIGN_VERBOSE */
/* 		    if(verbose) */
/* 		    	REprintf(" rm ind=%d\n", ind); */
/* #endif */
/* 		    // remove the 'ind'-th element from x@i and x@x : */
/* 		    nnz-- ; */
/* 		    for(k=ind; k < nnz; k++) { */
/* 			ri[k] = ri[k+1]; */
/* 			rx[k] = rx[k+1]; */
/* 		    } */
/* 		    for(k=j__ + 1; k <= ncol; k++) { */
/* 			rp[k] = rp[k] - 1; */
/* 		    } */
/* 		} */
/* 		else  */
  	        if(have_entry) {
		    // Case II ----- replace (non-empty) x[i,j] by v -------
#ifdef MATRIX_SUBASSIGN_VERBOSE
		    if(verbose)
		    	REprintf(" repl.  ind=%d\n", ind);
#endif
		    rx[ind] = v;
		} else {
		    // Case III ---- v != 0 : insert v into "empty" x[i,j] ----

		    // extend the  i  and  x  slot by one entry : ---------------------

		    if(nnz+1 > nnz_x) { // need to reallocate:
#ifdef MATRIX_SUBASSIGN_VERBOSE
			if(verbose) REprintf(" Realloc()ing: nnz_x=%d", nnz_x);
#endif
			// do it "only" 1x,..4x at the very most increasing by the
			// nnz-length of "value":
			nnz_x += (1 + nnz_val / 4);
#ifdef MATRIX_SUBASSIGN_VERBOSE
			if(verbose) REprintf("(nnz_v=%d) --> %d ", nnz_val, nnz_x);
#endif
			// C doc on realloc() says that the old content is *preserve*d
			ri = Realloc(ri, nnz_x, int);
			rx = Realloc(rx, nnz_x, double);
		    }

		    // 3) fill them ...

		    int i1 = ind;
#ifdef MATRIX_SUBASSIGN_VERBOSE
		    if(verbose)
		    	REprintf(" INSERT p12=(%d,%d) -> ind=%d -> i1 = %d\n",
		    		 p1,p2, ind, i1);
#endif

		    // shift the "upper values" *before* the insertion:
		    for(int l = nnz-1; l >= i1; l--) {
			ri[l+1] = ri[l];
			rx[l+1] = rx[l];
		    }
		    ri[i1] = i__;
		    rx[i1] = v;
		    nnz++;

		    // the columns j "right" of the current one :
		    for(k=j__ + 1; k <= ncol; k++)
			rp[k]++;
		}
	    }
#ifdef MATRIX_SUBASSIGN_VERBOSE
	    else if(verbose) REprintf("M_ij == v = %g\n", v);
#endif
	}// for( ii )
    }// for( jj )

    if(ctype_x == 1) { // triangularMatrix: copy the 'diag' and 'uplo' slots
	SET_SLOT(ans, Matrix_uploSym, duplicate(GET_SLOT(x, Matrix_uploSym)));
	SET_SLOT(ans, Matrix_diagSym, duplicate(GET_SLOT(x, Matrix_diagSym)));
    }
    // now assign the i- and x- slots,  free memory and return :
    Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym,  INTSXP, nnz)), ri, nnz);
    Memcpy(   REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), rx, nnz);
    Free(ri);
    Free(rx);
    UNPROTECT(4);
    return ans;
}

SEXP Csparse_MatrixMarket(SEXP x, SEXP fname)
{
    FILE *f = fopen(CHAR(asChar(fname)), "w");

    if (!f)
	error(_("failure to open file \"%s\" for writing"),
	      CHAR(asChar(fname)));
    if (!cholmod_write_sparse(f, AS_CHM_SP(x),
			      (CHM_SP)NULL, (char*) NULL, &c))
	error(_("cholmod_write_sparse returned error code"));
    fclose(f);
    return R_NilValue;
}


/**
 * Extract the diagonal entries from *triangular* Csparse matrix  __or__ a
 * cholmod_sparse factor (LDL = TRUE).
 *
 * @param n  dimension of the matrix.
 * @param x_p  'p' (column pointer) slot contents
 * @param x_x  'x' (non-zero entries) slot contents
 * @param perm 'perm' (= permutation vector) slot contents; only used for "diagBack"
 * @param resultKind a (SEXP) string indicating which kind of result is desired.
 *
 * @return  a SEXP, either a (double) number or a length n-vector of diagonal entries
 */
SEXP diag_tC_ptr(int n, int *x_p, double *x_x, int *perm, SEXP resultKind)
/*                                ^^^^^^ FIXME[Generalize] to int / ... */
{
    const char* res_ch = CHAR(STRING_ELT(resultKind,0));
    enum diag_kind { diag, diag_backpermuted, trace, prod, sum_log
    } res_kind = ((!strcmp(res_ch, "trace")) ? trace :
		  ((!strcmp(res_ch, "sumLog")) ? sum_log :
		   ((!strcmp(res_ch, "prod")) ? prod :
		    ((!strcmp(res_ch, "diag")) ? diag :
		     ((!strcmp(res_ch, "diagBack")) ? diag_backpermuted :
		      -1)))));
    int i, n_x, i_from = 0;
    SEXP ans = PROTECT(allocVector(REALSXP,
/*                                 ^^^^  FIXME[Generalize] */
				   (res_kind == diag ||
				    res_kind == diag_backpermuted) ? n : 1));
    double *v = REAL(ans);
/*  ^^^^^^      ^^^^  FIXME[Generalize] */

#define for_DIAG(v_ASSIGN)						\
    for(i = 0; i < n; i++, i_from += n_x) {				\
	/* looking at i-th column */					\
	n_x = x_p[i+1] - x_p[i];/* #{entries} in this column */	\
	v_ASSIGN;							\
    }

    /* NOTA BENE: we assume  -- uplo = "L" i.e. lower triangular matrix
     *            for uplo = "U" (makes sense with a "dtCMatrix" !),
     *            should use  x_x[i_from + (nx - 1)] instead of x_x[i_from],
     *            where nx = (x_p[i+1] - x_p[i])
     */

    switch(res_kind) {
    case trace:
	v[0] = 0.;
	for_DIAG(v[0] += x_x[i_from]);
	break;

    case sum_log:
	v[0] = 0.;
	for_DIAG(v[0] += log(x_x[i_from]));
	break;

    case prod:
	v[0] = 1.;
	for_DIAG(v[0] *= x_x[i_from]);
	break;

    case diag:
	for_DIAG(v[i] = x_x[i_from]);
	break;

    case diag_backpermuted:
	for_DIAG(v[i] = x_x[i_from]);

	warning(_("resultKind = 'diagBack' (back-permuted) is experimental"));
	/* now back_permute : */
	for(i = 0; i < n; i++) {
	    double tmp = v[i]; v[i] = v[perm[i]]; v[perm[i]] = tmp;
	    /*^^^^ FIXME[Generalize] */
	}
	break;

    default: /* -1 from above */
	error(_("diag_tC(): invalid 'resultKind'"));
	/* Wall: */ ans = R_NilValue; v = REAL(ans);
    }

    UNPROTECT(1);
    return ans;
}

/**
 * Extract the diagonal entries from *triangular* Csparse matrix  __or__ a
 * cholmod_sparse factor (LDL = TRUE).
 *
 * @param pslot  'p' (column pointer)   slot of Csparse matrix/factor
 * @param xslot  'x' (non-zero entries) slot of Csparse matrix/factor
 * @param perm_slot  'perm' (= permutation vector) slot of corresponding CHMfactor;
 *		     only used for "diagBack"
 * @param resultKind a (SEXP) string indicating which kind of result is desired.
 *
 * @return  a SEXP, either a (double) number or a length n-vector of diagonal entries
 */
SEXP diag_tC(SEXP pslot, SEXP xslot, SEXP perm_slot, SEXP resultKind)
{
    int n = length(pslot) - 1, /* n = ncol(.) = nrow(.) */
	*x_p  = INTEGER(pslot),
	*perm = INTEGER(perm_slot);
    double *x_x = REAL(xslot);
/*  ^^^^^^        ^^^^ FIXME[Generalize] to INTEGER(.) / LOGICAL(.) / ... xslot !*/

    return diag_tC_ptr(n, x_p, x_x, perm, resultKind);
}

/**
 * Create a Csparse matrix object from indices and/or pointers.
 *
 * @param cls name of actual class of object to create
 * @param i optional integer vector of length nnz of row indices
 * @param j optional integer vector of length nnz of column indices
 * @param p optional integer vector of length np of row or column pointers
 * @param np length of integer vector p.  Must be zero if p == (int*)NULL
 * @param x optional vector of values
 * @param nnz length of vectors i, j and/or x, whichever is to be used
 * @param dims optional integer vector of length 2 to be used as
 *     dimensions.  If dims == (int*)NULL then the maximum row and column
 *     index are used as the dimensions.
 * @param dimnames optional list of length 2 to be used as dimnames
 * @param index1 indicator of 1-based indices
 *
 * @return an SEXP of class cls inheriting from CsparseMatrix.
 */
SEXP create_Csparse(char* cls, int* i, int* j, int* p, int np,
		    void* x, int nnz, int* dims, SEXP dimnames,
		    int index1)
{
    SEXP ans;
    int *ij = (int*)NULL, *tri, *trj,
	mi, mj, mp, nrow = -1, ncol = -1;
    int xtype = -1;		/* -Wall */
    CHM_TR T;
    CHM_SP A;

    if (np < 0 || nnz < 0)
	error(_("negative vector lengths not allowed: np = %d, nnz = %d"),
	      np, nnz);
    if (1 != ((mi = (i == (int*)NULL)) +
	      (mj = (j == (int*)NULL)) +
	      (mp = (p == (int*)NULL))))
	error(_("exactly 1 of 'i', 'j' or 'p' must be NULL"));
    if (mp) {
	if (np) error(_("np = %d, must be zero when p is NULL"), np);
    } else {
	if (np) {		/* Expand p to form i or j */
	    if (!(p[0])) error(_("p[0] = %d, should be zero"), p[0]);
	    for (int ii = 0; ii < np; ii++)
		if (p[ii] > p[ii + 1])
		    error(_("p must be non-decreasing"));
	    if (p[np] != nnz)
		error("p[np] = %d != nnz = %d", p[np], nnz);
	    ij = Calloc(nnz, int);
	    if (mi) {
		i = ij;
		nrow = np;
	    } else {
		j = ij;
		ncol = np;
	    }
	    /* Expand p to 0-based indices */
	    for (int ii = 0; ii < np; ii++)
		for (int jj = p[ii]; jj < p[ii + 1]; jj++) ij[jj] = ii;
	} else {
	    if (nnz)
		error(_("Inconsistent dimensions: np = 0 and nnz = %d"),
		      nnz);
	}
    }
    /* calculate nrow and ncol */
    if (nrow < 0) {
	for (int ii = 0; ii < nnz; ii++) {
	    int i1 = i[ii] + (index1 ? 0 : 1); /* 1-based index */
	    if (i1 < 1) error(_("invalid row index at position %d"), ii);
	    if (i1 > nrow) nrow = i1;
	}
    }
    if (ncol < 0) {
	for (int jj = 0; jj < nnz; jj++) {
	    int j1 = j[jj] + (index1 ? 0 : 1);
	    if (j1 < 1) error(_("invalid column index at position %d"), jj);
	    if (j1 > ncol) ncol = j1;
	}
    }
    if (dims != (int*)NULL) {
	if (dims[0] > nrow) nrow = dims[0];
	if (dims[1] > ncol) ncol = dims[1];
    }
    /* check the class name */
    if (strlen(cls) != 8)
	error(_("strlen of cls argument = %d, should be 8"), strlen(cls));
    if (!strcmp(cls + 2, "CMatrix"))
	error(_("cls = \"%s\" does not end in \"CMatrix\""), cls);
    switch(cls[0]) {
    case 'd':
    case 'l':
	xtype = CHOLMOD_REAL;
    break;
    case 'n':
	xtype = CHOLMOD_PATTERN;
	break;
    default:
	error(_("cls = \"%s\" must begin with 'd', 'l' or 'n'"), cls);
    }
    if (cls[1] != 'g')
	error(_("Only 'g'eneral sparse matrix types allowed"));
    /* allocate and populate the triplet */
    T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0,
				 xtype, &c);
    T->x = x;
    tri = (int*)T->i;
    trj = (int*)T->j;
    for (int ii = 0; ii < nnz; ii++) {
	tri[ii] = i[ii] - ((!mi && index1) ? 1 : 0);
	trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0);
    }
    /* create the cholmod_sparse structure */
    A = cholmod_triplet_to_sparse(T, nnz, &c);
    cholmod_free_triplet(&T, &c);
    /* copy the information to the SEXP */
    ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls)));
/* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */
    /* allocate and copy common slots */
    nnz = cholmod_nnz(A, &c);
    dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2));
    dims[0] = A->nrow; dims[1] = A->ncol;
    Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1);
    Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), (int*)A->i, nnz);
    switch(cls[1]) {
    case 'd':
	Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), (double*)A->x, nnz);
	break;
    case 'l':
	error(_("code not yet written for cls = \"lgCMatrix\""));
    }
/* FIXME: dimnames are *NOT* put there yet (if non-NULL) */
    cholmod_free_sparse(&A, &c);
    UNPROTECT(1);
    return ans;
}

root@r-forge.r-project.org
ViewVC Help
Powered by ViewVC 1.0.0  
Thanks to:
Vienna University of Economics and Business Powered By FusionForge