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

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