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