SCM

SCM Repository

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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 2731 - (view) (download) (as text)

1 : bates 1218 /* Sparse matrices in compressed column-oriented form */
2 : mmaechler 2677
3 : bates 922 #include "Csparse.h"
4 : maechler 2120 #include "Tsparse.h"
5 : bates 922 #include "chm_common.h"
6 :    
7 : mmaechler 2279 /** "Cheap" C version of Csparse_validate() - *not* sorting : */
8 :     Rboolean isValid_Csparse(SEXP x)
9 :     {
10 :     /* NB: we do *NOT* check a potential 'x' slot here, at all */
11 :     SEXP pslot = GET_SLOT(x, Matrix_pSym),
12 :     islot = GET_SLOT(x, Matrix_iSym);
13 :     int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), j,
14 :     nrow = dims[0],
15 :     ncol = dims[1],
16 :     *xp = INTEGER(pslot),
17 :     *xi = INTEGER(islot);
18 :    
19 :     if (length(pslot) != dims[1] + 1)
20 :     return FALSE;
21 :     if (xp[0] != 0)
22 :     return FALSE;
23 :     if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
24 :     return FALSE;
25 :     for (j = 0; j < xp[ncol]; j++) {
26 :     if (xi[j] < 0 || xi[j] >= nrow)
27 :     return FALSE;
28 :     }
29 :     for (j = 0; j < ncol; j++) {
30 :     if (xp[j] > xp[j + 1])
31 :     return FALSE;
32 :     }
33 :     return TRUE;
34 :     }
35 :    
36 : mmaechler 2312 SEXP Csparse_validate(SEXP x) {
37 :     return Csparse_validate_(x, FALSE);
38 :     }
39 :    
40 :     SEXP Csparse_validate2(SEXP x, SEXP maybe_modify) {
41 :     return Csparse_validate_(x, asLogical(maybe_modify));
42 :     }
43 :    
44 :     SEXP Csparse_validate_(SEXP x, Rboolean maybe_modify)
45 : bates 922 {
46 : maechler 1575 /* NB: we do *NOT* check a potential 'x' slot here, at all */
47 : bates 922 SEXP pslot = GET_SLOT(x, Matrix_pSym),
48 :     islot = GET_SLOT(x, Matrix_iSym);
49 : maechler 1893 Rboolean sorted, strictly;
50 :     int j, k,
51 : bates 922 *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
52 : maechler 1660 nrow = dims[0],
53 :     ncol = dims[1],
54 : maechler 1654 *xp = INTEGER(pslot),
55 : bates 922 *xi = INTEGER(islot);
56 :    
57 : maechler 1654 if (length(pslot) != dims[1] + 1)
58 :     return mkString(_("slot p must have length = ncol(.) + 1"));
59 : bates 922 if (xp[0] != 0)
60 :     return mkString(_("first element of slot p must be zero"));
61 : maechler 1893 if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
62 : bates 1555 return
63 :     mkString(_("last element of slot p must match length of slots i and x"));
64 : mmaechler 2236 for (j = 0; j < xp[ncol]; j++) {
65 : bates 1555 if (xi[j] < 0 || xi[j] >= nrow)
66 :     return mkString(_("all row indices must be between 0 and nrow-1"));
67 :     }
68 : maechler 1893 sorted = TRUE; strictly = TRUE;
69 : bates 922 for (j = 0; j < ncol; j++) {
70 : mmaechler 2236 if (xp[j] > xp[j + 1])
71 : bates 922 return mkString(_("slot p must be non-decreasing"));
72 : mmaechler 2236 if(sorted) /* only act if >= 2 entries in column j : */
73 : maechler 1893 for (k = xp[j] + 1; k < xp[j + 1]; k++) {
74 :     if (xi[k] < xi[k - 1])
75 :     sorted = FALSE;
76 :     else if (xi[k] == xi[k - 1])
77 :     strictly = FALSE;
78 :     }
79 : bates 922 }
80 : maechler 1654 if (!sorted) {
81 : mmaechler 2312 if(maybe_modify) {
82 :     CHM_SP chx = (CHM_SP) alloca(sizeof(cholmod_sparse));
83 :     R_CheckStack();
84 :     as_cholmod_sparse(chx, x, FALSE, TRUE);/*-> cholmod_l_sort() ! */
85 :     /* as chx = AS_CHM_SP__(x) but ^^^^ sorting x in_place !!! */
86 :    
87 :     /* Now re-check that row indices are *strictly* increasing
88 :     * (and not just increasing) within each column : */
89 :     for (j = 0; j < ncol; j++) {
90 :     for (k = xp[j] + 1; k < xp[j + 1]; k++)
91 :     if (xi[k] == xi[k - 1])
92 :     return mkString(_("slot i is not *strictly* increasing inside a column (even after cholmod_l_sort)"));
93 :     }
94 :     } else { /* no modifying sorting : */
95 :     return mkString(_("row indices are not sorted within columns"));
96 :     }
97 : maechler 1893 } else if(!strictly) { /* sorted, but not strictly */
98 :     return mkString(_("slot i is not *strictly* increasing inside a column"));
99 : maechler 1654 }
100 : bates 922 return ScalarLogical(1);
101 :     }
102 :    
103 : maechler 1968 SEXP Rsparse_validate(SEXP x)
104 :     {
105 :     /* NB: we do *NOT* check a potential 'x' slot here, at all */
106 :     SEXP pslot = GET_SLOT(x, Matrix_pSym),
107 :     jslot = GET_SLOT(x, Matrix_jSym);
108 :     Rboolean sorted, strictly;
109 :     int i, k,
110 :     *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
111 :     nrow = dims[0],
112 :     ncol = dims[1],
113 :     *xp = INTEGER(pslot),
114 :     *xj = INTEGER(jslot);
115 :    
116 :     if (length(pslot) != dims[0] + 1)
117 :     return mkString(_("slot p must have length = nrow(.) + 1"));
118 :     if (xp[0] != 0)
119 :     return mkString(_("first element of slot p must be zero"));
120 :     if (length(jslot) < xp[nrow]) /* allow larger slots from over-allocation!*/
121 :     return
122 :     mkString(_("last element of slot p must match length of slots j and x"));
123 :     for (i = 0; i < length(jslot); i++) {
124 :     if (xj[i] < 0 || xj[i] >= ncol)
125 :     return mkString(_("all column indices must be between 0 and ncol-1"));
126 :     }
127 :     sorted = TRUE; strictly = TRUE;
128 :     for (i = 0; i < nrow; i++) {
129 :     if (xp[i] > xp[i+1])
130 :     return mkString(_("slot p must be non-decreasing"));
131 :     if(sorted)
132 :     for (k = xp[i] + 1; k < xp[i + 1]; k++) {
133 :     if (xj[k] < xj[k - 1])
134 :     sorted = FALSE;
135 :     else if (xj[k] == xj[k - 1])
136 :     strictly = FALSE;
137 :     }
138 :     }
139 :     if (!sorted)
140 : mmaechler 2661 /* cannot easily use cholmod_sort(.) ... -> "error out" :*/
141 : maechler 1968 return mkString(_("slot j is not increasing inside a column"));
142 :     else if(!strictly) /* sorted, but not strictly */
143 :     return mkString(_("slot j is not *strictly* increasing inside a column"));
144 :    
145 :     return ScalarLogical(1);
146 :     }
147 :    
148 :    
149 : maechler 1751 /* Called from ../R/Csparse.R : */
150 :     /* Can only return [dln]geMatrix (no symm/triang);
151 :     * FIXME: replace by non-CHOLMOD code ! */
152 : bates 1059 SEXP Csparse_to_dense(SEXP x)
153 :     {
154 : mmaechler 2223 CHM_SP chxs = AS_CHM_SP__(x);
155 : maechler 1751 /* This loses the symmetry property, since cholmod_dense has none,
156 :     * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices
157 :     * to numeric (CHOLMOD_REAL) ones : */
158 : mmaechler 2661 CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c);
159 : maechler 1751 int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x);
160 : maechler 1960 R_CheckStack();
161 : bates 1059
162 : maechler 1736 return chm_dense_to_SEXP(chxd, 1, Rkind, GET_SLOT(x, Matrix_DimNamesSym));
163 : bates 1059 }
164 :    
165 : mmaechler 2628 // FIXME: do not go via CHM (should not be too hard, to just *drop* the x-slot, right?
166 : maechler 1548 SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri)
167 : bates 1371 {
168 : mmaechler 2223 CHM_SP chxs = AS_CHM_SP__(x);
169 : mmaechler 2661 CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c);
170 : bates 1867 int tr = asLogical(tri);
171 : maechler 1960 R_CheckStack();
172 : bates 1371
173 : maechler 1960 return chm_sparse_to_SEXP(chxcp, 1/*do_free*/,
174 : bates 1867 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
175 :     0, tr ? diag_P(x) : "",
176 : maechler 1548 GET_SLOT(x, Matrix_DimNamesSym));
177 : bates 1371 }
178 :    
179 : mmaechler 2628 // n.CMatrix --> [dli].CMatrix (not going through CHM!)
180 :     SEXP nz_pattern_to_Csparse(SEXP x, SEXP res_kind)
181 :     {
182 :     return nz2Csparse(x, asInteger(res_kind));
183 :     }
184 :     // n.CMatrix --> [dli].CMatrix (not going through CHM!)
185 :     SEXP nz2Csparse(SEXP x, enum x_slot_kind r_kind)
186 :     {
187 :     const char *cl_x = class_P(x);
188 :     if(cl_x[0] != 'n') error(_("not a 'n.CMatrix'"));
189 :     if(cl_x[2] != 'C') error(_("not a CsparseMatrix"));
190 :     int nnz = LENGTH(GET_SLOT(x, Matrix_iSym));
191 :     SEXP ans;
192 :     char *ncl = strdup(cl_x);
193 :     double *dx_x; int *ix_x;
194 :     ncl[0] = (r_kind == x_double ? 'd' :
195 :     (r_kind == x_logical ? 'l' :
196 :     /* else (for now): r_kind == x_integer : */ 'i'));
197 :     PROTECT(ans = NEW_OBJECT(MAKE_CLASS(ncl)));
198 :     // create a correct 'x' slot:
199 :     switch(r_kind) {
200 :     int i;
201 :     case x_double: // 'd'
202 :     dx_x = REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz));
203 :     for (i=0; i < nnz; i++) dx_x[i] = 1.;
204 :     break;
205 :     case x_logical: // 'l'
206 :     ix_x = LOGICAL(ALLOC_SLOT(ans, Matrix_xSym, LGLSXP, nnz));
207 :     for (i=0; i < nnz; i++) ix_x[i] = TRUE;
208 :     break;
209 :     case x_integer: // 'i'
210 :     ix_x = INTEGER(ALLOC_SLOT(ans, Matrix_xSym, INTSXP, nnz));
211 :     for (i=0; i < nnz; i++) ix_x[i] = 1;
212 :     break;
213 :    
214 :     default:
215 :     error(_("nz2Csparse(): invalid/non-implemented r_kind = %d"),
216 :     r_kind);
217 :     }
218 :    
219 :     // now copy all other slots :
220 :     slot_dup(ans, x, Matrix_iSym);
221 :     slot_dup(ans, x, Matrix_pSym);
222 :     slot_dup(ans, x, Matrix_DimSym);
223 :     slot_dup(ans, x, Matrix_DimNamesSym);
224 :     if(ncl[1] != 'g') { // symmetric or triangular ...
225 :     slot_dup_if_has(ans, x, Matrix_uploSym);
226 :     slot_dup_if_has(ans, x, Matrix_diagSym);
227 :     }
228 :     UNPROTECT(1);
229 :     return ans;
230 :     }
231 :    
232 : bates 1366 SEXP Csparse_to_matrix(SEXP x)
233 : bates 922 {
234 : mmaechler 2661 return chm_dense_to_matrix(cholmod_sparse_to_dense(AS_CHM_SP__(x), &c),
235 : maechler 1960 1 /*do_free*/, GET_SLOT(x, Matrix_DimNamesSym));
236 : bates 1366 }
237 :    
238 :     SEXP Csparse_to_Tsparse(SEXP x, SEXP tri)
239 :     {
240 : mmaechler 2223 CHM_SP chxs = AS_CHM_SP__(x);
241 : mmaechler 2661 CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c);
242 : bates 1867 int tr = asLogical(tri);
243 : maechler 1736 int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
244 : maechler 1960 R_CheckStack();
245 : bates 922
246 : bates 1867 return chm_triplet_to_SEXP(chxt, 1,
247 :     tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
248 :     Rkind, tr ? diag_P(x) : "",
249 : bates 1371 GET_SLOT(x, Matrix_DimNamesSym));
250 : bates 922 }
251 :    
252 : bates 1448 /* this used to be called sCMatrix_to_gCMatrix(..) [in ./dsCMatrix.c ]: */
253 : bates 1371 SEXP Csparse_symmetric_to_general(SEXP x)
254 :     {
255 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chgx;
256 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
257 : maechler 1960 R_CheckStack();
258 : bates 1371
259 :     if (!(chx->stype))
260 : maechler 1548 error(_("Nonsymmetric matrix in Csparse_symmetric_to_general"));
261 : mmaechler 2661 chgx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
262 : maechler 1375 /* xtype: pattern, "real", complex or .. */
263 : maechler 1548 return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
264 : bates 1371 GET_SLOT(x, Matrix_DimNamesSym));
265 :     }
266 :    
267 : maechler 1618 SEXP Csparse_general_to_symmetric(SEXP x, SEXP uplo)
268 : maechler 1598 {
269 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chgx;
270 : maechler 2113 int uploT = (*CHAR(STRING_ELT(uplo,0)) == 'U') ? 1 : -1;
271 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
272 : maechler 1960 R_CheckStack();
273 : maechler 1598
274 : mmaechler 2661 chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c);
275 : maechler 1598 /* xtype: pattern, "real", complex or .. */
276 :     return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
277 :     GET_SLOT(x, Matrix_DimNamesSym));
278 :     }
279 :    
280 : bates 1369 SEXP Csparse_transpose(SEXP x, SEXP tri)
281 : bates 922 {
282 : maechler 1921 /* TODO: lgCMatrix & igC* currently go via double prec. cholmod -
283 :     * since cholmod (& cs) lacks sparse 'int' matrices */
284 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
285 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
286 : mmaechler 2661 CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c);
287 : bates 1366 SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp;
288 : bates 1867 int tr = asLogical(tri);
289 : maechler 1960 R_CheckStack();
290 : bates 1369
291 : bates 1366 tmp = VECTOR_ELT(dn, 0); /* swap the dimnames */
292 :     SET_VECTOR_ELT(dn, 0, VECTOR_ELT(dn, 1));
293 :     SET_VECTOR_ELT(dn, 1, tmp);
294 :     UNPROTECT(1);
295 : bates 1867 return chm_sparse_to_SEXP(chxt, 1, /* SWAP 'uplo' for triangular */
296 :     tr ? ((*uplo_P(x) == 'U') ? -1 : 1) : 0,
297 :     Rkind, tr ? diag_P(x) : "", dn);
298 : bates 922 }
299 :    
300 :     SEXP Csparse_Csparse_prod(SEXP a, SEXP b)
301 :     {
302 : maechler 2120 CHM_SP
303 : mmaechler 2223 cha = AS_CHM_SP(a),
304 :     chb = AS_CHM_SP(b),
305 : mmaechler 2661 chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0,
306 : mmaechler 2490 /* values:= is_numeric (T/F) */ cha->xtype > 0,
307 :     /*out sorted:*/ 1, &c);
308 : maechler 2125 const char *cl_a = class_P(a), *cl_b = class_P(b);
309 :     char diag[] = {'\0', '\0'};
310 :     int uploT = 0;
311 : mmaechler 2495 SEXP dn = PROTECT(allocVector(VECSXP, 2));
312 : maechler 1960 R_CheckStack();
313 : bates 922
314 : mmaechler 2494 #ifdef DEBUG_Matrix_verbose
315 : mmaechler 2490 Rprintf("DBG Csparse_C*_prod(%s, %s)\n", cl_a, cl_b);
316 :     #endif
317 :    
318 : maechler 2125 /* Preserve triangularity and even unit-triangularity if appropriate.
319 :     * Note that in that case, the multiplication itself should happen
320 :     * faster. But there's no support for that in CHOLMOD */
321 :    
322 :     /* UGLY hack -- rather should have (fast!) C-level version of
323 :     * is(a, "triangularMatrix") etc */
324 :     if (cl_a[1] == 't' && cl_b[1] == 't')
325 :     /* FIXME: fails for "Cholesky","BunchKaufmann"..*/
326 :     if(*uplo_P(a) == *uplo_P(b)) { /* both upper, or both lower tri. */
327 :     uploT = (*uplo_P(a) == 'U') ? 1 : -1;
328 :     if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
329 :     /* "remove the diagonal entries": */
330 :     chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
331 :     diag[0]= 'U';
332 :     }
333 :     else diag[0]= 'N';
334 :     }
335 : bates 1366 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
336 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
337 :     SET_VECTOR_ELT(dn, 1,
338 :     duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), 1)));
339 : mmaechler 2495 UNPROTECT(1);
340 : maechler 2125 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
341 : bates 922 }
342 :    
343 : maechler 1659 SEXP Csparse_Csparse_crossprod(SEXP a, SEXP b, SEXP trans)
344 : bates 1657 {
345 : maechler 1659 int tr = asLogical(trans);
346 : maechler 2120 CHM_SP
347 : mmaechler 2223 cha = AS_CHM_SP(a),
348 :     chb = AS_CHM_SP(b),
349 : maechler 2120 chTr, chc;
350 : maechler 2125 const char *cl_a = class_P(a), *cl_b = class_P(b);
351 :     char diag[] = {'\0', '\0'};
352 :     int uploT = 0;
353 : mmaechler 2495 SEXP dn = PROTECT(allocVector(VECSXP, 2));
354 : maechler 1960 R_CheckStack();
355 : bates 1657
356 : mmaechler 2661 chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c);
357 :     chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb,
358 : maechler 2125 /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c);
359 : mmaechler 2661 cholmod_free_sparse(&chTr, &c);
360 : maechler 1659
361 : maechler 2125 /* Preserve triangularity and unit-triangularity if appropriate;
362 :     * see Csparse_Csparse_prod() for comments */
363 :     if (cl_a[1] == 't' && cl_b[1] == 't')
364 :     if(*uplo_P(a) != *uplo_P(b)) { /* one 'U', the other 'L' */
365 :     uploT = (*uplo_P(b) == 'U') ? 1 : -1;
366 :     if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
367 :     chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
368 :     diag[0]= 'U';
369 :     }
370 :     else diag[0]= 'N';
371 :     }
372 : bates 1657 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
373 : maechler 1659 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), (tr) ? 0 : 1)));
374 : bates 1657 SET_VECTOR_ELT(dn, 1,
375 : maechler 1659 duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), (tr) ? 0 : 1)));
376 : mmaechler 2495 UNPROTECT(1);
377 : maechler 2125 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
378 : bates 1657 }
379 :    
380 : bates 922 SEXP Csparse_dense_prod(SEXP a, SEXP b)
381 :     {
382 : mmaechler 2223 CHM_SP cha = AS_CHM_SP(a);
383 : maechler 1660 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
384 : maechler 1960 CHM_DN chb = AS_CHM_DN(b_M);
385 : mmaechler 2661 CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow,
386 : maechler 1960 chb->xtype, &c);
387 :     SEXP dn = PROTECT(allocVector(VECSXP, 2));
388 :     double one[] = {1,0}, zero[] = {0,0};
389 : mmaechler 2628 int nprot = 2;
390 : maechler 1960 R_CheckStack();
391 : mmaechler 2628 /* Tim Davis, please FIXME: currently (2010-11) *fails* when a is a pattern matrix:*/
392 :     if(cha->xtype == CHOLMOD_PATTERN) {
393 :     /* warning(_("Csparse_dense_prod(): cholmod_sdmult() not yet implemented for pattern./ ngCMatrix" */
394 :     /* " --> slightly inefficient coercion")); */
395 : bates 922
396 : mmaechler 2628 // This *fails* to produce a CHOLMOD_REAL ..
397 :     // CHM_SP chd = cholmod_l_copy(cha, cha->stype, CHOLMOD_REAL, &c);
398 :     // --> use our Matrix-classes
399 :     SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
400 :     cha = AS_CHM_SP(da);
401 :     }
402 : mmaechler 2661 cholmod_sdmult(cha, 0, one, zero, chb, chc, &c);
403 : maechler 1660 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
404 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
405 :     SET_VECTOR_ELT(dn, 1,
406 :     duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
407 : mmaechler 2628 UNPROTECT(nprot);
408 : maechler 1660 return chm_dense_to_SEXP(chc, 1, 0, dn);
409 : bates 922 }
410 : maechler 925
411 : bates 1067 SEXP Csparse_dense_crossprod(SEXP a, SEXP b)
412 :     {
413 : mmaechler 2223 CHM_SP cha = AS_CHM_SP(a);
414 : maechler 1660 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
415 : maechler 1960 CHM_DN chb = AS_CHM_DN(b_M);
416 : mmaechler 2661 CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol,
417 : maechler 1960 chb->xtype, &c);
418 : mmaechler 2629 SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2;
419 : maechler 1960 double one[] = {1,0}, zero[] = {0,0};
420 :     R_CheckStack();
421 : mmaechler 2629 // -- see Csparse_dense_prod() above :
422 :     if(cha->xtype == CHOLMOD_PATTERN) {
423 :     SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
424 :     cha = AS_CHM_SP(da);
425 :     }
426 : mmaechler 2661 cholmod_sdmult(cha, 1, one, zero, chb, chc, &c);
427 : maechler 1660 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
428 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1)));
429 :     SET_VECTOR_ELT(dn, 1,
430 :     duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
431 : mmaechler 2629 UNPROTECT(nprot);
432 : maechler 1660 return chm_dense_to_SEXP(chc, 1, 0, dn);
433 : bates 1067 }
434 :    
435 : maechler 2125 /* Computes x'x or x x' -- *also* for Tsparse (triplet = TRUE)
436 :     see Csparse_Csparse_crossprod above for x'y and x y' */
437 : bates 928 SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet)
438 : bates 922 {
439 : maechler 957 int trip = asLogical(triplet),
440 :     tr = asLogical(trans); /* gets reversed because _aat is tcrossprod */
441 : mmaechler 2491 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
442 : mmaechler 2223 CHM_TR cht = trip ? AS_CHM_TR(x) : (CHM_TR) NULL;
443 : mmaechler 2491 #else /* workaround needed:*/
444 : mmaechler 2516 SEXP xx = PROTECT(Tsparse_diagU2N(x));
445 :     CHM_TR cht = trip ? AS_CHM_TR__(xx) : (CHM_TR) NULL;
446 : mmaechler 2491 #endif
447 : maechler 1960 CHM_SP chcp, chxt,
448 : maechler 2120 chx = (trip ?
449 : mmaechler 2661 cholmod_triplet_to_sparse(cht, cht->nnz, &c) :
450 : mmaechler 2223 AS_CHM_SP(x));
451 : bates 1366 SEXP dn = PROTECT(allocVector(VECSXP, 2));
452 : maechler 1960 R_CheckStack();
453 : bates 922
454 : mmaechler 2661 if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c);
455 :     chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c);
456 : maechler 2120 if(!chcp) {
457 :     UNPROTECT(1);
458 : mmaechler 2661 error(_("Csparse_crossprod(): error return from cholmod_aat()"));
459 : maechler 2120 }
460 : mmaechler 2661 cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c);
461 : bates 1360 chcp->stype = 1;
462 : mmaechler 2661 if (trip) cholmod_free_sparse(&chx, &c);
463 :     if (!tr) cholmod_free_sparse(&chxt, &c);
464 : maechler 1960 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
465 : bates 1366 duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym),
466 : maechler 1660 (tr) ? 0 : 1)));
467 : bates 1366 SET_VECTOR_ELT(dn, 1, duplicate(VECTOR_ELT(dn, 0)));
468 : mmaechler 2516 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
469 : bates 1366 UNPROTECT(1);
470 : mmaechler 2516 #else
471 :     UNPROTECT(2);
472 :     #endif
473 : maechler 1548 return chm_sparse_to_SEXP(chcp, 1, 0, 0, "", dn);
474 : bates 922 }
475 : bates 923
476 : mmaechler 2646 /* Csparse_drop(x, tol): drop entries with absolute value < tol, i.e,
477 :     * at least all "explicit" zeros */
478 : maechler 1618 SEXP Csparse_drop(SEXP x, SEXP tol)
479 :     {
480 : mmaechler 2175 const char *cl = class_P(x);
481 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
482 :     int tr = (cl[1] == 't');
483 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
484 : mmaechler 2661 CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c);
485 : maechler 1618 double dtol = asReal(tol);
486 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
487 : maechler 1960 R_CheckStack();
488 : maechler 1618
489 : mmaechler 2661 if(!cholmod_drop(dtol, ans, &c))
490 :     error(_("cholmod_drop() failed"));
491 : mmaechler 2725 return chm_sparse_to_SEXP(ans, 1,
492 : mmaechler 2175 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
493 :     Rkind, tr ? diag_P(x) : "",
494 : maechler 1736 GET_SLOT(x, Matrix_DimNamesSym));
495 : maechler 1618 }
496 :    
497 : bates 1218 SEXP Csparse_horzcat(SEXP x, SEXP y)
498 :     {
499 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
500 : mmaechler 2540 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
501 :     Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
502 :     Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
503 : maechler 1960 R_CheckStack();
504 : maechler 1375
505 : mmaechler 2540 /* TODO: currently drops dimnames - and we fix at R level */
506 : mmaechler 2661 return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c),
507 : maechler 1960 1, 0, Rkind, "", R_NilValue);
508 : bates 1218 }
509 :    
510 :     SEXP Csparse_vertcat(SEXP x, SEXP y)
511 :     {
512 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
513 : mmaechler 2540 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
514 :     Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
515 :     Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
516 : maechler 1960 R_CheckStack();
517 : maechler 1375
518 : mmaechler 2540 /* TODO: currently drops dimnames - and we fix at R level */
519 : mmaechler 2661 return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c),
520 : maechler 1960 1, 0, Rkind, "", R_NilValue);
521 : bates 1218 }
522 : bates 1265
523 :     SEXP Csparse_band(SEXP x, SEXP k1, SEXP k2)
524 :     {
525 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
526 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
527 : mmaechler 2661 CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c);
528 : maechler 1960 R_CheckStack();
529 : bates 1265
530 : maechler 1736 return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "",
531 :     GET_SLOT(x, Matrix_DimNamesSym));
532 : bates 1265 }
533 : bates 1366
534 :     SEXP Csparse_diagU2N(SEXP x)
535 :     {
536 : maechler 2120 const char *cl = class_P(x);
537 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
538 :     if (cl[1] != 't' || *diag_P(x) != 'U') {
539 : maechler 2125 /* "trivially fast" when not triangular (<==> no 'diag' slot),
540 :     or not *unit* triangular */
541 : maechler 1708 return (x);
542 :     }
543 : maechler 2125 else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */
544 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
545 : mmaechler 2661 CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c);
546 : maechler 1708 double one[] = {1, 0};
547 : mmaechler 2661 CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c);
548 : maechler 1710 int uploT = (*uplo_P(x) == 'U') ? 1 : -1;
549 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
550 : bates 1366
551 : maechler 1960 R_CheckStack();
552 : mmaechler 2661 cholmod_free_sparse(&eye, &c);
553 : maechler 1708 return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N",
554 : maechler 1736 GET_SLOT(x, Matrix_DimNamesSym));
555 : maechler 1708 }
556 : bates 1366 }
557 :    
558 : maechler 2125 SEXP Csparse_diagN2U(SEXP x)
559 :     {
560 :     const char *cl = class_P(x);
561 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
562 :     if (cl[1] != 't' || *diag_P(x) != 'N') {
563 :     /* "trivially fast" when not triangular (<==> no 'diag' slot),
564 :     or already *unit* triangular */
565 :     return (x);
566 :     }
567 :     else { /* triangular with diag='N'): now drop the diagonal */
568 :     /* duplicate, since chx will be modified: */
569 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(duplicate(x));
570 : maechler 2125 int uploT = (*uplo_P(x) == 'U') ? 1 : -1,
571 :     Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
572 :     R_CheckStack();
573 :    
574 :     chm_diagN2U(chx, uploT, /* do_realloc */ FALSE);
575 :    
576 :     return chm_sparse_to_SEXP(chx, /*dofree*/ 0/* or 1 ?? */,
577 :     uploT, Rkind, "U",
578 :     GET_SLOT(x, Matrix_DimNamesSym));
579 :     }
580 :     }
581 :    
582 : mmaechler 2519 /**
583 :     * "Indexing" aka subsetting : Compute x[i,j], also for vectors i and j
584 :     * Working via CHOLMOD_submatrix, see ./CHOLMOD/MatrixOps/cholmod_submatrix.c
585 :     * @param x CsparseMatrix
586 :     * @param i row indices (0-origin), or NULL (R's)
587 :     * @param j columns indices (0-origin), or NULL
588 :     *
589 :     * @return x[i,j] still CsparseMatrix --- currently, this loses dimnames
590 :     */
591 : bates 1366 SEXP Csparse_submatrix(SEXP x, SEXP i, SEXP j)
592 :     {
593 : mmaechler 2519 CHM_SP chx = AS_CHM_SP(x); /* << does diagU2N() when needed */
594 : bates 1366 int rsize = (isNull(i)) ? -1 : LENGTH(i),
595 :     csize = (isNull(j)) ? -1 : LENGTH(j);
596 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
597 : maechler 1960 R_CheckStack();
598 : bates 1366
599 :     if (rsize >= 0 && !isInteger(i))
600 :     error(_("Index i must be NULL or integer"));
601 :     if (csize >= 0 && !isInteger(j))
602 :     error(_("Index j must be NULL or integer"));
603 : maechler 1736
604 : dmbates 2731 if (!chx->stype) {/* non-symmetric Matrix */
605 :     return chm_sparse_to_SEXP(cholmod_submatrix(chx,
606 :     (rsize < 0) ? NULL : INTEGER(i), rsize,
607 :     (csize < 0) ? NULL : INTEGER(j), csize,
608 :     TRUE, TRUE, &c),
609 :     1, 0, Rkind, "",
610 :     /* FIXME: drops dimnames */ R_NilValue);
611 :     }
612 :     /* for now, cholmod_submatrix() only accepts "generalMatrix" */
613 :     CHM_SP tmp = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
614 :     CHM_SP ans = cholmod_submatrix(tmp,
615 :     (rsize < 0) ? NULL : INTEGER(i), rsize,
616 :     (csize < 0) ? NULL : INTEGER(j), csize,
617 :     TRUE, TRUE, &c);
618 :     cholmod_free_sparse(&tmp, &c);
619 :     return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "", R_NilValue);
620 : bates 1366 }
621 : bates 2049
622 : mmaechler 2684 #define _d_Csp_
623 :     #include "t_Csparse_subassign.c"
624 : mmaechler 2661
625 : mmaechler 2684 #define _l_Csp_
626 :     #include "t_Csparse_subassign.c"
627 : mmaechler 2661
628 : mmaechler 2684 #define _i_Csp_
629 :     #include "t_Csparse_subassign.c"
630 : mmaechler 2661
631 : mmaechler 2684 #define _n_Csp_
632 :     #include "t_Csparse_subassign.c"
633 : mmaechler 2677
634 : mmaechler 2684 #define _z_Csp_
635 :     #include "t_Csparse_subassign.c"
636 : mmaechler 2677
637 : mmaechler 2661
638 : mmaechler 2677
639 : bates 2049 SEXP Csparse_MatrixMarket(SEXP x, SEXP fname)
640 :     {
641 :     FILE *f = fopen(CHAR(asChar(fname)), "w");
642 :    
643 :     if (!f)
644 :     error(_("failure to open file \"%s\" for writing"),
645 :     CHAR(asChar(fname)));
646 : mmaechler 2661 if (!cholmod_write_sparse(f, AS_CHM_SP(x),
647 : maechler 2120 (CHM_SP)NULL, (char*) NULL, &c))
648 : mmaechler 2661 error(_("cholmod_write_sparse returned error code"));
649 : bates 2049 fclose(f);
650 :     return R_NilValue;
651 :     }
652 : maechler 2137
653 :    
654 :     /**
655 :     * Extract the diagonal entries from *triangular* Csparse matrix __or__ a
656 :     * cholmod_sparse factor (LDL = TRUE).
657 :     *
658 :     * @param n dimension of the matrix.
659 :     * @param x_p 'p' (column pointer) slot contents
660 :     * @param x_x 'x' (non-zero entries) slot contents
661 : mmaechler 2175 * @param perm 'perm' (= permutation vector) slot contents; only used for "diagBack"
662 : maechler 2137 * @param resultKind a (SEXP) string indicating which kind of result is desired.
663 :     *
664 :     * @return a SEXP, either a (double) number or a length n-vector of diagonal entries
665 :     */
666 :     SEXP diag_tC_ptr(int n, int *x_p, double *x_x, int *perm, SEXP resultKind)
667 :     /* ^^^^^^ FIXME[Generalize] to int / ... */
668 :     {
669 :     const char* res_ch = CHAR(STRING_ELT(resultKind,0));
670 :     enum diag_kind { diag, diag_backpermuted, trace, prod, sum_log
671 :     } res_kind = ((!strcmp(res_ch, "trace")) ? trace :
672 :     ((!strcmp(res_ch, "sumLog")) ? sum_log :
673 :     ((!strcmp(res_ch, "prod")) ? prod :
674 :     ((!strcmp(res_ch, "diag")) ? diag :
675 :     ((!strcmp(res_ch, "diagBack")) ? diag_backpermuted :
676 :     -1)))));
677 :     int i, n_x, i_from = 0;
678 :     SEXP ans = PROTECT(allocVector(REALSXP,
679 :     /* ^^^^ FIXME[Generalize] */
680 :     (res_kind == diag ||
681 :     res_kind == diag_backpermuted) ? n : 1));
682 :     double *v = REAL(ans);
683 :     /* ^^^^^^ ^^^^ FIXME[Generalize] */
684 :    
685 :     #define for_DIAG(v_ASSIGN) \
686 :     for(i = 0; i < n; i++, i_from += n_x) { \
687 :     /* looking at i-th column */ \
688 :     n_x = x_p[i+1] - x_p[i];/* #{entries} in this column */ \
689 :     v_ASSIGN; \
690 :     }
691 :    
692 :     /* NOTA BENE: we assume -- uplo = "L" i.e. lower triangular matrix
693 :     * for uplo = "U" (makes sense with a "dtCMatrix" !),
694 :     * should use x_x[i_from + (nx - 1)] instead of x_x[i_from],
695 :     * where nx = (x_p[i+1] - x_p[i])
696 :     */
697 :    
698 :     switch(res_kind) {
699 :     case trace:
700 :     v[0] = 0.;
701 :     for_DIAG(v[0] += x_x[i_from]);
702 :     break;
703 :    
704 :     case sum_log:
705 :     v[0] = 0.;
706 :     for_DIAG(v[0] += log(x_x[i_from]));
707 :     break;
708 :    
709 :     case prod:
710 :     v[0] = 1.;
711 :     for_DIAG(v[0] *= x_x[i_from]);
712 :     break;
713 :    
714 :     case diag:
715 :     for_DIAG(v[i] = x_x[i_from]);
716 :     break;
717 :    
718 :     case diag_backpermuted:
719 :     for_DIAG(v[i] = x_x[i_from]);
720 :    
721 : mmaechler 2175 warning(_("resultKind = 'diagBack' (back-permuted) is experimental"));
722 : maechler 2144 /* now back_permute : */
723 :     for(i = 0; i < n; i++) {
724 :     double tmp = v[i]; v[i] = v[perm[i]]; v[perm[i]] = tmp;
725 :     /*^^^^ FIXME[Generalize] */
726 :     }
727 : maechler 2137 break;
728 :    
729 :     default: /* -1 from above */
730 : mmaechler 2387 error(_("diag_tC(): invalid 'resultKind'"));
731 : maechler 2137 /* Wall: */ ans = R_NilValue; v = REAL(ans);
732 :     }
733 :    
734 :     UNPROTECT(1);
735 :     return ans;
736 :     }
737 :    
738 :     /**
739 :     * Extract the diagonal entries from *triangular* Csparse matrix __or__ a
740 :     * cholmod_sparse factor (LDL = TRUE).
741 :     *
742 :     * @param pslot 'p' (column pointer) slot of Csparse matrix/factor
743 :     * @param xslot 'x' (non-zero entries) slot of Csparse matrix/factor
744 : mmaechler 2175 * @param perm_slot 'perm' (= permutation vector) slot of corresponding CHMfactor;
745 :     * only used for "diagBack"
746 : maechler 2137 * @param resultKind a (SEXP) string indicating which kind of result is desired.
747 :     *
748 :     * @return a SEXP, either a (double) number or a length n-vector of diagonal entries
749 :     */
750 :     SEXP diag_tC(SEXP pslot, SEXP xslot, SEXP perm_slot, SEXP resultKind)
751 :     {
752 :     int n = length(pslot) - 1, /* n = ncol(.) = nrow(.) */
753 :     *x_p = INTEGER(pslot),
754 :     *perm = INTEGER(perm_slot);
755 :     double *x_x = REAL(xslot);
756 :     /* ^^^^^^ ^^^^ FIXME[Generalize] to INTEGER(.) / LOGICAL(.) / ... xslot !*/
757 :    
758 :     return diag_tC_ptr(n, x_p, x_x, perm, resultKind);
759 :     }
760 : dmbates 2319
761 :     /**
762 :     * Create a Csparse matrix object from indices and/or pointers.
763 :     *
764 :     * @param cls name of actual class of object to create
765 :     * @param i optional integer vector of length nnz of row indices
766 :     * @param j optional integer vector of length nnz of column indices
767 :     * @param p optional integer vector of length np of row or column pointers
768 :     * @param np length of integer vector p. Must be zero if p == (int*)NULL
769 :     * @param x optional vector of values
770 :     * @param nnz length of vectors i, j and/or x, whichever is to be used
771 :     * @param dims optional integer vector of length 2 to be used as
772 :     * dimensions. If dims == (int*)NULL then the maximum row and column
773 :     * index are used as the dimensions.
774 :     * @param dimnames optional list of length 2 to be used as dimnames
775 :     * @param index1 indicator of 1-based indices
776 :     *
777 :     * @return an SEXP of class cls inheriting from CsparseMatrix.
778 :     */
779 :     SEXP create_Csparse(char* cls, int* i, int* j, int* p, int np,
780 :     void* x, int nnz, int* dims, SEXP dimnames,
781 :     int index1)
782 :     {
783 :     SEXP ans;
784 :     int *ij = (int*)NULL, *tri, *trj,
785 :     mi, mj, mp, nrow = -1, ncol = -1;
786 :     int xtype = -1; /* -Wall */
787 :     CHM_TR T;
788 :     CHM_SP A;
789 :    
790 :     if (np < 0 || nnz < 0)
791 :     error(_("negative vector lengths not allowed: np = %d, nnz = %d"),
792 :     np, nnz);
793 :     if (1 != ((mi = (i == (int*)NULL)) +
794 :     (mj = (j == (int*)NULL)) +
795 :     (mp = (p == (int*)NULL))))
796 :     error(_("exactly 1 of 'i', 'j' or 'p' must be NULL"));
797 :     if (mp) {
798 :     if (np) error(_("np = %d, must be zero when p is NULL"), np);
799 :     } else {
800 :     if (np) { /* Expand p to form i or j */
801 :     if (!(p[0])) error(_("p[0] = %d, should be zero"), p[0]);
802 :     for (int ii = 0; ii < np; ii++)
803 :     if (p[ii] > p[ii + 1])
804 :     error(_("p must be non-decreasing"));
805 :     if (p[np] != nnz)
806 : mmaechler 2387 error("p[np] = %d != nnz = %d", p[np], nnz);
807 : dmbates 2319 ij = Calloc(nnz, int);
808 :     if (mi) {
809 :     i = ij;
810 :     nrow = np;
811 :     } else {
812 :     j = ij;
813 :     ncol = np;
814 :     }
815 : mmaechler 2661 /* Expand p to 0-based indices */
816 : dmbates 2319 for (int ii = 0; ii < np; ii++)
817 :     for (int jj = p[ii]; jj < p[ii + 1]; jj++) ij[jj] = ii;
818 :     } else {
819 :     if (nnz)
820 :     error(_("Inconsistent dimensions: np = 0 and nnz = %d"),
821 :     nnz);
822 :     }
823 :     }
824 : mmaechler 2661 /* calculate nrow and ncol */
825 : dmbates 2319 if (nrow < 0) {
826 :     for (int ii = 0; ii < nnz; ii++) {
827 :     int i1 = i[ii] + (index1 ? 0 : 1); /* 1-based index */
828 :     if (i1 < 1) error(_("invalid row index at position %d"), ii);
829 :     if (i1 > nrow) nrow = i1;
830 :     }
831 :     }
832 :     if (ncol < 0) {
833 :     for (int jj = 0; jj < nnz; jj++) {
834 :     int j1 = j[jj] + (index1 ? 0 : 1);
835 : mmaechler 2387 if (j1 < 1) error(_("invalid column index at position %d"), jj);
836 : dmbates 2319 if (j1 > ncol) ncol = j1;
837 :     }
838 :     }
839 :     if (dims != (int*)NULL) {
840 :     if (dims[0] > nrow) nrow = dims[0];
841 :     if (dims[1] > ncol) ncol = dims[1];
842 :     }
843 : mmaechler 2661 /* check the class name */
844 : dmbates 2319 if (strlen(cls) != 8)
845 :     error(_("strlen of cls argument = %d, should be 8"), strlen(cls));
846 :     if (!strcmp(cls + 2, "CMatrix"))
847 :     error(_("cls = \"%s\" does not end in \"CMatrix\""), cls);
848 :     switch(cls[0]) {
849 :     case 'd':
850 :     case 'l':
851 : mmaechler 2661 xtype = CHOLMOD_REAL;
852 :     break;
853 : dmbates 2319 case 'n':
854 : mmaechler 2661 xtype = CHOLMOD_PATTERN;
855 :     break;
856 : dmbates 2319 default:
857 : mmaechler 2661 error(_("cls = \"%s\" must begin with 'd', 'l' or 'n'"), cls);
858 : dmbates 2319 }
859 :     if (cls[1] != 'g')
860 :     error(_("Only 'g'eneral sparse matrix types allowed"));
861 : mmaechler 2661 /* allocate and populate the triplet */
862 :     T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0,
863 :     xtype, &c);
864 : dmbates 2319 T->x = x;
865 :     tri = (int*)T->i;
866 :     trj = (int*)T->j;
867 :     for (int ii = 0; ii < nnz; ii++) {
868 :     tri[ii] = i[ii] - ((!mi && index1) ? 1 : 0);
869 :     trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0);
870 :     }
871 : mmaechler 2661 /* create the cholmod_sparse structure */
872 :     A = cholmod_triplet_to_sparse(T, nnz, &c);
873 :     cholmod_free_triplet(&T, &c);
874 :     /* copy the information to the SEXP */
875 : dmbates 2319 ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls)));
876 :     /* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */
877 : mmaechler 2661 /* allocate and copy common slots */
878 :     nnz = cholmod_nnz(A, &c);
879 : dmbates 2319 dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2));
880 :     dims[0] = A->nrow; dims[1] = A->ncol;
881 :     Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1);
882 :     Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), (int*)A->i, nnz);
883 :     switch(cls[1]) {
884 :     case 'd':
885 :     Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), (double*)A->x, nnz);
886 :     break;
887 :     case 'l':
888 :     error(_("code not yet written for cls = \"lgCMatrix\""));
889 :     }
890 : mmaechler 2646 /* FIXME: dimnames are *NOT* put there yet (if non-NULL) */
891 : mmaechler 2661 cholmod_free_sparse(&A, &c);
892 : dmbates 2319 UNPROTECT(1);
893 :     return ans;
894 :     }

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