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

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