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

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

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