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 2681 - (view) (download) (as text)

1 : bates 1218 /* Sparse matrices in compressed column-oriented form */
2 : mmaechler 2677
3 :     #include <stdint.h> // C99 for int64_t
4 : bates 922 #include "Csparse.h"
5 : maechler 2120 #include "Tsparse.h"
6 : bates 922 #include "chm_common.h"
7 :    
8 : mmaechler 2279 /** "Cheap" C version of Csparse_validate() - *not* sorting : */
9 :     Rboolean isValid_Csparse(SEXP x)
10 :     {
11 :     /* NB: we do *NOT* check a potential 'x' slot here, at all */
12 :     SEXP pslot = GET_SLOT(x, Matrix_pSym),
13 :     islot = GET_SLOT(x, Matrix_iSym);
14 :     int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), j,
15 :     nrow = dims[0],
16 :     ncol = dims[1],
17 :     *xp = INTEGER(pslot),
18 :     *xi = INTEGER(islot);
19 :    
20 :     if (length(pslot) != dims[1] + 1)
21 :     return FALSE;
22 :     if (xp[0] != 0)
23 :     return FALSE;
24 :     if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
25 :     return FALSE;
26 :     for (j = 0; j < xp[ncol]; j++) {
27 :     if (xi[j] < 0 || xi[j] >= nrow)
28 :     return FALSE;
29 :     }
30 :     for (j = 0; j < ncol; j++) {
31 :     if (xp[j] > xp[j + 1])
32 :     return FALSE;
33 :     }
34 :     return TRUE;
35 :     }
36 :    
37 : mmaechler 2312 SEXP Csparse_validate(SEXP x) {
38 :     return Csparse_validate_(x, FALSE);
39 :     }
40 :    
41 :     SEXP Csparse_validate2(SEXP x, SEXP maybe_modify) {
42 :     return Csparse_validate_(x, asLogical(maybe_modify));
43 :     }
44 :    
45 :     SEXP Csparse_validate_(SEXP x, Rboolean maybe_modify)
46 : bates 922 {
47 : maechler 1575 /* NB: we do *NOT* check a potential 'x' slot here, at all */
48 : bates 922 SEXP pslot = GET_SLOT(x, Matrix_pSym),
49 :     islot = GET_SLOT(x, Matrix_iSym);
50 : maechler 1893 Rboolean sorted, strictly;
51 :     int j, k,
52 : bates 922 *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
53 : maechler 1660 nrow = dims[0],
54 :     ncol = dims[1],
55 : maechler 1654 *xp = INTEGER(pslot),
56 : bates 922 *xi = INTEGER(islot);
57 :    
58 : maechler 1654 if (length(pslot) != dims[1] + 1)
59 :     return mkString(_("slot p must have length = ncol(.) + 1"));
60 : bates 922 if (xp[0] != 0)
61 :     return mkString(_("first element of slot p must be zero"));
62 : maechler 1893 if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/
63 : bates 1555 return
64 :     mkString(_("last element of slot p must match length of slots i and x"));
65 : mmaechler 2236 for (j = 0; j < xp[ncol]; j++) {
66 : bates 1555 if (xi[j] < 0 || xi[j] >= nrow)
67 :     return mkString(_("all row indices must be between 0 and nrow-1"));
68 :     }
69 : maechler 1893 sorted = TRUE; strictly = TRUE;
70 : bates 922 for (j = 0; j < ncol; j++) {
71 : mmaechler 2236 if (xp[j] > xp[j + 1])
72 : bates 922 return mkString(_("slot p must be non-decreasing"));
73 : mmaechler 2236 if(sorted) /* only act if >= 2 entries in column j : */
74 : maechler 1893 for (k = xp[j] + 1; k < xp[j + 1]; k++) {
75 :     if (xi[k] < xi[k - 1])
76 :     sorted = FALSE;
77 :     else if (xi[k] == xi[k - 1])
78 :     strictly = FALSE;
79 :     }
80 : bates 922 }
81 : maechler 1654 if (!sorted) {
82 : mmaechler 2312 if(maybe_modify) {
83 :     CHM_SP chx = (CHM_SP) alloca(sizeof(cholmod_sparse));
84 :     R_CheckStack();
85 :     as_cholmod_sparse(chx, x, FALSE, TRUE);/*-> cholmod_l_sort() ! */
86 :     /* as chx = AS_CHM_SP__(x) but ^^^^ sorting x in_place !!! */
87 :    
88 :     /* Now re-check that row indices are *strictly* increasing
89 :     * (and not just increasing) within each column : */
90 :     for (j = 0; j < ncol; j++) {
91 :     for (k = xp[j] + 1; k < xp[j + 1]; k++)
92 :     if (xi[k] == xi[k - 1])
93 :     return mkString(_("slot i is not *strictly* increasing inside a column (even after cholmod_l_sort)"));
94 :     }
95 :     } else { /* no modifying sorting : */
96 :     return mkString(_("row indices are not sorted within columns"));
97 :     }
98 : maechler 1893 } else if(!strictly) { /* sorted, but not strictly */
99 :     return mkString(_("slot i is not *strictly* increasing inside a column"));
100 : maechler 1654 }
101 : bates 922 return ScalarLogical(1);
102 :     }
103 :    
104 : maechler 1968 SEXP Rsparse_validate(SEXP x)
105 :     {
106 :     /* NB: we do *NOT* check a potential 'x' slot here, at all */
107 :     SEXP pslot = GET_SLOT(x, Matrix_pSym),
108 :     jslot = GET_SLOT(x, Matrix_jSym);
109 :     Rboolean sorted, strictly;
110 :     int i, k,
111 :     *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)),
112 :     nrow = dims[0],
113 :     ncol = dims[1],
114 :     *xp = INTEGER(pslot),
115 :     *xj = INTEGER(jslot);
116 :    
117 :     if (length(pslot) != dims[0] + 1)
118 :     return mkString(_("slot p must have length = nrow(.) + 1"));
119 :     if (xp[0] != 0)
120 :     return mkString(_("first element of slot p must be zero"));
121 :     if (length(jslot) < xp[nrow]) /* allow larger slots from over-allocation!*/
122 :     return
123 :     mkString(_("last element of slot p must match length of slots j and x"));
124 :     for (i = 0; i < length(jslot); i++) {
125 :     if (xj[i] < 0 || xj[i] >= ncol)
126 :     return mkString(_("all column indices must be between 0 and ncol-1"));
127 :     }
128 :     sorted = TRUE; strictly = TRUE;
129 :     for (i = 0; i < nrow; i++) {
130 :     if (xp[i] > xp[i+1])
131 :     return mkString(_("slot p must be non-decreasing"));
132 :     if(sorted)
133 :     for (k = xp[i] + 1; k < xp[i + 1]; k++) {
134 :     if (xj[k] < xj[k - 1])
135 :     sorted = FALSE;
136 :     else if (xj[k] == xj[k - 1])
137 :     strictly = FALSE;
138 :     }
139 :     }
140 :     if (!sorted)
141 : mmaechler 2661 /* cannot easily use cholmod_sort(.) ... -> "error out" :*/
142 : maechler 1968 return mkString(_("slot j is not increasing inside a column"));
143 :     else if(!strictly) /* sorted, but not strictly */
144 :     return mkString(_("slot j is not *strictly* increasing inside a column"));
145 :    
146 :     return ScalarLogical(1);
147 :     }
148 :    
149 :    
150 : maechler 1751 /* Called from ../R/Csparse.R : */
151 :     /* Can only return [dln]geMatrix (no symm/triang);
152 :     * FIXME: replace by non-CHOLMOD code ! */
153 : bates 1059 SEXP Csparse_to_dense(SEXP x)
154 :     {
155 : mmaechler 2223 CHM_SP chxs = AS_CHM_SP__(x);
156 : maechler 1751 /* This loses the symmetry property, since cholmod_dense has none,
157 :     * BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices
158 :     * to numeric (CHOLMOD_REAL) ones : */
159 : mmaechler 2661 CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c);
160 : maechler 1751 int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x);
161 : maechler 1960 R_CheckStack();
162 : bates 1059
163 : maechler 1736 return chm_dense_to_SEXP(chxd, 1, Rkind, GET_SLOT(x, Matrix_DimNamesSym));
164 : bates 1059 }
165 :    
166 : mmaechler 2628 // FIXME: do not go via CHM (should not be too hard, to just *drop* the x-slot, right?
167 : maechler 1548 SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri)
168 : bates 1371 {
169 : mmaechler 2223 CHM_SP chxs = AS_CHM_SP__(x);
170 : mmaechler 2661 CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c);
171 : bates 1867 int tr = asLogical(tri);
172 : maechler 1960 R_CheckStack();
173 : bates 1371
174 : maechler 1960 return chm_sparse_to_SEXP(chxcp, 1/*do_free*/,
175 : bates 1867 tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
176 :     0, tr ? diag_P(x) : "",
177 : maechler 1548 GET_SLOT(x, Matrix_DimNamesSym));
178 : bates 1371 }
179 :    
180 : mmaechler 2628 // n.CMatrix --> [dli].CMatrix (not going through CHM!)
181 :     SEXP nz_pattern_to_Csparse(SEXP x, SEXP res_kind)
182 :     {
183 :     return nz2Csparse(x, asInteger(res_kind));
184 :     }
185 :     // n.CMatrix --> [dli].CMatrix (not going through CHM!)
186 :     SEXP nz2Csparse(SEXP x, enum x_slot_kind r_kind)
187 :     {
188 :     const char *cl_x = class_P(x);
189 :     if(cl_x[0] != 'n') error(_("not a 'n.CMatrix'"));
190 :     if(cl_x[2] != 'C') error(_("not a CsparseMatrix"));
191 :     int nnz = LENGTH(GET_SLOT(x, Matrix_iSym));
192 :     SEXP ans;
193 :     char *ncl = strdup(cl_x);
194 :     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 2223 CHM_SP chx = AS_CHM_SP__(x), chgx;
271 : maechler 2113 int uploT = (*CHAR(STRING_ELT(uplo,0)) == 'U') ? 1 : -1;
272 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
273 : maechler 1960 R_CheckStack();
274 : maechler 1598
275 : mmaechler 2661 chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c);
276 : maechler 1598 /* xtype: pattern, "real", complex or .. */
277 :     return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "",
278 :     GET_SLOT(x, Matrix_DimNamesSym));
279 :     }
280 :    
281 : bates 1369 SEXP Csparse_transpose(SEXP x, SEXP tri)
282 : bates 922 {
283 : maechler 1921 /* TODO: lgCMatrix & igC* currently go via double prec. cholmod -
284 :     * since cholmod (& cs) lacks sparse 'int' matrices */
285 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
286 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
287 : mmaechler 2661 CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c);
288 : bates 1366 SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp;
289 : bates 1867 int tr = asLogical(tri);
290 : maechler 1960 R_CheckStack();
291 : bates 1369
292 : bates 1366 tmp = VECTOR_ELT(dn, 0); /* swap the dimnames */
293 :     SET_VECTOR_ELT(dn, 0, VECTOR_ELT(dn, 1));
294 :     SET_VECTOR_ELT(dn, 1, tmp);
295 :     UNPROTECT(1);
296 : bates 1867 return chm_sparse_to_SEXP(chxt, 1, /* SWAP 'uplo' for triangular */
297 :     tr ? ((*uplo_P(x) == 'U') ? -1 : 1) : 0,
298 :     Rkind, tr ? diag_P(x) : "", dn);
299 : bates 922 }
300 :    
301 :     SEXP Csparse_Csparse_prod(SEXP a, SEXP b)
302 :     {
303 : maechler 2120 CHM_SP
304 : mmaechler 2223 cha = AS_CHM_SP(a),
305 :     chb = AS_CHM_SP(b),
306 : mmaechler 2661 chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0,
307 : mmaechler 2490 /* values:= is_numeric (T/F) */ cha->xtype > 0,
308 :     /*out sorted:*/ 1, &c);
309 : maechler 2125 const char *cl_a = class_P(a), *cl_b = class_P(b);
310 :     char diag[] = {'\0', '\0'};
311 :     int uploT = 0;
312 : mmaechler 2495 SEXP dn = PROTECT(allocVector(VECSXP, 2));
313 : maechler 1960 R_CheckStack();
314 : bates 922
315 : mmaechler 2494 #ifdef DEBUG_Matrix_verbose
316 : mmaechler 2490 Rprintf("DBG Csparse_C*_prod(%s, %s)\n", cl_a, cl_b);
317 :     #endif
318 :    
319 : maechler 2125 /* Preserve triangularity and even unit-triangularity if appropriate.
320 :     * Note that in that case, the multiplication itself should happen
321 :     * faster. But there's no support for that in CHOLMOD */
322 :    
323 :     /* UGLY hack -- rather should have (fast!) C-level version of
324 :     * is(a, "triangularMatrix") etc */
325 :     if (cl_a[1] == 't' && cl_b[1] == 't')
326 :     /* FIXME: fails for "Cholesky","BunchKaufmann"..*/
327 :     if(*uplo_P(a) == *uplo_P(b)) { /* both upper, or both lower tri. */
328 :     uploT = (*uplo_P(a) == 'U') ? 1 : -1;
329 :     if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
330 :     /* "remove the diagonal entries": */
331 :     chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
332 :     diag[0]= 'U';
333 :     }
334 :     else diag[0]= 'N';
335 :     }
336 : bates 1366 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
337 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
338 :     SET_VECTOR_ELT(dn, 1,
339 :     duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), 1)));
340 : mmaechler 2495 UNPROTECT(1);
341 : maechler 2125 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
342 : bates 922 }
343 :    
344 : maechler 1659 SEXP Csparse_Csparse_crossprod(SEXP a, SEXP b, SEXP trans)
345 : bates 1657 {
346 : maechler 1659 int tr = asLogical(trans);
347 : maechler 2120 CHM_SP
348 : mmaechler 2223 cha = AS_CHM_SP(a),
349 :     chb = AS_CHM_SP(b),
350 : maechler 2120 chTr, chc;
351 : maechler 2125 const char *cl_a = class_P(a), *cl_b = class_P(b);
352 :     char diag[] = {'\0', '\0'};
353 :     int uploT = 0;
354 : mmaechler 2495 SEXP dn = PROTECT(allocVector(VECSXP, 2));
355 : maechler 1960 R_CheckStack();
356 : bates 1657
357 : mmaechler 2661 chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c);
358 :     chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb,
359 : maechler 2125 /*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c);
360 : mmaechler 2661 cholmod_free_sparse(&chTr, &c);
361 : maechler 1659
362 : maechler 2125 /* Preserve triangularity and unit-triangularity if appropriate;
363 :     * see Csparse_Csparse_prod() for comments */
364 :     if (cl_a[1] == 't' && cl_b[1] == 't')
365 :     if(*uplo_P(a) != *uplo_P(b)) { /* one 'U', the other 'L' */
366 :     uploT = (*uplo_P(b) == 'U') ? 1 : -1;
367 :     if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */
368 :     chm_diagN2U(chc, uploT, /* do_realloc */ FALSE);
369 :     diag[0]= 'U';
370 :     }
371 :     else diag[0]= 'N';
372 :     }
373 : bates 1657 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
374 : maechler 1659 duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), (tr) ? 0 : 1)));
375 : bates 1657 SET_VECTOR_ELT(dn, 1,
376 : maechler 1659 duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), (tr) ? 0 : 1)));
377 : mmaechler 2495 UNPROTECT(1);
378 : maechler 2125 return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn);
379 : bates 1657 }
380 :    
381 : bates 922 SEXP Csparse_dense_prod(SEXP a, SEXP b)
382 :     {
383 : mmaechler 2223 CHM_SP cha = AS_CHM_SP(a);
384 : maechler 1660 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
385 : maechler 1960 CHM_DN chb = AS_CHM_DN(b_M);
386 : mmaechler 2661 CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow,
387 : maechler 1960 chb->xtype, &c);
388 :     SEXP dn = PROTECT(allocVector(VECSXP, 2));
389 :     double one[] = {1,0}, zero[] = {0,0};
390 : mmaechler 2628 int nprot = 2;
391 : maechler 1960 R_CheckStack();
392 : mmaechler 2628 /* Tim Davis, please FIXME: currently (2010-11) *fails* when a is a pattern matrix:*/
393 :     if(cha->xtype == CHOLMOD_PATTERN) {
394 :     /* warning(_("Csparse_dense_prod(): cholmod_sdmult() not yet implemented for pattern./ ngCMatrix" */
395 :     /* " --> slightly inefficient coercion")); */
396 : bates 922
397 : mmaechler 2628 // This *fails* to produce a CHOLMOD_REAL ..
398 :     // CHM_SP chd = cholmod_l_copy(cha, cha->stype, CHOLMOD_REAL, &c);
399 :     // --> use our Matrix-classes
400 :     SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
401 :     cha = AS_CHM_SP(da);
402 :     }
403 : mmaechler 2661 cholmod_sdmult(cha, 0, one, zero, chb, chc, &c);
404 : maechler 1660 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
405 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0)));
406 :     SET_VECTOR_ELT(dn, 1,
407 :     duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
408 : mmaechler 2628 UNPROTECT(nprot);
409 : maechler 1660 return chm_dense_to_SEXP(chc, 1, 0, dn);
410 : bates 922 }
411 : maechler 925
412 : bates 1067 SEXP Csparse_dense_crossprod(SEXP a, SEXP b)
413 :     {
414 : mmaechler 2223 CHM_SP cha = AS_CHM_SP(a);
415 : maechler 1660 SEXP b_M = PROTECT(mMatrix_as_dgeMatrix(b));
416 : maechler 1960 CHM_DN chb = AS_CHM_DN(b_M);
417 : mmaechler 2661 CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol,
418 : maechler 1960 chb->xtype, &c);
419 : mmaechler 2629 SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2;
420 : maechler 1960 double one[] = {1,0}, zero[] = {0,0};
421 :     R_CheckStack();
422 : mmaechler 2629 // -- see Csparse_dense_prod() above :
423 :     if(cha->xtype == CHOLMOD_PATTERN) {
424 :     SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++;
425 :     cha = AS_CHM_SP(da);
426 :     }
427 : mmaechler 2661 cholmod_sdmult(cha, 1, one, zero, chb, chc, &c);
428 : maechler 1660 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
429 :     duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1)));
430 :     SET_VECTOR_ELT(dn, 1,
431 :     duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1)));
432 : mmaechler 2629 UNPROTECT(nprot);
433 : maechler 1660 return chm_dense_to_SEXP(chc, 1, 0, dn);
434 : bates 1067 }
435 :    
436 : maechler 2125 /* Computes x'x or x x' -- *also* for Tsparse (triplet = TRUE)
437 :     see Csparse_Csparse_crossprod above for x'y and x y' */
438 : bates 928 SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet)
439 : bates 922 {
440 : maechler 957 int trip = asLogical(triplet),
441 :     tr = asLogical(trans); /* gets reversed because _aat is tcrossprod */
442 : mmaechler 2491 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
443 : mmaechler 2223 CHM_TR cht = trip ? AS_CHM_TR(x) : (CHM_TR) NULL;
444 : mmaechler 2491 #else /* workaround needed:*/
445 : mmaechler 2516 SEXP xx = PROTECT(Tsparse_diagU2N(x));
446 :     CHM_TR cht = trip ? AS_CHM_TR__(xx) : (CHM_TR) NULL;
447 : mmaechler 2491 #endif
448 : maechler 1960 CHM_SP chcp, chxt,
449 : maechler 2120 chx = (trip ?
450 : mmaechler 2661 cholmod_triplet_to_sparse(cht, cht->nnz, &c) :
451 : mmaechler 2223 AS_CHM_SP(x));
452 : bates 1366 SEXP dn = PROTECT(allocVector(VECSXP, 2));
453 : maechler 1960 R_CheckStack();
454 : bates 922
455 : mmaechler 2661 if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c);
456 :     chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c);
457 : maechler 2120 if(!chcp) {
458 :     UNPROTECT(1);
459 : mmaechler 2661 error(_("Csparse_crossprod(): error return from cholmod_aat()"));
460 : maechler 2120 }
461 : mmaechler 2661 cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c);
462 : bates 1360 chcp->stype = 1;
463 : mmaechler 2661 if (trip) cholmod_free_sparse(&chx, &c);
464 :     if (!tr) cholmod_free_sparse(&chxt, &c);
465 : maechler 1960 SET_VECTOR_ELT(dn, 0, /* establish dimnames */
466 : bates 1366 duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym),
467 : maechler 1660 (tr) ? 0 : 1)));
468 : bates 1366 SET_VECTOR_ELT(dn, 1, duplicate(VECTOR_ELT(dn, 0)));
469 : mmaechler 2516 #ifdef AS_CHM_DIAGU2N_FIXED_FINALLY
470 : bates 1366 UNPROTECT(1);
471 : mmaechler 2516 #else
472 :     UNPROTECT(2);
473 :     #endif
474 : maechler 1548 return chm_sparse_to_SEXP(chcp, 1, 0, 0, "", dn);
475 : bates 922 }
476 : bates 923
477 : mmaechler 2646 /* Csparse_drop(x, tol): drop entries with absolute value < tol, i.e,
478 :     * at least all "explicit" zeros */
479 : maechler 1618 SEXP Csparse_drop(SEXP x, SEXP tol)
480 :     {
481 : mmaechler 2175 const char *cl = class_P(x);
482 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
483 :     int tr = (cl[1] == 't');
484 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
485 : mmaechler 2661 CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c);
486 : maechler 1618 double dtol = asReal(tol);
487 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
488 : maechler 1960 R_CheckStack();
489 : maechler 1618
490 : mmaechler 2661 if(!cholmod_drop(dtol, ans, &c))
491 :     error(_("cholmod_drop() failed"));
492 : mmaechler 2175 return chm_sparse_to_SEXP(ans, 1,
493 :     tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0,
494 :     Rkind, tr ? diag_P(x) : "",
495 : maechler 1736 GET_SLOT(x, Matrix_DimNamesSym));
496 : maechler 1618 }
497 :    
498 : bates 1218 SEXP Csparse_horzcat(SEXP x, SEXP y)
499 :     {
500 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
501 : mmaechler 2540 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
502 :     Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
503 :     Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
504 : maechler 1960 R_CheckStack();
505 : maechler 1375
506 : mmaechler 2540 /* TODO: currently drops dimnames - and we fix at R level */
507 : mmaechler 2661 return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c),
508 : maechler 1960 1, 0, Rkind, "", R_NilValue);
509 : bates 1218 }
510 :    
511 :     SEXP Csparse_vertcat(SEXP x, SEXP y)
512 :     {
513 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y);
514 : mmaechler 2540 int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0,
515 :     Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : 0,
516 :     Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0;
517 : maechler 1960 R_CheckStack();
518 : maechler 1375
519 : mmaechler 2540 /* TODO: currently drops dimnames - and we fix at R level */
520 : mmaechler 2661 return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c),
521 : maechler 1960 1, 0, Rkind, "", R_NilValue);
522 : bates 1218 }
523 : bates 1265
524 :     SEXP Csparse_band(SEXP x, SEXP k1, SEXP k2)
525 :     {
526 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
527 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
528 : mmaechler 2661 CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c);
529 : maechler 1960 R_CheckStack();
530 : bates 1265
531 : maechler 1736 return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "",
532 :     GET_SLOT(x, Matrix_DimNamesSym));
533 : bates 1265 }
534 : bates 1366
535 :     SEXP Csparse_diagU2N(SEXP x)
536 :     {
537 : maechler 2120 const char *cl = class_P(x);
538 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
539 :     if (cl[1] != 't' || *diag_P(x) != 'U') {
540 : maechler 2125 /* "trivially fast" when not triangular (<==> no 'diag' slot),
541 :     or not *unit* triangular */
542 : maechler 1708 return (x);
543 :     }
544 : maechler 2125 else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */
545 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(x);
546 : mmaechler 2661 CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c);
547 : maechler 1708 double one[] = {1, 0};
548 : mmaechler 2661 CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c);
549 : maechler 1710 int uploT = (*uplo_P(x) == 'U') ? 1 : -1;
550 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
551 : bates 1366
552 : maechler 1960 R_CheckStack();
553 : mmaechler 2661 cholmod_free_sparse(&eye, &c);
554 : maechler 1708 return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N",
555 : maechler 1736 GET_SLOT(x, Matrix_DimNamesSym));
556 : maechler 1708 }
557 : bates 1366 }
558 :    
559 : maechler 2125 SEXP Csparse_diagN2U(SEXP x)
560 :     {
561 :     const char *cl = class_P(x);
562 :     /* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */
563 :     if (cl[1] != 't' || *diag_P(x) != 'N') {
564 :     /* "trivially fast" when not triangular (<==> no 'diag' slot),
565 :     or already *unit* triangular */
566 :     return (x);
567 :     }
568 :     else { /* triangular with diag='N'): now drop the diagonal */
569 :     /* duplicate, since chx will be modified: */
570 : mmaechler 2223 CHM_SP chx = AS_CHM_SP__(duplicate(x));
571 : maechler 2125 int uploT = (*uplo_P(x) == 'U') ? 1 : -1,
572 :     Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
573 :     R_CheckStack();
574 :    
575 :     chm_diagN2U(chx, uploT, /* do_realloc */ FALSE);
576 :    
577 :     return chm_sparse_to_SEXP(chx, /*dofree*/ 0/* or 1 ?? */,
578 :     uploT, Rkind, "U",
579 :     GET_SLOT(x, Matrix_DimNamesSym));
580 :     }
581 :     }
582 :    
583 : mmaechler 2519 /**
584 :     * "Indexing" aka subsetting : Compute x[i,j], also for vectors i and j
585 :     * Working via CHOLMOD_submatrix, see ./CHOLMOD/MatrixOps/cholmod_submatrix.c
586 :     * @param x CsparseMatrix
587 :     * @param i row indices (0-origin), or NULL (R's)
588 :     * @param j columns indices (0-origin), or NULL
589 :     *
590 :     * @return x[i,j] still CsparseMatrix --- currently, this loses dimnames
591 :     */
592 : bates 1366 SEXP Csparse_submatrix(SEXP x, SEXP i, SEXP j)
593 :     {
594 : mmaechler 2519 CHM_SP chx = AS_CHM_SP(x); /* << does diagU2N() when needed */
595 : bates 1366 int rsize = (isNull(i)) ? -1 : LENGTH(i),
596 :     csize = (isNull(j)) ? -1 : LENGTH(j);
597 : maechler 1736 int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0;
598 : maechler 1960 R_CheckStack();
599 : bates 1366
600 :     if (rsize >= 0 && !isInteger(i))
601 :     error(_("Index i must be NULL or integer"));
602 :     if (csize >= 0 && !isInteger(j))
603 :     error(_("Index j must be NULL or integer"));
604 : maechler 1736
605 : mmaechler 2519 if (chx->stype) /* symmetricMatrix */
606 :     /* for now, cholmod_submatrix() only accepts "generalMatrix" */
607 : mmaechler 2661 chx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c);
608 : mmaechler 2519
609 : mmaechler 2661 return chm_sparse_to_SEXP(cholmod_submatrix(chx,
610 : mmaechler 2519 (rsize < 0) ? NULL : INTEGER(i), rsize,
611 :     (csize < 0) ? NULL : INTEGER(j), csize,
612 :     TRUE, TRUE, &c),
613 : maechler 1736 1, 0, Rkind, "",
614 :     /* FIXME: drops dimnames */ R_NilValue);
615 : bates 1366 }
616 : bates 2049
617 : mmaechler 2661 /**
618 :     * Subassignment: x[i,j] <- value
619 :     *
620 :     * @param x
621 :     * @param i_ integer row index 0-origin vector (as returned from R .ind.prep2())
622 :     * @param j_ integer column index 0-origin vector
623 :     * @param value currently must be a dsparseVector {which is recycled if needed}
624 :     *
625 :     * @return a Csparse matrix like x, but with the values replaced
626 :     */
627 :     SEXP Csparse_subassign(SEXP x, SEXP i_, SEXP j_, SEXP value)
628 :     {
629 : mmaechler 2677 // TODO: for other classes consider using a trick as RallocedReal() in ./chm_common.c
630 : mmaechler 2661 static const char
631 : mmaechler 2681 *valid_cM [] = { // the only ones, for "the moment". FIXME: extend (!)
632 :     "dgCMatrix",// 0
633 :     "dtCMatrix",// 1
634 :     ""},
635 :     // value: assume a "dsparseVector" for now -- slots: (i, length, x)
636 : mmaechler 2661 *valid_spv[] = {"dsparseVector",
637 :     ""};
638 :    
639 : mmaechler 2681 int ctype_x = Matrix_check_class_etc(x, valid_cM),
640 :     ctype_v = Matrix_check_class_etc(value, valid_spv);
641 :     if (ctype_x < 0)
642 : mmaechler 2661 error(_("invalid class of 'x' in Csparse_subassign()"));
643 : mmaechler 2681 if (ctype_v < 0)
644 : mmaechler 2661 error(_("invalid class of 'value' in Csparse_subassign()"));
645 :    
646 : mmaechler 2677 SEXP
647 :     islot = GET_SLOT(x, Matrix_iSym),
648 :     dimslot = GET_SLOT(x, Matrix_DimSym),
649 :     i_cp = PROTECT(coerceVector(i_, INTSXP)),
650 :     j_cp = PROTECT(coerceVector(j_, INTSXP));
651 : mmaechler 2673 // for d.CMatrix and l.CMatrix but not n.CMatrix:
652 : mmaechler 2661
653 : mmaechler 2677 int *dims = INTEGER(dimslot),
654 :     ncol = dims[1], /* nrow = dims[0], */
655 :     *i = INTEGER(i_cp), len_i = LENGTH(i_cp),
656 :     *j = INTEGER(j_cp), len_j = LENGTH(j_cp),
657 :     k,
658 :     nnz_x = LENGTH(islot);
659 :     int nnz = nnz_x;
660 :    
661 :     #define MATRIX_SUBASSIGN_VERBOSE
662 :     // Temporary hack for debugging --- remove eventually -- FIXME
663 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
664 :     Rboolean verbose = i[0] < 0;
665 :     if(verbose) i[0] = -i[0];
666 :     #endif
667 :    
668 :     SEXP val_i_slot;
669 :     PROTECT(val_i_slot = coerceVector(GET_SLOT(value, Matrix_iSym), REALSXP));
670 :     double *val_i = REAL(val_i_slot);
671 :     int nnz_val = LENGTH(GET_SLOT(value, Matrix_iSym));
672 : mmaechler 2673 // for dsparseVector only:
673 : mmaechler 2661 double *val_x = REAL (GET_SLOT(value, Matrix_xSym));
674 : mmaechler 2677 int64_t len_val = (int64_t) asReal(GET_SLOT(value, Matrix_lengthSym));
675 :     /* llen_i = (int64_t) len_i; */
676 : mmaechler 2661
677 : mmaechler 2677 double z_ans = 0.;
678 :    
679 :     SEXP ans;
680 :     /* Instead of simple "duplicate": PROTECT(ans = duplicate(x)) , build up: */
681 : mmaechler 2681 // Assuming that ans will have the same basic Matrix type as x :
682 :     ans = PROTECT(NEW_OBJECT(MAKE_CLASS(valid_cM[ctype_x])));
683 : mmaechler 2677 SET_SLOT(ans, Matrix_DimSym, duplicate(dimslot));
684 :     SET_SLOT(ans, Matrix_DimNamesSym, duplicate(GET_SLOT(x, Matrix_DimNamesSym)));
685 :     SET_SLOT(ans, Matrix_pSym, duplicate(GET_SLOT(x, Matrix_pSym)));
686 :     SEXP r_pslot = GET_SLOT(ans, Matrix_pSym);
687 :     // and assign the i- and x- slots at the end, as they are potentially modified
688 :     // not just in content, but also in their *length*
689 :     int *rp = INTEGER(r_pslot),
690 :     *ri = Calloc(nnz_x, int); // to contain the final i - slot
691 : mmaechler 2661 // for d.CMatrix only:
692 : mmaechler 2677 double *rx = Calloc(nnz_x, double); // to contain the final x - slot
693 :     Memcpy(ri, INTEGER(islot), nnz_x);
694 :     Memcpy(rx, REAL(GET_SLOT(x, Matrix_xSym)), nnz_x);
695 :     // NB: nnz_x : will always be the "current allocated length" of (i, x) slots
696 :     // -- nnz : the current *used* length; always nnz <= nnz_x
697 : mmaechler 2661
698 : mmaechler 2677 int jj, j_val = 0; // in "running" conceptionally through all value[i+ jj*len_i]
699 :     // values, we are "below"/"before" the (j_val)-th non-zero one.
700 :     // e.g. if value = (0,0,...,0), have nnz_val == 0, j_val must remain == 0
701 :     int64_t ii_val;// == "running" index (i + jj*len_i) % len_val for value[]
702 :     for(jj = 0, ii_val=0; jj < len_j; jj++) {
703 :     int j__ = j[jj];
704 :     /* int64_t j_l = jj * llen_i; */
705 :     R_CheckUserInterrupt();
706 :     for(int ii = 0; ii < len_i; ii++, ii_val++) {
707 :     int i__ = i[ii], p1, p2;
708 :     if(nnz_val && ii_val >= len_val) { // "recycle" indexing into value[]
709 :     ii_val -= len_val; // = (ii + jj*len_i) % len_val
710 :     j_val = 0;
711 :     }
712 :     int64_t ii_v1;//= ii_val + 1;
713 :     double v, /* := value[(ii + j_l) % len_val]
714 :     = dsparseVector_sub((ii + j_l) % len_val,
715 :     nnz_val, val_i, val_x, len_val)
716 :     */
717 :     M_ij;
718 :     int ind;
719 :     Rboolean have_entry = FALSE;
720 : mmaechler 2673
721 : mmaechler 2677 // note that rp[]'s may have *changed* even when 'j' remained!
722 :     // "FIXME": do this only *when* rp[] has changed
723 :     p1 = rp[j__], p2 = rp[j__ + 1];
724 : mmaechler 2673
725 : mmaechler 2677 // v := value[(ii + j_l) % len_val] = value[ii_val]
726 :     v = z_ans;
727 :     if(j_val < nnz_val) { // maybe find v := non-zero value[ii_val]
728 :     ii_v1 = ii_val + 1;
729 :     if(ii_v1 < val_i[j_val]) { // typical case: are still in zero-stretch
730 :     v = z_ans; // v = 0
731 :     } else if(ii_v1 == val_i[j_val]) { // have a match
732 :     v = val_x[j_val];
733 :     j_val++;// from now on, look at the next non-zero entry
734 :     } else { // ii_v1 > val_i[j_val]
735 :     REprintf("programming thinko in Csparse_subassign(*, i=%d,j=%d): ii_v=%d, v@i[j_val=%ld]=%g\n",
736 :     i__,j__, ii_v1, j_val, val_i[j_val]);
737 :     j_val++;// from now on, look at the next non-zero entry
738 :     }
739 :     }
740 :     // --------------- M_ij := getM(i., j.) --------------------------------
741 :     M_ij = z_ans; // as in ./t_sparseVector.c
742 :     for(ind = p1; ind < p2; ind++) {
743 :     if(ri[ind] >= i__) {
744 :     if(ri[ind] == i__) {
745 :     M_ij = rx[ind];
746 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
747 :     if(verbose) REprintf("have entry x[%d, %d] = %g\n",
748 :     i__, j__, M_ij);
749 :     #endif
750 :     have_entry = TRUE;
751 :     } else { // ri[ind] > i__
752 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
753 :     if(verbose)
754 :     REprintf("@i > i__ = %d --> ind-- = %d\n", i__, ind);
755 :     #endif
756 :     }
757 :     break;
758 :     }
759 :     }
760 : mmaechler 2673
761 : mmaechler 2677 //-- R: if(getM(i., j.) != (v <- getV(ii, jj)))
762 : mmaechler 2673
763 : mmaechler 2677 if(M_ij != v) { // contents differ ==> value needs to be changed
764 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
765 :     if(verbose)
766 :     REprintf("setting x[%d, %d] <- %g", i__,j__, v);
767 :     #endif
768 :     // (otherwise: nothing to do):
769 :     // setM(i__, j__, v)
770 :     // ----------------------------------------------------------
771 : mmaechler 2673
772 : mmaechler 2677 // Case I --------------------------------------------
773 :     /* if(v == z_ans) { // remove x[i, j] = M_ij which we know is *non*-zero */
774 : mmaechler 2681 //-------- Better (memory-management) *NOT* to remove, but rather at the very end
775 :     // currently using drop0() in R code
776 : mmaechler 2677 /* // we know : have_entry = TRUE ; */
777 :     /* // ri[ind] == i__; M_ij = rx[ind]; */
778 :     /* #ifdef MATRIX_SUBASSIGN_VERBOSE */
779 :     /* if(verbose) */
780 :     /* REprintf(" rm ind=%d\n", ind); */
781 :     /* #endif */
782 :     /* // remove the 'ind'-th element from x@i and x@x : */
783 :     /* nnz-- ; */
784 :     /* for(k=ind; k < nnz; k++) { */
785 :     /* ri[k] = ri[k+1]; */
786 :     /* rx[k] = rx[k+1]; */
787 :     /* } */
788 :     /* for(k=j__ + 1; k <= ncol; k++) { */
789 :     /* rp[k] = rp[k] - 1; */
790 :     /* } */
791 :     /* } */
792 :     /* else */
793 :     if(have_entry) {
794 :     // Case II ----- replace (non-empty) x[i,j] by v -------
795 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
796 :     if(verbose)
797 :     REprintf(" repl. ind=%d\n", ind);
798 :     #endif
799 :     rx[ind] = v;
800 :     } else {
801 :     // Case III ---- v != 0 : insert v into "empty" x[i,j] ----
802 : mmaechler 2673
803 : mmaechler 2677 // extend the i and x slot by one entry : ---------------------
804 : mmaechler 2673
805 : mmaechler 2677 if(nnz+1 > nnz_x) { // need to reallocate:
806 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
807 :     if(verbose) REprintf(" Realloc()ing: nnz_x=%d", nnz_x);
808 :     #endif
809 :     // do it "only" 1x,..4x at the very most increasing by the
810 :     // nnz-length of "value":
811 :     nnz_x += (1 + nnz_val / 4);
812 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
813 :     if(verbose) REprintf("(nnz_v=%d) --> %d ", nnz_val, nnz_x);
814 :     #endif
815 :     // C doc on realloc() says that the old content is *preserve*d
816 :     ri = Realloc(ri, nnz_x, int);
817 :     rx = Realloc(rx, nnz_x, double);
818 :     }
819 :    
820 :     // 3) fill them ...
821 :    
822 :     int i1 = ind;
823 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
824 :     if(verbose)
825 :     REprintf(" INSERT p12=(%d,%d) -> ind=%d -> i1 = %d\n",
826 :     p1,p2, ind, i1);
827 :     #endif
828 :    
829 :     // shift the "upper values" *before* the insertion:
830 :     for(int l = nnz-1; l >= i1; l--) {
831 :     ri[l+1] = ri[l];
832 :     rx[l+1] = rx[l];
833 :     }
834 :     ri[i1] = i__;
835 :     rx[i1] = v;
836 :     nnz++;
837 :    
838 :     // the columns j "right" of the current one :
839 :     for(k=j__ + 1; k <= ncol; k++)
840 :     rp[k]++;
841 :     }
842 :     }
843 :     #ifdef MATRIX_SUBASSIGN_VERBOSE
844 :     else if(verbose) REprintf("M_ij == v = %g\n", v);
845 :     #endif
846 :     }// for( ii )
847 :     }// for( jj )
848 :    
849 : mmaechler 2681 if(ctype_x == 1) { // triangularMatrix: copy the 'diag' and 'uplo' slots
850 :     SET_SLOT(ans, Matrix_uploSym, duplicate(GET_SLOT(x, Matrix_uploSym)));
851 :     SET_SLOT(ans, Matrix_diagSym, duplicate(GET_SLOT(x, Matrix_diagSym)));
852 :     }
853 : mmaechler 2677 // now assign the i- and x- slots, free memory and return :
854 :     Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), ri, nnz);
855 :     Memcpy( REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), rx, nnz);
856 :     Free(ri);
857 :     Free(rx);
858 :     UNPROTECT(4);
859 : mmaechler 2661 return ans;
860 :     }
861 :    
862 : bates 2049 SEXP Csparse_MatrixMarket(SEXP x, SEXP fname)
863 :     {
864 :     FILE *f = fopen(CHAR(asChar(fname)), "w");
865 :    
866 :     if (!f)
867 :     error(_("failure to open file \"%s\" for writing"),
868 :     CHAR(asChar(fname)));
869 : mmaechler 2661 if (!cholmod_write_sparse(f, AS_CHM_SP(x),
870 : maechler 2120 (CHM_SP)NULL, (char*) NULL, &c))
871 : mmaechler 2661 error(_("cholmod_write_sparse returned error code"));
872 : bates 2049 fclose(f);
873 :     return R_NilValue;
874 :     }
875 : maechler 2137
876 :    
877 :     /**
878 :     * Extract the diagonal entries from *triangular* Csparse matrix __or__ a
879 :     * cholmod_sparse factor (LDL = TRUE).
880 :     *
881 :     * @param n dimension of the matrix.
882 :     * @param x_p 'p' (column pointer) slot contents
883 :     * @param x_x 'x' (non-zero entries) slot contents
884 : mmaechler 2175 * @param perm 'perm' (= permutation vector) slot contents; only used for "diagBack"
885 : maechler 2137 * @param resultKind a (SEXP) string indicating which kind of result is desired.
886 :     *
887 :     * @return a SEXP, either a (double) number or a length n-vector of diagonal entries
888 :     */
889 :     SEXP diag_tC_ptr(int n, int *x_p, double *x_x, int *perm, SEXP resultKind)
890 :     /* ^^^^^^ FIXME[Generalize] to int / ... */
891 :     {
892 :     const char* res_ch = CHAR(STRING_ELT(resultKind,0));
893 :     enum diag_kind { diag, diag_backpermuted, trace, prod, sum_log
894 :     } res_kind = ((!strcmp(res_ch, "trace")) ? trace :
895 :     ((!strcmp(res_ch, "sumLog")) ? sum_log :
896 :     ((!strcmp(res_ch, "prod")) ? prod :
897 :     ((!strcmp(res_ch, "diag")) ? diag :
898 :     ((!strcmp(res_ch, "diagBack")) ? diag_backpermuted :
899 :     -1)))));
900 :     int i, n_x, i_from = 0;
901 :     SEXP ans = PROTECT(allocVector(REALSXP,
902 :     /* ^^^^ FIXME[Generalize] */
903 :     (res_kind == diag ||
904 :     res_kind == diag_backpermuted) ? n : 1));
905 :     double *v = REAL(ans);
906 :     /* ^^^^^^ ^^^^ FIXME[Generalize] */
907 :    
908 :     #define for_DIAG(v_ASSIGN) \
909 :     for(i = 0; i < n; i++, i_from += n_x) { \
910 :     /* looking at i-th column */ \
911 :     n_x = x_p[i+1] - x_p[i];/* #{entries} in this column */ \
912 :     v_ASSIGN; \
913 :     }
914 :    
915 :     /* NOTA BENE: we assume -- uplo = "L" i.e. lower triangular matrix
916 :     * for uplo = "U" (makes sense with a "dtCMatrix" !),
917 :     * should use x_x[i_from + (nx - 1)] instead of x_x[i_from],
918 :     * where nx = (x_p[i+1] - x_p[i])
919 :     */
920 :    
921 :     switch(res_kind) {
922 :     case trace:
923 :     v[0] = 0.;
924 :     for_DIAG(v[0] += x_x[i_from]);
925 :     break;
926 :    
927 :     case sum_log:
928 :     v[0] = 0.;
929 :     for_DIAG(v[0] += log(x_x[i_from]));
930 :     break;
931 :    
932 :     case prod:
933 :     v[0] = 1.;
934 :     for_DIAG(v[0] *= x_x[i_from]);
935 :     break;
936 :    
937 :     case diag:
938 :     for_DIAG(v[i] = x_x[i_from]);
939 :     break;
940 :    
941 :     case diag_backpermuted:
942 :     for_DIAG(v[i] = x_x[i_from]);
943 :    
944 : mmaechler 2175 warning(_("resultKind = 'diagBack' (back-permuted) is experimental"));
945 : maechler 2144 /* now back_permute : */
946 :     for(i = 0; i < n; i++) {
947 :     double tmp = v[i]; v[i] = v[perm[i]]; v[perm[i]] = tmp;
948 :     /*^^^^ FIXME[Generalize] */
949 :     }
950 : maechler 2137 break;
951 :    
952 :     default: /* -1 from above */
953 : mmaechler 2387 error(_("diag_tC(): invalid 'resultKind'"));
954 : maechler 2137 /* Wall: */ ans = R_NilValue; v = REAL(ans);
955 :     }
956 :    
957 :     UNPROTECT(1);
958 :     return ans;
959 :     }
960 :    
961 :     /**
962 :     * Extract the diagonal entries from *triangular* Csparse matrix __or__ a
963 :     * cholmod_sparse factor (LDL = TRUE).
964 :     *
965 :     * @param pslot 'p' (column pointer) slot of Csparse matrix/factor
966 :     * @param xslot 'x' (non-zero entries) slot of Csparse matrix/factor
967 : mmaechler 2175 * @param perm_slot 'perm' (= permutation vector) slot of corresponding CHMfactor;
968 :     * only used for "diagBack"
969 : maechler 2137 * @param resultKind a (SEXP) string indicating which kind of result is desired.
970 :     *
971 :     * @return a SEXP, either a (double) number or a length n-vector of diagonal entries
972 :     */
973 :     SEXP diag_tC(SEXP pslot, SEXP xslot, SEXP perm_slot, SEXP resultKind)
974 :     {
975 :     int n = length(pslot) - 1, /* n = ncol(.) = nrow(.) */
976 :     *x_p = INTEGER(pslot),
977 :     *perm = INTEGER(perm_slot);
978 :     double *x_x = REAL(xslot);
979 :     /* ^^^^^^ ^^^^ FIXME[Generalize] to INTEGER(.) / LOGICAL(.) / ... xslot !*/
980 :    
981 :     return diag_tC_ptr(n, x_p, x_x, perm, resultKind);
982 :     }
983 : dmbates 2319
984 :     /**
985 :     * Create a Csparse matrix object from indices and/or pointers.
986 :     *
987 :     * @param cls name of actual class of object to create
988 :     * @param i optional integer vector of length nnz of row indices
989 :     * @param j optional integer vector of length nnz of column indices
990 :     * @param p optional integer vector of length np of row or column pointers
991 :     * @param np length of integer vector p. Must be zero if p == (int*)NULL
992 :     * @param x optional vector of values
993 :     * @param nnz length of vectors i, j and/or x, whichever is to be used
994 :     * @param dims optional integer vector of length 2 to be used as
995 :     * dimensions. If dims == (int*)NULL then the maximum row and column
996 :     * index are used as the dimensions.
997 :     * @param dimnames optional list of length 2 to be used as dimnames
998 :     * @param index1 indicator of 1-based indices
999 :     *
1000 :     * @return an SEXP of class cls inheriting from CsparseMatrix.
1001 :     */
1002 :     SEXP create_Csparse(char* cls, int* i, int* j, int* p, int np,
1003 :     void* x, int nnz, int* dims, SEXP dimnames,
1004 :     int index1)
1005 :     {
1006 :     SEXP ans;
1007 :     int *ij = (int*)NULL, *tri, *trj,
1008 :     mi, mj, mp, nrow = -1, ncol = -1;
1009 :     int xtype = -1; /* -Wall */
1010 :     CHM_TR T;
1011 :     CHM_SP A;
1012 :    
1013 :     if (np < 0 || nnz < 0)
1014 :     error(_("negative vector lengths not allowed: np = %d, nnz = %d"),
1015 :     np, nnz);
1016 :     if (1 != ((mi = (i == (int*)NULL)) +
1017 :     (mj = (j == (int*)NULL)) +
1018 :     (mp = (p == (int*)NULL))))
1019 :     error(_("exactly 1 of 'i', 'j' or 'p' must be NULL"));
1020 :     if (mp) {
1021 :     if (np) error(_("np = %d, must be zero when p is NULL"), np);
1022 :     } else {
1023 :     if (np) { /* Expand p to form i or j */
1024 :     if (!(p[0])) error(_("p[0] = %d, should be zero"), p[0]);
1025 :     for (int ii = 0; ii < np; ii++)
1026 :     if (p[ii] > p[ii + 1])
1027 :     error(_("p must be non-decreasing"));
1028 :     if (p[np] != nnz)
1029 : mmaechler 2387 error("p[np] = %d != nnz = %d", p[np], nnz);
1030 : dmbates 2319 ij = Calloc(nnz, int);
1031 :     if (mi) {
1032 :     i = ij;
1033 :     nrow = np;
1034 :     } else {
1035 :     j = ij;
1036 :     ncol = np;
1037 :     }
1038 : mmaechler 2661 /* Expand p to 0-based indices */
1039 : dmbates 2319 for (int ii = 0; ii < np; ii++)
1040 :     for (int jj = p[ii]; jj < p[ii + 1]; jj++) ij[jj] = ii;
1041 :     } else {
1042 :     if (nnz)
1043 :     error(_("Inconsistent dimensions: np = 0 and nnz = %d"),
1044 :     nnz);
1045 :     }
1046 :     }
1047 : mmaechler 2661 /* calculate nrow and ncol */
1048 : dmbates 2319 if (nrow < 0) {
1049 :     for (int ii = 0; ii < nnz; ii++) {
1050 :     int i1 = i[ii] + (index1 ? 0 : 1); /* 1-based index */
1051 :     if (i1 < 1) error(_("invalid row index at position %d"), ii);
1052 :     if (i1 > nrow) nrow = i1;
1053 :     }
1054 :     }
1055 :     if (ncol < 0) {
1056 :     for (int jj = 0; jj < nnz; jj++) {
1057 :     int j1 = j[jj] + (index1 ? 0 : 1);
1058 : mmaechler 2387 if (j1 < 1) error(_("invalid column index at position %d"), jj);
1059 : dmbates 2319 if (j1 > ncol) ncol = j1;
1060 :     }
1061 :     }
1062 :     if (dims != (int*)NULL) {
1063 :     if (dims[0] > nrow) nrow = dims[0];
1064 :     if (dims[1] > ncol) ncol = dims[1];
1065 :     }
1066 : mmaechler 2661 /* check the class name */
1067 : dmbates 2319 if (strlen(cls) != 8)
1068 :     error(_("strlen of cls argument = %d, should be 8"), strlen(cls));
1069 :     if (!strcmp(cls + 2, "CMatrix"))
1070 :     error(_("cls = \"%s\" does not end in \"CMatrix\""), cls);
1071 :     switch(cls[0]) {
1072 :     case 'd':
1073 :     case 'l':
1074 : mmaechler 2661 xtype = CHOLMOD_REAL;
1075 :     break;
1076 : dmbates 2319 case 'n':
1077 : mmaechler 2661 xtype = CHOLMOD_PATTERN;
1078 :     break;
1079 : dmbates 2319 default:
1080 : mmaechler 2661 error(_("cls = \"%s\" must begin with 'd', 'l' or 'n'"), cls);
1081 : dmbates 2319 }
1082 :     if (cls[1] != 'g')
1083 :     error(_("Only 'g'eneral sparse matrix types allowed"));
1084 : mmaechler 2661 /* allocate and populate the triplet */
1085 :     T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0,
1086 :     xtype, &c);
1087 : dmbates 2319 T->x = x;
1088 :     tri = (int*)T->i;
1089 :     trj = (int*)T->j;
1090 :     for (int ii = 0; ii < nnz; ii++) {
1091 :     tri[ii] = i[ii] - ((!mi && index1) ? 1 : 0);
1092 :     trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0);
1093 :     }
1094 : mmaechler 2661 /* create the cholmod_sparse structure */
1095 :     A = cholmod_triplet_to_sparse(T, nnz, &c);
1096 :     cholmod_free_triplet(&T, &c);
1097 :     /* copy the information to the SEXP */
1098 : dmbates 2319 ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls)));
1099 :     /* FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c */
1100 : mmaechler 2661 /* allocate and copy common slots */
1101 :     nnz = cholmod_nnz(A, &c);
1102 : dmbates 2319 dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2));
1103 :     dims[0] = A->nrow; dims[1] = A->ncol;
1104 :     Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1);
1105 :     Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), (int*)A->i, nnz);
1106 :     switch(cls[1]) {
1107 :     case 'd':
1108 :     Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), (double*)A->x, nnz);
1109 :     break;
1110 :     case 'l':
1111 :     error(_("code not yet written for cls = \"lgCMatrix\""));
1112 :     }
1113 : mmaechler 2646 /* FIXME: dimnames are *NOT* put there yet (if non-NULL) */
1114 : mmaechler 2661 cholmod_free_sparse(&A, &c);
1115 : dmbates 2319 UNPROTECT(1);
1116 :     return ans;
1117 :     }

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