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

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