1 |
/* Sparse matrices in compress column-oriented form */ |
/* Sparse matrices in compressed column-oriented form */ |
2 |
|
|
3 |
#include "Csparse.h" |
#include "Csparse.h" |
4 |
#ifdef USE_CHOLMOD |
#include "Tsparse.h" |
5 |
#include "chm_common.h" |
#include "chm_common.h" |
|
#endif /* USE_CHOLMOD */ |
|
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 */ |
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, ncol = length(pslot) - 1, |
int *dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), j, |
14 |
*dims = INTEGER(GET_SLOT(x, Matrix_DimSym)), |
nrow = dims[0], |
15 |
nrow, *xp = INTEGER(pslot), |
ncol = dims[1], |
16 |
|
*xp = INTEGER(pslot), |
17 |
*xi = INTEGER(islot); |
*xi = INTEGER(islot); |
18 |
|
|
19 |
nrow = dims[0]; |
if (length(pslot) != dims[1] + 1) |
20 |
if (length(pslot) <= 0) |
return FALSE; |
|
return mkString(_("slot p must have length > 0")); |
|
21 |
if (xp[0] != 0) |
if (xp[0] != 0) |
22 |
return mkString(_("first element of slot p must be zero")); |
return FALSE; |
23 |
if (length(islot) != xp[ncol]) |
if (length(islot) < xp[ncol]) /* allow larger slots from over-allocation!*/ |
24 |
return mkString(_("last element of slot p must match length of slots i and x")); |
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++) { |
for (j = 0; j < ncol; j++) { |
30 |
if (xp[j] > xp[j+1]) |
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 |
|
|
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 |
|
SEXP Csparse_validate2(SEXP x, SEXP maybe_modify) { |
49 |
|
return Csparse_validate_(x, asLogical(maybe_modify)); |
50 |
|
} |
51 |
|
|
52 |
|
// 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 |
|
} |
58 |
|
|
59 |
|
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")); |
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 |
} |
} |
|
for (j = 0; j < length(islot); j++) { |
|
|
if (xi[j] < 0 || xi[j] >= nrow) |
|
|
return mkString(_("all row indices must be between 0 and nrow-1")); |
|
94 |
} |
} |
95 |
|
if (!sorted) |
96 |
|
/* cannot easily use cholmod_sort(.) ... -> "error out" :*/ |
97 |
|
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); |
return ScalarLogical(1); |
102 |
} |
} |
103 |
|
|
104 |
SEXP Csparse_to_Tsparse(SEXP x) |
/** |
105 |
|
* From a CsparseMatrix, produce a dense one. |
106 |
|
* Directly deals with symmetric, triangular and general. |
107 |
|
* Called from ../R/Csparse.R's C2dense() |
108 |
|
* |
109 |
|
* @param x a CsparseMatrix: currently all 9 of "[dln][gst]CMatrix" |
110 |
|
* @param symm_or_tri integer (NA, < 0, > 0, = 0) specifying the knowledge of the caller about x: |
111 |
|
* NA : unknown => will be determined |
112 |
|
* = 0 : "generalMatrix" (not symm or tri); |
113 |
|
* < 0 : "triangularMatrix" |
114 |
|
* > 0 : "symmetricMatrix" |
115 |
|
* |
116 |
|
* @return a "denseMatrix" |
117 |
|
*/ |
118 |
|
SEXP Csparse_to_dense(SEXP x, SEXP symm_or_tri) |
119 |
{ |
{ |
120 |
#ifdef USE_CHOLMOD |
Rboolean is_sym, is_tri; |
121 |
cholmod_sparse *chxs = as_cholmod_sparse(x); |
int is_sym_or_tri = asInteger(symm_or_tri), |
122 |
cholmod_triplet *chxt = cholmod_sparse_to_triplet(chxs, &c); |
ctype = 0; // <- default = "dgC" |
123 |
|
static const char *valid[] = { MATRIX_VALID_Csparse, ""}; |
124 |
|
if(is_sym_or_tri == NA_INTEGER) { // find if is(x, "symmetricMatrix") : |
125 |
|
ctype = Matrix_check_class_etc(x, valid); |
126 |
|
is_sym = (ctype % 3 == 1); |
127 |
|
is_tri = (ctype % 3 == 2); |
128 |
|
} else { |
129 |
|
is_sym = is_sym_or_tri > 0; |
130 |
|
is_tri = is_sym_or_tri < 0; |
131 |
|
// => both are FALSE iff is_.. == 0 |
132 |
|
if(is_sym || is_tri) |
133 |
|
ctype = Matrix_check_class_etc(x, valid); |
134 |
|
} |
135 |
|
CHM_SP chxs = AS_CHM_SP__(x);// -> chxs->stype = +- 1 <==> symmetric |
136 |
|
R_CheckStack(); |
137 |
|
if(is_tri && *diag_P(x) == 'U') { // ==> x := diagU2N(x), directly for chxs |
138 |
|
CHM_SP eye = cholmod_speye(chxs->nrow, chxs->ncol, chxs->xtype, &c); |
139 |
|
double one[] = {1, 0}; |
140 |
|
CHM_SP ans = cholmod_add(chxs, eye, one, one, |
141 |
|
/* values: */ ((ctype / 3) != 2), // TRUE iff not "nMatrix" |
142 |
|
TRUE, &c); |
143 |
|
cholmod_free_sparse(&eye, &c); |
144 |
|
chxs = cholmod_copy_sparse(ans, &c); |
145 |
|
cholmod_free_sparse(&ans, &c); |
146 |
|
} |
147 |
|
/* The following loses the symmetry property, since cholmod_dense has none, |
148 |
|
* BUT, much worse (FIXME!), it also transforms CHOLMOD_PATTERN ("n") matrices |
149 |
|
* to numeric (CHOLMOD_REAL) ones {and we "revert" via chm_dense_to_SEXP()}: */ |
150 |
|
CHM_DN chxd = cholmod_sparse_to_dense(chxs, &c); |
151 |
|
int Rkind = (chxs->xtype == CHOLMOD_PATTERN)? -1 : Real_kind(x); |
152 |
|
|
153 |
free(chxs); |
SEXP ans = chm_dense_to_SEXP(chxd, 1, Rkind,GET_SLOT(x, Matrix_DimNamesSym)); |
154 |
return chm_triplet_to_SEXP(chxt, 1); |
// -> a [dln]geMatrix |
155 |
#else |
if(is_sym) { // ==> want [dln]syMatrix |
156 |
error("General conversion requires CHOLMOD"); |
const char cl1 = class_P(ans)[0]; |
157 |
return R_NilValue; /* -Wall */ |
PROTECT(ans); |
158 |
#endif /* USE_CHOLMOD */ |
SEXP aa = PROTECT(NEW_OBJECT(MAKE_CLASS((cl1 == 'd') ? "dsyMatrix" : |
159 |
|
((cl1 == 'l') ? "lsyMatrix" : "nsyMatrix")))); |
160 |
|
// No need to duplicate() as slots of ans are freshly allocated and ans will not be used |
161 |
|
SET_SLOT(aa, Matrix_xSym, GET_SLOT(ans, Matrix_xSym)); |
162 |
|
SET_SLOT(aa, Matrix_DimSym, GET_SLOT(ans, Matrix_DimSym)); |
163 |
|
SET_SLOT(aa, Matrix_DimNamesSym,GET_SLOT(ans, Matrix_DimNamesSym)); |
164 |
|
SET_SLOT(aa, Matrix_uploSym, mkString((chxs->stype > 0) ? "U" : "L")); |
165 |
|
UNPROTECT(2); |
166 |
|
return aa; |
167 |
|
} |
168 |
|
else if(is_tri) { // ==> want [dln]trMatrix |
169 |
|
const char cl1 = class_P(ans)[0]; |
170 |
|
PROTECT(ans); |
171 |
|
SEXP aa = PROTECT(NEW_OBJECT(MAKE_CLASS((cl1 == 'd') ? "dtrMatrix" : |
172 |
|
((cl1 == 'l') ? "ltrMatrix" : "ntrMatrix")))); |
173 |
|
// No need to duplicate() as slots of ans are freshly allocated and ans will not be used |
174 |
|
SET_SLOT(aa, Matrix_xSym, GET_SLOT(ans, Matrix_xSym)); |
175 |
|
SET_SLOT(aa, Matrix_DimSym, GET_SLOT(ans, Matrix_DimSym)); |
176 |
|
SET_SLOT(aa, Matrix_DimNamesSym,GET_SLOT(ans, Matrix_DimNamesSym)); |
177 |
|
slot_dup(aa, x, Matrix_uploSym); |
178 |
|
/* already by NEW_OBJECT(..) above: |
179 |
|
SET_SLOT(aa, Matrix_diagSym, mkString("N")); */ |
180 |
|
UNPROTECT(2); |
181 |
|
return aa; |
182 |
|
} |
183 |
|
else |
184 |
|
return ans; |
185 |
} |
} |
186 |
|
|
187 |
SEXP Csparse_transpose(SEXP x) |
// FIXME: do not go via CHM (should not be too hard, to just *drop* the x-slot, right? |
188 |
|
SEXP Csparse_to_nz_pattern(SEXP x, SEXP tri) |
189 |
{ |
{ |
190 |
#ifdef USE_CHOLMOD |
CHM_SP chxs = AS_CHM_SP__(x); |
191 |
cholmod_sparse *chx = as_cholmod_sparse(x); |
CHM_SP chxcp = cholmod_copy(chxs, chxs->stype, CHOLMOD_PATTERN, &c); |
192 |
cholmod_sparse *chxt = cholmod_transpose(chx, (int) chx->xtype, &c); |
int tr = asLogical(tri); |
193 |
|
R_CheckStack(); |
194 |
|
|
195 |
free(chx); |
return chm_sparse_to_SEXP(chxcp, 1/*do_free*/, |
196 |
return chm_sparse_to_SEXP(chxt, 1); |
tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, |
197 |
#else |
0, tr ? diag_P(x) : "", |
198 |
error("General conversion requires CHOLMOD"); |
GET_SLOT(x, Matrix_DimNamesSym)); |
199 |
|
} |
200 |
|
|
201 |
|
// n.CMatrix --> [dli].CMatrix (not going through CHM!) |
202 |
|
SEXP nz_pattern_to_Csparse(SEXP x, SEXP res_kind) |
203 |
|
{ |
204 |
|
return nz2Csparse(x, asInteger(res_kind)); |
205 |
|
} |
206 |
|
|
207 |
|
// n.CMatrix --> [dli].CMatrix (not going through CHM!) |
208 |
|
// NOTE: use chm_MOD_xtype(() to change type of 'cholmod_sparse' matrix |
209 |
|
SEXP nz2Csparse(SEXP x, enum x_slot_kind r_kind) |
210 |
|
{ |
211 |
|
const char *cl_x = class_P(x); |
212 |
|
if(cl_x[0] != 'n') error(_("not a 'n.CMatrix'")); |
213 |
|
if(cl_x[2] != 'C') error(_("not a CsparseMatrix")); |
214 |
|
int nnz = LENGTH(GET_SLOT(x, Matrix_iSym)); |
215 |
|
SEXP ans; |
216 |
|
char *ncl = alloca(strlen(cl_x) + 1); /* not much memory required */ |
217 |
|
strcpy(ncl, cl_x); |
218 |
|
double *dx_x; int *ix_x; |
219 |
|
ncl[0] = (r_kind == x_double ? 'd' : |
220 |
|
(r_kind == x_logical ? 'l' : |
221 |
|
/* else (for now): r_kind == x_integer : */ 'i')); |
222 |
|
PROTECT(ans = NEW_OBJECT(MAKE_CLASS(ncl))); |
223 |
|
// create a correct 'x' slot: |
224 |
|
switch(r_kind) { |
225 |
|
int i; |
226 |
|
case x_double: // 'd' |
227 |
|
dx_x = REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)); |
228 |
|
for (i=0; i < nnz; i++) dx_x[i] = 1.; |
229 |
|
break; |
230 |
|
case x_logical: // 'l' |
231 |
|
ix_x = LOGICAL(ALLOC_SLOT(ans, Matrix_xSym, LGLSXP, nnz)); |
232 |
|
for (i=0; i < nnz; i++) ix_x[i] = TRUE; |
233 |
|
break; |
234 |
|
case x_integer: // 'i' |
235 |
|
ix_x = INTEGER(ALLOC_SLOT(ans, Matrix_xSym, INTSXP, nnz)); |
236 |
|
for (i=0; i < nnz; i++) ix_x[i] = 1; |
237 |
|
break; |
238 |
|
|
239 |
|
default: |
240 |
|
error(_("nz2Csparse(): invalid/non-implemented r_kind = %d"), |
241 |
|
r_kind); |
242 |
|
} |
243 |
|
|
244 |
|
// now copy all other slots : |
245 |
|
slot_dup(ans, x, Matrix_iSym); |
246 |
|
slot_dup(ans, x, Matrix_pSym); |
247 |
|
slot_dup(ans, x, Matrix_DimSym); |
248 |
|
slot_dup(ans, x, Matrix_DimNamesSym); |
249 |
|
if(ncl[1] != 'g') { // symmetric or triangular ... |
250 |
|
slot_dup_if_has(ans, x, Matrix_uploSym); |
251 |
|
slot_dup_if_has(ans, x, Matrix_diagSym); |
252 |
|
} |
253 |
|
UNPROTECT(1); |
254 |
|
return ans; |
255 |
|
} |
256 |
|
|
257 |
|
SEXP Csparse_to_matrix(SEXP x, SEXP chk, SEXP symm) |
258 |
|
{ |
259 |
|
int is_sym = asLogical(symm); |
260 |
|
if(is_sym == NA_LOGICAL) { // find if is(x, "symmetricMatrix") : |
261 |
|
static const char *valid[] = { MATRIX_VALID_Csparse, ""}; |
262 |
|
int ctype = Matrix_check_class_etc(x, valid); |
263 |
|
is_sym = (ctype % 3 == 1); |
264 |
|
} |
265 |
|
return chm_dense_to_matrix( |
266 |
|
cholmod_sparse_to_dense(AS_CHM_SP2(x, asLogical(chk)), &c), |
267 |
|
1 /*do_free*/, |
268 |
|
(is_sym |
269 |
|
? symmetric_DimNames(GET_SLOT(x, Matrix_DimNamesSym)) |
270 |
|
: GET_SLOT(x, Matrix_DimNamesSym))); |
271 |
|
} |
272 |
|
|
273 |
|
SEXP Csparse_to_vector(SEXP x) |
274 |
|
{ |
275 |
|
return chm_dense_to_vector(cholmod_sparse_to_dense(AS_CHM_SP__(x), &c), 1); |
276 |
|
} |
277 |
|
|
278 |
|
SEXP Csparse_to_Tsparse(SEXP x, SEXP tri) |
279 |
|
{ |
280 |
|
CHM_SP chxs = AS_CHM_SP__(x); |
281 |
|
CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c); |
282 |
|
int tr = asLogical(tri); |
283 |
|
int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
284 |
|
R_CheckStack(); |
285 |
|
|
286 |
|
return chm_triplet_to_SEXP(chxt, 1, |
287 |
|
tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, |
288 |
|
Rkind, tr ? diag_P(x) : "", |
289 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
290 |
|
} |
291 |
|
|
292 |
|
SEXP Csparse_to_tCsparse(SEXP x, SEXP uplo, SEXP diag) |
293 |
|
{ |
294 |
|
CHM_SP chxs = AS_CHM_SP__(x); |
295 |
|
int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
296 |
|
R_CheckStack(); |
297 |
|
return chm_sparse_to_SEXP(chxs, /* dofree = */ 0, |
298 |
|
/* uploT = */ (*CHAR(asChar(uplo)) == 'U')? 1: -1, |
299 |
|
Rkind, /* diag = */ CHAR(STRING_ELT(diag, 0)), |
300 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
301 |
|
} |
302 |
|
|
303 |
|
SEXP Csparse_to_tTsparse(SEXP x, SEXP uplo, SEXP diag) |
304 |
|
{ |
305 |
|
CHM_SP chxs = AS_CHM_SP__(x); |
306 |
|
CHM_TR chxt = cholmod_sparse_to_triplet(chxs, &c); |
307 |
|
int Rkind = (chxs->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
308 |
|
R_CheckStack(); |
309 |
|
return chm_triplet_to_SEXP(chxt, 1, |
310 |
|
/* uploT = */ (*CHAR(asChar(uplo)) == 'U')? 1: -1, |
311 |
|
Rkind, /* diag = */ CHAR(STRING_ELT(diag, 0)), |
312 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
313 |
|
} |
314 |
|
|
315 |
|
|
316 |
|
SEXP Csparse_symmetric_to_general(SEXP x) |
317 |
|
{ |
318 |
|
CHM_SP chx = AS_CHM_SP__(x), chgx; |
319 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
320 |
|
R_CheckStack(); |
321 |
|
|
322 |
|
if (!(chx->stype)) |
323 |
|
error(_("Nonsymmetric matrix in Csparse_symmetric_to_general")); |
324 |
|
chgx = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c); |
325 |
|
/* xtype: pattern, "real", complex or .. */ |
326 |
|
return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "", |
327 |
|
symmetric_DimNames(GET_SLOT(x, Matrix_DimNamesSym))); |
328 |
|
} |
329 |
|
|
330 |
|
SEXP Csparse_general_to_symmetric(SEXP x, SEXP uplo, SEXP sym_dmns) |
331 |
|
{ |
332 |
|
int *adims = INTEGER(GET_SLOT(x, Matrix_DimSym)), n = adims[0]; |
333 |
|
if(n != adims[1]) { |
334 |
|
error(_("Csparse_general_to_symmetric(): matrix is not square!")); |
335 |
return R_NilValue; /* -Wall */ |
return R_NilValue; /* -Wall */ |
|
#endif /* USE_CHOLMOD */ |
|
336 |
} |
} |
337 |
|
CHM_SP chx = AS_CHM_SP__(x), chgx; |
338 |
|
int uploT = (*CHAR(asChar(uplo)) == 'U') ? 1 : -1; |
339 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
340 |
|
R_CheckStack(); |
341 |
|
chgx = cholmod_copy(chx, /* stype: */ uploT, chx->xtype, &c); |
342 |
|
|
343 |
|
SEXP dns = GET_SLOT(x, Matrix_DimNamesSym); |
344 |
|
if(asLogical(sym_dmns)) |
345 |
|
dns = symmetric_DimNames(dns); |
346 |
|
else if((!isNull(VECTOR_ELT(dns, 0)) && |
347 |
|
!isNull(VECTOR_ELT(dns, 1))) || |
348 |
|
!isNull(getAttrib(dns, R_NamesSymbol))) { |
349 |
|
/* symmetrize them if both are not NULL |
350 |
|
* or names(dimnames(.)) is asymmetric : */ |
351 |
|
dns = PROTECT(duplicate(dns)); |
352 |
|
if(!equal_string_vectors(VECTOR_ELT(dns, 0), |
353 |
|
VECTOR_ELT(dns, 1))) { |
354 |
|
if(uploT == 1) |
355 |
|
SET_VECTOR_ELT(dns, 0, VECTOR_ELT(dns,1)); |
356 |
|
else |
357 |
|
SET_VECTOR_ELT(dns, 1, VECTOR_ELT(dns,0)); |
358 |
|
} |
359 |
|
SEXP nms_dns = getAttrib(dns, R_NamesSymbol); |
360 |
|
if(!isNull(nms_dns) && // names(dimnames(.)) : |
361 |
|
!R_compute_identical(STRING_ELT(nms_dns, 0), |
362 |
|
STRING_ELT(nms_dns, 1), 16)) { |
363 |
|
if(uploT == 1) |
364 |
|
SET_STRING_ELT(nms_dns, 0, STRING_ELT(nms_dns,1)); |
365 |
|
else |
366 |
|
SET_STRING_ELT(nms_dns, 1, STRING_ELT(nms_dns,0)); |
367 |
|
setAttrib(dns, R_NamesSymbol, nms_dns); |
368 |
|
} |
369 |
|
UNPROTECT(1); |
370 |
|
} |
371 |
|
/* xtype: pattern, "real", complex or .. */ |
372 |
|
return chm_sparse_to_SEXP(chgx, 1, 0, Rkind, "", dns); |
373 |
|
} |
374 |
|
|
375 |
|
SEXP Csparse_transpose(SEXP x, SEXP tri) |
376 |
|
{ |
377 |
|
/* TODO: lgCMatrix & igC* currently go via double prec. cholmod - |
378 |
|
* since cholmod (& cs) lacks sparse 'int' matrices */ |
379 |
|
CHM_SP chx = AS_CHM_SP__(x); |
380 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
381 |
|
CHM_SP chxt = cholmod_transpose(chx, chx->xtype, &c); |
382 |
|
SEXP dn = PROTECT(duplicate(GET_SLOT(x, Matrix_DimNamesSym))), tmp; |
383 |
|
int tr = asLogical(tri); |
384 |
|
R_CheckStack(); |
385 |
|
|
386 |
|
tmp = VECTOR_ELT(dn, 0); /* swap the dimnames */ |
387 |
|
SET_VECTOR_ELT(dn, 0, VECTOR_ELT(dn, 1)); |
388 |
|
SET_VECTOR_ELT(dn, 1, tmp); |
389 |
|
if(!isNull(tmp = getAttrib(dn, R_NamesSymbol))) { // swap names(dimnames(.)): |
390 |
|
SEXP nms_dns = PROTECT(allocVector(VECSXP, 2)); |
391 |
|
SET_VECTOR_ELT(nms_dns, 1, STRING_ELT(tmp, 0)); |
392 |
|
SET_VECTOR_ELT(nms_dns, 0, STRING_ELT(tmp, 1)); |
393 |
|
setAttrib(dn, R_NamesSymbol, nms_dns); |
394 |
|
UNPROTECT(1); |
395 |
|
} |
396 |
|
UNPROTECT(1); |
397 |
|
return chm_sparse_to_SEXP(chxt, 1, /* SWAP 'uplo' for triangular */ |
398 |
|
tr ? ((*uplo_P(x) == 'U') ? -1 : 1) : 0, |
399 |
|
Rkind, tr ? diag_P(x) : "", dn); |
400 |
|
} |
401 |
|
|
402 |
SEXP Csparse_Csparse_prod(SEXP a, SEXP b) |
SEXP Csparse_Csparse_prod(SEXP a, SEXP b) |
403 |
{ |
{ |
404 |
#ifdef USE_CHOLMOD |
CHM_SP |
405 |
cholmod_sparse *cha = as_cholmod_sparse(a), *chb = as_cholmod_sparse(b); |
cha = AS_CHM_SP(a), |
406 |
cholmod_sparse *chc = cholmod_ssmult(cha, chb, 0, (int) cha->xtype, 1, &c); |
chb = AS_CHM_SP(b), |
407 |
|
chc = cholmod_ssmult(cha, chb, /*out_stype:*/ 0, |
408 |
|
/* values:= is_numeric (T/F) */ cha->xtype > 0, |
409 |
|
/*out sorted:*/ 1, &c); |
410 |
|
const char *cl_a = class_P(a), *cl_b = class_P(b); |
411 |
|
char diag[] = {'\0', '\0'}; |
412 |
|
int uploT = 0; |
413 |
|
SEXP dn = PROTECT(allocVector(VECSXP, 2)); |
414 |
|
R_CheckStack(); |
415 |
|
|
416 |
free(cha); free(chb); |
#ifdef DEBUG_Matrix_verbose |
417 |
return chm_sparse_to_SEXP(chc, 1); |
Rprintf("DBG Csparse_C*_prod(%s, %s)\n", cl_a, cl_b); |
418 |
#else |
#endif |
419 |
error("General multiplication requires CHOLMOD"); |
|
420 |
return R_NilValue; /* -Wall */ |
/* Preserve triangularity and even unit-triangularity if appropriate. |
421 |
#endif /* USE_CHOLMOD */ |
* Note that in that case, the multiplication itself should happen |
422 |
|
* faster. But there's no support for that in CHOLMOD */ |
423 |
|
|
424 |
|
/* UGLY hack -- rather should have (fast!) C-level version of |
425 |
|
* is(a, "triangularMatrix") etc */ |
426 |
|
if (cl_a[1] == 't' && cl_b[1] == 't') |
427 |
|
/* FIXME: fails for "Cholesky","BunchKaufmann"..*/ |
428 |
|
if(*uplo_P(a) == *uplo_P(b)) { /* both upper, or both lower tri. */ |
429 |
|
uploT = (*uplo_P(a) == 'U') ? 1 : -1; |
430 |
|
if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */ |
431 |
|
/* "remove the diagonal entries": */ |
432 |
|
chm_diagN2U(chc, uploT, /* do_realloc */ FALSE); |
433 |
|
diag[0]= 'U'; |
434 |
|
} |
435 |
|
else diag[0]= 'N'; |
436 |
|
} |
437 |
|
SET_VECTOR_ELT(dn, 0, /* establish dimnames */ |
438 |
|
duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0))); |
439 |
|
SET_VECTOR_ELT(dn, 1, |
440 |
|
duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), 1))); |
441 |
|
UNPROTECT(1); |
442 |
|
return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn); |
443 |
|
} |
444 |
|
|
445 |
|
SEXP Csparse_Csparse_crossprod(SEXP a, SEXP b, SEXP trans) |
446 |
|
{ |
447 |
|
int tr = asLogical(trans); |
448 |
|
CHM_SP |
449 |
|
cha = AS_CHM_SP(a), |
450 |
|
chb = AS_CHM_SP(b), |
451 |
|
chTr, chc; |
452 |
|
const char *cl_a = class_P(a), *cl_b = class_P(b); |
453 |
|
char diag[] = {'\0', '\0'}; |
454 |
|
int uploT = 0; |
455 |
|
SEXP dn = PROTECT(allocVector(VECSXP, 2)); |
456 |
|
R_CheckStack(); |
457 |
|
|
458 |
|
chTr = cholmod_transpose((tr) ? chb : cha, chb->xtype, &c); |
459 |
|
chc = cholmod_ssmult((tr) ? cha : chTr, (tr) ? chTr : chb, |
460 |
|
/*out_stype:*/ 0, cha->xtype, /*out sorted:*/ 1, &c); |
461 |
|
cholmod_free_sparse(&chTr, &c); |
462 |
|
|
463 |
|
/* Preserve triangularity and unit-triangularity if appropriate; |
464 |
|
* see Csparse_Csparse_prod() for comments */ |
465 |
|
if (cl_a[1] == 't' && cl_b[1] == 't') |
466 |
|
if(*uplo_P(a) != *uplo_P(b)) { /* one 'U', the other 'L' */ |
467 |
|
uploT = (*uplo_P(b) == 'U') ? 1 : -1; |
468 |
|
if(*diag_P(a) == 'U' && *diag_P(b) == 'U') { /* return UNIT-triag. */ |
469 |
|
chm_diagN2U(chc, uploT, /* do_realloc */ FALSE); |
470 |
|
diag[0]= 'U'; |
471 |
|
} |
472 |
|
else diag[0]= 'N'; |
473 |
|
} |
474 |
|
SET_VECTOR_ELT(dn, 0, /* establish dimnames */ |
475 |
|
duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), (tr) ? 0 : 1))); |
476 |
|
SET_VECTOR_ELT(dn, 1, |
477 |
|
duplicate(VECTOR_ELT(GET_SLOT(b, Matrix_DimNamesSym), (tr) ? 0 : 1))); |
478 |
|
UNPROTECT(1); |
479 |
|
return chm_sparse_to_SEXP(chc, 1, uploT, /*Rkind*/0, diag, dn); |
480 |
} |
} |
481 |
|
|
482 |
SEXP Csparse_dense_prod(SEXP a, SEXP b) |
SEXP Csparse_dense_prod(SEXP a, SEXP b) |
483 |
{ |
{ |
484 |
#ifdef USE_CHOLMOD |
CHM_SP cha = AS_CHM_SP(a); |
485 |
cholmod_sparse *cha = as_cholmod_sparse(a); |
SEXP b_M = PROTECT(mMatrix_as_dgeMatrix2(b, // transpose_if_vector = |
486 |
cholmod_dense *chb = as_cholmod_dense(b); |
cha->ncol == 1)); |
487 |
cholmod_dense *chc = cholmod_allocate_dense(cha->nrow, chb->ncol, |
CHM_DN chb = AS_CHM_DN(b_M); |
488 |
cha->nrow, chb->xtype, &c); |
CHM_DN chc = cholmod_allocate_dense(cha->nrow, chb->ncol, cha->nrow, |
489 |
double alpha = 1, beta = 0; |
chb->xtype, &c); |
490 |
|
SEXP dn = PROTECT(allocVector(VECSXP, 2)); |
491 |
cholmod_sdmult(cha, 0, &alpha, &beta, chb, chc, &c); |
double one[] = {1,0}, zero[] = {0,0}; |
492 |
free(cha); free(chb); |
int nprot = 2; |
493 |
return chm_dense_to_SEXP(chc, 1); |
R_CheckStack(); |
494 |
#else |
/* Tim Davis, please FIXME: currently (2010-11) *fails* when a is a pattern matrix:*/ |
495 |
error("General multiplication requires CHOLMOD"); |
if(cha->xtype == CHOLMOD_PATTERN) { |
496 |
return R_NilValue; /* -Wall */ |
/* warning(_("Csparse_dense_prod(): cholmod_sdmult() not yet implemented for pattern./ ngCMatrix" */ |
497 |
#endif /* USE_CHOLMOD */ |
/* " --> slightly inefficient coercion")); */ |
498 |
|
|
499 |
|
// This *fails* to produce a CHOLMOD_REAL .. |
500 |
|
// CHM_SP chd = cholmod_l_copy(cha, cha->stype, CHOLMOD_REAL, &c); |
501 |
|
// --> use our Matrix-classes |
502 |
|
SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++; |
503 |
|
cha = AS_CHM_SP(da); |
504 |
|
} |
505 |
|
cholmod_sdmult(cha, 0, one, zero, chb, chc, &c); |
506 |
|
SET_VECTOR_ELT(dn, 0, /* establish dimnames */ |
507 |
|
duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 0))); |
508 |
|
SET_VECTOR_ELT(dn, 1, |
509 |
|
duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1))); |
510 |
|
UNPROTECT(nprot); |
511 |
|
return chm_dense_to_SEXP(chc, 1, 0, dn); |
512 |
|
} |
513 |
|
|
514 |
|
SEXP Csparse_dense_crossprod(SEXP a, SEXP b) |
515 |
|
{ |
516 |
|
CHM_SP cha = AS_CHM_SP(a); |
517 |
|
SEXP b_M = PROTECT(mMatrix_as_dgeMatrix2(b, // transpose_if_vector = |
518 |
|
cha->nrow == 1)); |
519 |
|
CHM_DN chb = AS_CHM_DN(b_M); |
520 |
|
CHM_DN chc = cholmod_allocate_dense(cha->ncol, chb->ncol, cha->ncol, |
521 |
|
chb->xtype, &c); |
522 |
|
SEXP dn = PROTECT(allocVector(VECSXP, 2)); int nprot = 2; |
523 |
|
double one[] = {1,0}, zero[] = {0,0}; |
524 |
|
R_CheckStack(); |
525 |
|
// -- see Csparse_dense_prod() above : |
526 |
|
if(cha->xtype == CHOLMOD_PATTERN) { |
527 |
|
SEXP da = PROTECT(nz2Csparse(a, x_double)); nprot++; |
528 |
|
cha = AS_CHM_SP(da); |
529 |
|
} |
530 |
|
cholmod_sdmult(cha, 1, one, zero, chb, chc, &c); |
531 |
|
SET_VECTOR_ELT(dn, 0, /* establish dimnames */ |
532 |
|
duplicate(VECTOR_ELT(GET_SLOT(a, Matrix_DimNamesSym), 1))); |
533 |
|
SET_VECTOR_ELT(dn, 1, |
534 |
|
duplicate(VECTOR_ELT(GET_SLOT(b_M, Matrix_DimNamesSym), 1))); |
535 |
|
UNPROTECT(nprot); |
536 |
|
return chm_dense_to_SEXP(chc, 1, 0, dn); |
537 |
} |
} |
538 |
|
|
539 |
|
/* Computes x'x or x x' -- *also* for Tsparse (triplet = TRUE) |
540 |
|
see Csparse_Csparse_crossprod above for x'y and x y' */ |
541 |
SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet) |
SEXP Csparse_crossprod(SEXP x, SEXP trans, SEXP triplet) |
542 |
{ |
{ |
|
#ifdef USE_CHOLMOD |
|
543 |
int trip = asLogical(triplet), |
int trip = asLogical(triplet), |
544 |
tr = asLogical(trans); /* gets reversed because _aat is tcrossprod */ |
tr = asLogical(trans); /* gets reversed because _aat is tcrossprod */ |
545 |
cholmod_triplet |
#ifdef AS_CHM_DIAGU2N_FIXED_FINALLY |
546 |
*cht = trip ? as_cholmod_triplet(x) : (cholmod_triplet*) NULL; |
CHM_TR cht = trip ? AS_CHM_TR(x) : (CHM_TR) NULL; |
547 |
cholmod_sparse *chcp, *chxt, |
#else /* workaround needed:*/ |
548 |
*chx = trip ? cholmod_triplet_to_sparse(cht, cht->nnz, &c) |
SEXP xx = PROTECT(Tsparse_diagU2N(x)); |
549 |
: as_cholmod_sparse(x); |
CHM_TR cht = trip ? AS_CHM_TR__(xx) : (CHM_TR) NULL; |
550 |
|
#endif |
551 |
|
CHM_SP chcp, chxt, |
552 |
|
chx = (trip ? |
553 |
|
cholmod_triplet_to_sparse(cht, cht->nnz, &c) : |
554 |
|
AS_CHM_SP(x)); |
555 |
|
SEXP dn = PROTECT(allocVector(VECSXP, 2)); |
556 |
|
R_CheckStack(); |
557 |
|
|
558 |
if (!tr) |
if (!tr) chxt = cholmod_transpose(chx, chx->xtype, &c); |
|
chxt = cholmod_transpose(chx, (int) chx->xtype, &c); |
|
559 |
chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c); |
chcp = cholmod_aat((!tr) ? chxt : chx, (int *) NULL, 0, chx->xtype, &c); |
560 |
if(!chcp) |
if(!chcp) { |
561 |
error("Csparse_crossprod(): error return from cholmod_aat()"); |
UNPROTECT(1); |
562 |
|
error(_("Csparse_crossprod(): error return from cholmod_aat()")); |
|
if (trip) { |
|
|
cholmod_free_sparse(&chx, &c); |
|
|
free(cht); |
|
|
} else { |
|
|
free(chx); |
|
563 |
} |
} |
564 |
|
cholmod_band_inplace(0, chcp->ncol, chcp->xtype, chcp, &c); |
565 |
|
chcp->stype = 1; |
566 |
|
if (trip) cholmod_free_sparse(&chx, &c); |
567 |
if (!tr) cholmod_free_sparse(&chxt, &c); |
if (!tr) cholmod_free_sparse(&chxt, &c); |
568 |
return chm_sparse_to_SEXP(chcp, 1); |
SET_VECTOR_ELT(dn, 0, /* establish dimnames */ |
569 |
|
duplicate(VECTOR_ELT(GET_SLOT(x, Matrix_DimNamesSym), |
570 |
|
(tr) ? 0 : 1))); |
571 |
|
SET_VECTOR_ELT(dn, 1, duplicate(VECTOR_ELT(dn, 0))); |
572 |
|
#ifdef AS_CHM_DIAGU2N_FIXED_FINALLY |
573 |
|
UNPROTECT(1); |
574 |
#else |
#else |
575 |
error("General crossproduct requires CHOLMOD"); |
UNPROTECT(2); |
576 |
return R_NilValue; /* -Wall */ |
#endif |
577 |
#endif /* USE_CHOLMOD */ |
return chm_sparse_to_SEXP(chcp, 1, 0, 0, "", dn); |
578 |
|
} |
579 |
|
|
580 |
|
/* Csparse_drop(x, tol): drop entries with absolute value < tol, i.e, |
581 |
|
* at least all "explicit" zeros */ |
582 |
|
SEXP Csparse_drop(SEXP x, SEXP tol) |
583 |
|
{ |
584 |
|
const char *cl = class_P(x); |
585 |
|
/* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ |
586 |
|
int tr = (cl[1] == 't'); |
587 |
|
CHM_SP chx = AS_CHM_SP__(x); |
588 |
|
CHM_SP ans = cholmod_copy(chx, chx->stype, chx->xtype, &c); |
589 |
|
double dtol = asReal(tol); |
590 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
591 |
|
R_CheckStack(); |
592 |
|
|
593 |
|
if(!cholmod_drop(dtol, ans, &c)) |
594 |
|
error(_("cholmod_drop() failed")); |
595 |
|
return chm_sparse_to_SEXP(ans, 1, |
596 |
|
tr ? ((*uplo_P(x) == 'U') ? 1 : -1) : 0, |
597 |
|
Rkind, tr ? diag_P(x) : "", |
598 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
599 |
|
} |
600 |
|
|
601 |
|
SEXP Csparse_horzcat(SEXP x, SEXP y) |
602 |
|
{ |
603 |
|
#define CSPARSE_CAT(_KIND_) \ |
604 |
|
CHM_SP chx = AS_CHM_SP__(x), chy = AS_CHM_SP__(y); \ |
605 |
|
R_CheckStack(); \ |
606 |
|
int Rk_x = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : -3, \ |
607 |
|
Rk_y = (chy->xtype != CHOLMOD_PATTERN) ? Real_kind(y) : -3, Rkind; \ |
608 |
|
if(Rk_x == -3 || Rk_y == -3) { /* at least one of them is patter"n" */ \ |
609 |
|
if(Rk_x == -3 && Rk_y == -3) { /* fine */ \ |
610 |
|
} else { /* only one is a patter"n" \ |
611 |
|
* "Bug" in cholmod_horzcat()/vertcat(): returns patter"n" matrix if one of them is */ \ |
612 |
|
Rboolean ok; \ |
613 |
|
if(Rk_x == -3) { \ |
614 |
|
ok = chm_MOD_xtype(CHOLMOD_REAL, chx, &c); Rk_x = 0; \ |
615 |
|
} else if(Rk_y == -3) { \ |
616 |
|
ok = chm_MOD_xtype(CHOLMOD_REAL, chy, &c); Rk_y = 0; \ |
617 |
|
} else \ |
618 |
|
error(_("Impossible Rk_x/Rk_y in Csparse_%s(), please report"), _KIND_); \ |
619 |
|
if(!ok) \ |
620 |
|
error(_("chm_MOD_xtype() was not successful in Csparse_%s(), please report"), \ |
621 |
|
_KIND_); \ |
622 |
|
} \ |
623 |
|
} \ |
624 |
|
Rkind = /* logical if both x and y are */ (Rk_x == 1 && Rk_y == 1) ? 1 : 0 |
625 |
|
|
626 |
|
CSPARSE_CAT("horzcat"); |
627 |
|
// TODO: currently drops dimnames - and we fix at R level; |
628 |
|
|
629 |
|
return chm_sparse_to_SEXP(cholmod_horzcat(chx, chy, 1, &c), |
630 |
|
1, 0, Rkind, "", R_NilValue); |
631 |
|
} |
632 |
|
|
633 |
|
SEXP Csparse_vertcat(SEXP x, SEXP y) |
634 |
|
{ |
635 |
|
CSPARSE_CAT("vertcat"); |
636 |
|
// TODO: currently drops dimnames - and we fix at R level; |
637 |
|
|
638 |
|
return chm_sparse_to_SEXP(cholmod_vertcat(chx, chy, 1, &c), |
639 |
|
1, 0, Rkind, "", R_NilValue); |
640 |
|
} |
641 |
|
|
642 |
|
SEXP Csparse_band(SEXP x, SEXP k1, SEXP k2) |
643 |
|
{ |
644 |
|
CHM_SP chx = AS_CHM_SP__(x); |
645 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
646 |
|
CHM_SP ans = cholmod_band(chx, asInteger(k1), asInteger(k2), chx->xtype, &c); |
647 |
|
R_CheckStack(); |
648 |
|
|
649 |
|
return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "", |
650 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
651 |
|
} |
652 |
|
|
653 |
|
SEXP Csparse_diagU2N(SEXP x) |
654 |
|
{ |
655 |
|
const char *cl = class_P(x); |
656 |
|
/* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ |
657 |
|
if (cl[1] != 't' || *diag_P(x) != 'U') { |
658 |
|
/* "trivially fast" when not triangular (<==> no 'diag' slot), |
659 |
|
or not *unit* triangular */ |
660 |
|
return (x); |
661 |
|
} |
662 |
|
else { /* unit triangular (diag='U'): "fill the diagonal" & diag:= "N" */ |
663 |
|
CHM_SP chx = AS_CHM_SP__(x); |
664 |
|
CHM_SP eye = cholmod_speye(chx->nrow, chx->ncol, chx->xtype, &c); |
665 |
|
double one[] = {1, 0}; |
666 |
|
CHM_SP ans = cholmod_add(chx, eye, one, one, TRUE, TRUE, &c); |
667 |
|
int uploT = (*uplo_P(x) == 'U') ? 1 : -1; |
668 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
669 |
|
|
670 |
|
R_CheckStack(); |
671 |
|
cholmod_free_sparse(&eye, &c); |
672 |
|
return chm_sparse_to_SEXP(ans, 1, uploT, Rkind, "N", |
673 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
674 |
|
} |
675 |
|
} |
676 |
|
|
677 |
|
SEXP Csparse_diagN2U(SEXP x) |
678 |
|
{ |
679 |
|
const char *cl = class_P(x); |
680 |
|
/* dtCMatrix, etc; [1] = the second character =?= 't' for triangular */ |
681 |
|
if (cl[1] != 't' || *diag_P(x) != 'N') { |
682 |
|
/* "trivially fast" when not triangular (<==> no 'diag' slot), |
683 |
|
or already *unit* triangular */ |
684 |
|
return (x); |
685 |
|
} |
686 |
|
else { /* triangular with diag='N'): now drop the diagonal */ |
687 |
|
/* duplicate, since chx will be modified: */ |
688 |
|
SEXP xx = PROTECT(duplicate(x)); |
689 |
|
CHM_SP chx = AS_CHM_SP__(xx); |
690 |
|
int uploT = (*uplo_P(x) == 'U') ? 1 : -1, |
691 |
|
Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
692 |
|
R_CheckStack(); |
693 |
|
|
694 |
|
chm_diagN2U(chx, uploT, /* do_realloc */ FALSE); |
695 |
|
|
696 |
|
SEXP ans = chm_sparse_to_SEXP(chx, /*dofree*/ 0/* or 1 ?? */, |
697 |
|
uploT, Rkind, "U", |
698 |
|
GET_SLOT(x, Matrix_DimNamesSym)); |
699 |
|
UNPROTECT(1);// only now ! |
700 |
|
return ans; |
701 |
|
} |
702 |
|
} |
703 |
|
|
704 |
|
/** |
705 |
|
* "Indexing" aka subsetting : Compute x[i,j], also for vectors i and j |
706 |
|
* Working via CHOLMOD_submatrix, see ./CHOLMOD/MatrixOps/cholmod_submatrix.c |
707 |
|
* @param x CsparseMatrix |
708 |
|
* @param i row indices (0-origin), or NULL (R's) |
709 |
|
* @param j columns indices (0-origin), or NULL |
710 |
|
* |
711 |
|
* @return x[i,j] still CsparseMatrix --- currently, this loses dimnames |
712 |
|
*/ |
713 |
|
SEXP Csparse_submatrix(SEXP x, SEXP i, SEXP j) |
714 |
|
{ |
715 |
|
CHM_SP chx = AS_CHM_SP(x); /* << does diagU2N() when needed */ |
716 |
|
int rsize = (isNull(i)) ? -1 : LENGTH(i), |
717 |
|
csize = (isNull(j)) ? -1 : LENGTH(j); |
718 |
|
int Rkind = (chx->xtype != CHOLMOD_PATTERN) ? Real_kind(x) : 0; |
719 |
|
R_CheckStack(); |
720 |
|
|
721 |
|
if (rsize >= 0 && !isInteger(i)) |
722 |
|
error(_("Index i must be NULL or integer")); |
723 |
|
if (csize >= 0 && !isInteger(j)) |
724 |
|
error(_("Index j must be NULL or integer")); |
725 |
|
|
726 |
|
#define CHM_SUB(_M_, _i_, _j_) \ |
727 |
|
cholmod_submatrix(_M_, \ |
728 |
|
(rsize < 0) ? NULL : INTEGER(_i_), rsize, \ |
729 |
|
(csize < 0) ? NULL : INTEGER(_j_), csize, \ |
730 |
|
TRUE, TRUE, &c) |
731 |
|
CHM_SP ans; |
732 |
|
if (!chx->stype) {/* non-symmetric Matrix */ |
733 |
|
ans = CHM_SUB(chx, i, j); |
734 |
|
} |
735 |
|
else { |
736 |
|
/* for now, cholmod_submatrix() only accepts "generalMatrix" */ |
737 |
|
CHM_SP tmp = cholmod_copy(chx, /* stype: */ 0, chx->xtype, &c); |
738 |
|
ans = CHM_SUB(tmp, i, j); |
739 |
|
cholmod_free_sparse(&tmp, &c); |
740 |
|
} |
741 |
|
|
742 |
|
// "FIXME": currently dropping dimnames, and adding them afterwards in R : |
743 |
|
/* // dimnames: */ |
744 |
|
/* SEXP x_dns = GET_SLOT(x, Matrix_DimNamesSym), */ |
745 |
|
/* dn = PROTECT(allocVector(VECSXP, 2)); */ |
746 |
|
return chm_sparse_to_SEXP(ans, 1, 0, Rkind, "", /* dimnames: */ R_NilValue); |
747 |
|
} |
748 |
|
|
749 |
|
#define _d_Csp_ |
750 |
|
#include "t_Csparse_subassign.c" |
751 |
|
|
752 |
|
#define _l_Csp_ |
753 |
|
#include "t_Csparse_subassign.c" |
754 |
|
|
755 |
|
#define _i_Csp_ |
756 |
|
#include "t_Csparse_subassign.c" |
757 |
|
|
758 |
|
#define _n_Csp_ |
759 |
|
#include "t_Csparse_subassign.c" |
760 |
|
|
761 |
|
#define _z_Csp_ |
762 |
|
#include "t_Csparse_subassign.c" |
763 |
|
|
764 |
|
|
765 |
|
|
766 |
|
SEXP Csparse_MatrixMarket(SEXP x, SEXP fname) |
767 |
|
{ |
768 |
|
FILE *f = fopen(CHAR(asChar(fname)), "w"); |
769 |
|
|
770 |
|
if (!f) |
771 |
|
error(_("failure to open file \"%s\" for writing"), |
772 |
|
CHAR(asChar(fname))); |
773 |
|
if (!cholmod_write_sparse(f, AS_CHM_SP(x), |
774 |
|
(CHM_SP)NULL, (char*) NULL, &c)) |
775 |
|
error(_("cholmod_write_sparse returned error code")); |
776 |
|
fclose(f); |
777 |
|
return R_NilValue; |
778 |
|
} |
779 |
|
|
780 |
|
|
781 |
|
/** |
782 |
|
* Extract the diagonal entries from *triangular* Csparse matrix __or__ a |
783 |
|
* cholmod_sparse factor (LDL = TRUE). |
784 |
|
* |
785 |
|
* @param n dimension of the matrix. |
786 |
|
* @param x_p 'p' (column pointer) slot contents |
787 |
|
* @param x_x 'x' (non-zero entries) slot contents |
788 |
|
* @param perm 'perm' (= permutation vector) slot contents; only used for "diagBack" |
789 |
|
* @param resultKind a (SEXP) string indicating which kind of result is desired. |
790 |
|
* |
791 |
|
* @return a SEXP, either a (double) number or a length n-vector of diagonal entries |
792 |
|
*/ |
793 |
|
SEXP diag_tC_ptr(int n, int *x_p, double *x_x, Rboolean is_U, int *perm, |
794 |
|
/* ^^^^^^ FIXME[Generalize] to int / ... */ |
795 |
|
SEXP resultKind) |
796 |
|
{ |
797 |
|
const char* res_ch = CHAR(STRING_ELT(resultKind,0)); |
798 |
|
enum diag_kind { diag, diag_backpermuted, trace, prod, sum_log, min, max, range |
799 |
|
} res_kind = ((!strcmp(res_ch, "trace")) ? trace : |
800 |
|
((!strcmp(res_ch, "sumLog")) ? sum_log : |
801 |
|
((!strcmp(res_ch, "prod")) ? prod : |
802 |
|
((!strcmp(res_ch, "min")) ? min : |
803 |
|
((!strcmp(res_ch, "max")) ? max : |
804 |
|
((!strcmp(res_ch, "range")) ? range : |
805 |
|
((!strcmp(res_ch, "diag")) ? diag : |
806 |
|
((!strcmp(res_ch, "diagBack")) ? diag_backpermuted : |
807 |
|
-1)))))))); |
808 |
|
int i, n_x, i_from; |
809 |
|
SEXP ans = PROTECT(allocVector(REALSXP, |
810 |
|
/* ^^^^ FIXME[Generalize] */ |
811 |
|
(res_kind == diag || |
812 |
|
res_kind == diag_backpermuted) ? n : |
813 |
|
(res_kind == range ? 2 : 1))); |
814 |
|
double *v = REAL(ans); |
815 |
|
/* ^^^^^^ ^^^^ FIXME[Generalize] */ |
816 |
|
|
817 |
|
i_from = (is_U ? -1 : 0); |
818 |
|
|
819 |
|
#define for_DIAG(v_ASSIGN) \ |
820 |
|
for(i = 0; i < n; i++) { \ |
821 |
|
/* looking at i-th column */ \ |
822 |
|
n_x = x_p[i+1] - x_p[i];/* #{entries} in this column */ \ |
823 |
|
if( is_U) i_from += n_x; \ |
824 |
|
v_ASSIGN; \ |
825 |
|
if(!is_U) i_from += n_x; \ |
826 |
|
} |
827 |
|
|
828 |
|
/* NOTA BENE: we assume -- uplo = "L" i.e. lower triangular matrix |
829 |
|
* for uplo = "U" (makes sense with a "dtCMatrix" !), |
830 |
|
* should use x_x[i_from + (n_x - 1)] instead of x_x[i_from], |
831 |
|
* where n_x = (x_p[i+1] - x_p[i]) |
832 |
|
*/ |
833 |
|
|
834 |
|
switch(res_kind) { |
835 |
|
case trace: // = sum |
836 |
|
v[0] = 0.; |
837 |
|
for_DIAG(v[0] += x_x[i_from]); |
838 |
|
break; |
839 |
|
|
840 |
|
case sum_log: |
841 |
|
v[0] = 0.; |
842 |
|
for_DIAG(v[0] += log(x_x[i_from])); |
843 |
|
break; |
844 |
|
|
845 |
|
case prod: |
846 |
|
v[0] = 1.; |
847 |
|
for_DIAG(v[0] *= x_x[i_from]); |
848 |
|
break; |
849 |
|
|
850 |
|
case min: |
851 |
|
v[0] = R_PosInf; |
852 |
|
for_DIAG(if(v[0] > x_x[i_from]) v[0] = x_x[i_from]); |
853 |
|
break; |
854 |
|
|
855 |
|
case max: |
856 |
|
v[0] = R_NegInf; |
857 |
|
for_DIAG(if(v[0] < x_x[i_from]) v[0] = x_x[i_from]); |
858 |
|
break; |
859 |
|
|
860 |
|
case range: |
861 |
|
v[0] = R_PosInf; |
862 |
|
v[1] = R_NegInf; |
863 |
|
for_DIAG(if(v[0] > x_x[i_from]) v[0] = x_x[i_from]; |
864 |
|
if(v[1] < x_x[i_from]) v[1] = x_x[i_from]); |
865 |
|
break; |
866 |
|
|
867 |
|
case diag: |
868 |
|
for_DIAG(v[i] = x_x[i_from]); |
869 |
|
break; |
870 |
|
|
871 |
|
case diag_backpermuted: |
872 |
|
for_DIAG(v[i] = x_x[i_from]); |
873 |
|
|
874 |
|
warning(_("%s = '%s' (back-permuted) is experimental"), |
875 |
|
"resultKind", "diagBack"); |
876 |
|
/* now back_permute : */ |
877 |
|
for(i = 0; i < n; i++) { |
878 |
|
double tmp = v[i]; v[i] = v[perm[i]]; v[perm[i]] = tmp; |
879 |
|
/*^^^^ FIXME[Generalize] */ |
880 |
|
} |
881 |
|
break; |
882 |
|
|
883 |
|
default: /* -1 from above */ |
884 |
|
error(_("diag_tC(): invalid 'resultKind'")); |
885 |
|
/* Wall: */ ans = R_NilValue; v = REAL(ans); |
886 |
|
} |
887 |
|
|
888 |
|
UNPROTECT(1); |
889 |
|
return ans; |
890 |
} |
} |
891 |
|
|
892 |
|
/** |
893 |
|
* Extract the diagonal entries from *triangular* Csparse matrix __or__ a |
894 |
|
* cholmod_sparse factor (LDL = TRUE). |
895 |
|
* |
896 |
|
* @param obj -- now a cholmod_sparse factor or a dtCMatrix |
897 |
|
* @param pslot 'p' (column pointer) slot of Csparse matrix/factor |
898 |
|
* @param xslot 'x' (non-zero entries) slot of Csparse matrix/factor |
899 |
|
* @param perm_slot 'perm' (= permutation vector) slot of corresponding CHMfactor; |
900 |
|
* only used for "diagBack" |
901 |
|
* @param resultKind a (SEXP) string indicating which kind of result is desired. |
902 |
|
* |
903 |
|
* @return a SEXP, either a (double) number or a length n-vector of diagonal entries |
904 |
|
*/ |
905 |
|
SEXP diag_tC(SEXP obj, SEXP resultKind) |
906 |
|
{ |
907 |
|
|
908 |
|
SEXP |
909 |
|
pslot = GET_SLOT(obj, Matrix_pSym), |
910 |
|
xslot = GET_SLOT(obj, Matrix_xSym); |
911 |
|
Rboolean is_U = (R_has_slot(obj, Matrix_uploSym) && |
912 |
|
*CHAR(asChar(GET_SLOT(obj, Matrix_uploSym))) == 'U'); |
913 |
|
int n = length(pslot) - 1, /* n = ncol(.) = nrow(.) */ |
914 |
|
*x_p = INTEGER(pslot), pp = -1, *perm; |
915 |
|
double *x_x = REAL(xslot); |
916 |
|
/* ^^^^^^ ^^^^ FIXME[Generalize] to INTEGER(.) / LOGICAL(.) / ... xslot !*/ |
917 |
|
|
918 |
|
if(R_has_slot(obj, Matrix_permSym)) |
919 |
|
perm = INTEGER(GET_SLOT(obj, Matrix_permSym)); |
920 |
|
else perm = &pp; |
921 |
|
|
922 |
|
return diag_tC_ptr(n, x_p, x_x, is_U, perm, resultKind); |
923 |
|
} |
924 |
|
|
925 |
|
|
926 |
|
/** |
927 |
|
* Create a Csparse matrix object from indices and/or pointers. |
928 |
|
* |
929 |
|
* @param cls name of actual class of object to create |
930 |
|
* @param i optional integer vector of length nnz of row indices |
931 |
|
* @param j optional integer vector of length nnz of column indices |
932 |
|
* @param p optional integer vector of length np of row or column pointers |
933 |
|
* @param np length of integer vector p. Must be zero if p == (int*)NULL |
934 |
|
* @param x optional vector of values |
935 |
|
* @param nnz length of vectors i, j and/or x, whichever is to be used |
936 |
|
* @param dims optional integer vector of length 2 to be used as |
937 |
|
* dimensions. If dims == (int*)NULL then the maximum row and column |
938 |
|
* index are used as the dimensions. |
939 |
|
* @param dimnames optional list of length 2 to be used as dimnames |
940 |
|
* @param index1 indicator of 1-based indices |
941 |
|
* |
942 |
|
* @return an SEXP of class cls inheriting from CsparseMatrix. |
943 |
|
*/ |
944 |
|
SEXP create_Csparse(char* cls, int* i, int* j, int* p, int np, |
945 |
|
void* x, int nnz, int* dims, SEXP dimnames, |
946 |
|
int index1) |
947 |
|
{ |
948 |
|
SEXP ans; |
949 |
|
int *ij = (int*)NULL, *tri, *trj, |
950 |
|
mi, mj, mp, nrow = -1, ncol = -1; |
951 |
|
int xtype = -1; /* -Wall */ |
952 |
|
CHM_TR T; |
953 |
|
CHM_SP A; |
954 |
|
|
955 |
|
if (np < 0 || nnz < 0) |
956 |
|
error(_("negative vector lengths not allowed: np = %d, nnz = %d"), |
957 |
|
np, nnz); |
958 |
|
if (1 != ((mi = (i == (int*)NULL)) + |
959 |
|
(mj = (j == (int*)NULL)) + |
960 |
|
(mp = (p == (int*)NULL)))) |
961 |
|
error(_("exactly 1 of 'i', 'j' or 'p' must be NULL")); |
962 |
|
if (mp) { |
963 |
|
if (np) error(_("np = %d, must be zero when p is NULL"), np); |
964 |
|
} else { |
965 |
|
if (np) { /* Expand p to form i or j */ |
966 |
|
if (!(p[0])) error(_("p[0] = %d, should be zero"), p[0]); |
967 |
|
for (int ii = 0; ii < np; ii++) |
968 |
|
if (p[ii] > p[ii + 1]) |
969 |
|
error(_("p must be non-decreasing")); |
970 |
|
if (p[np] != nnz) |
971 |
|
error("p[np] = %d != nnz = %d", p[np], nnz); |
972 |
|
ij = Calloc(nnz, int); |
973 |
|
if (mi) { |
974 |
|
i = ij; |
975 |
|
nrow = np; |
976 |
|
} else { |
977 |
|
j = ij; |
978 |
|
ncol = np; |
979 |
|
} |
980 |
|
/* Expand p to 0-based indices */ |
981 |
|
for (int ii = 0; ii < np; ii++) |
982 |
|
for (int jj = p[ii]; jj < p[ii + 1]; jj++) ij[jj] = ii; |
983 |
|
} else { |
984 |
|
if (nnz) |
985 |
|
error(_("Inconsistent dimensions: np = 0 and nnz = %d"), |
986 |
|
nnz); |
987 |
|
} |
988 |
|
} |
989 |
|
/* calculate nrow and ncol */ |
990 |
|
if (nrow < 0) { |
991 |
|
for (int ii = 0; ii < nnz; ii++) { |
992 |
|
int i1 = i[ii] + (index1 ? 0 : 1); /* 1-based index */ |
993 |
|
if (i1 < 1) error(_("invalid row index at position %d"), ii); |
994 |
|
if (i1 > nrow) nrow = i1; |
995 |
|
} |
996 |
|
} |
997 |
|
if (ncol < 0) { |
998 |
|
for (int jj = 0; jj < nnz; jj++) { |
999 |
|
int j1 = j[jj] + (index1 ? 0 : 1); |
1000 |
|
if (j1 < 1) error(_("invalid column index at position %d"), jj); |
1001 |
|
if (j1 > ncol) ncol = j1; |
1002 |
|
} |
1003 |
|
} |
1004 |
|
if (dims != (int*)NULL) { |
1005 |
|
if (dims[0] > nrow) nrow = dims[0]; |
1006 |
|
if (dims[1] > ncol) ncol = dims[1]; |
1007 |
|
} |
1008 |
|
/* check the class name */ |
1009 |
|
if (strlen(cls) != 8) |
1010 |
|
error(_("strlen of cls argument = %d, should be 8"), strlen(cls)); |
1011 |
|
if (!strcmp(cls + 2, "CMatrix")) |
1012 |
|
error(_("cls = \"%s\" does not end in \"CMatrix\""), cls); |
1013 |
|
switch(cls[0]) { |
1014 |
|
case 'd': |
1015 |
|
case 'l': |
1016 |
|
xtype = CHOLMOD_REAL; |
1017 |
|
break; |
1018 |
|
case 'n': |
1019 |
|
xtype = CHOLMOD_PATTERN; |
1020 |
|
break; |
1021 |
|
default: |
1022 |
|
error(_("cls = \"%s\" must begin with 'd', 'l' or 'n'"), cls); |
1023 |
|
} |
1024 |
|
if (cls[1] != 'g') |
1025 |
|
error(_("Only 'g'eneral sparse matrix types allowed")); |
1026 |
|
/* allocate and populate the triplet */ |
1027 |
|
T = cholmod_allocate_triplet((size_t)nrow, (size_t)ncol, (size_t)nnz, 0, |
1028 |
|
xtype, &c); |
1029 |
|
T->x = x; |
1030 |
|
tri = (int*)T->i; |
1031 |
|
trj = (int*)T->j; |
1032 |
|
for (int ii = 0; ii < nnz; ii++) { |
1033 |
|
tri[ii] = i[ii] - ((!mi && index1) ? 1 : 0); |
1034 |
|
trj[ii] = j[ii] - ((!mj && index1) ? 1 : 0); |
1035 |
|
} |
1036 |
|
/* create the cholmod_sparse structure */ |
1037 |
|
A = cholmod_triplet_to_sparse(T, nnz, &c); |
1038 |
|
cholmod_free_triplet(&T, &c); |
1039 |
|
/* copy the information to the SEXP */ |
1040 |
|
ans = PROTECT(NEW_OBJECT(MAKE_CLASS(cls))); |
1041 |
|
// FIXME: This has been copied from chm_sparse_to_SEXP in chm_common.c |
1042 |
|
/* allocate and copy common slots */ |
1043 |
|
nnz = cholmod_nnz(A, &c); |
1044 |
|
dims = INTEGER(ALLOC_SLOT(ans, Matrix_DimSym, INTSXP, 2)); |
1045 |
|
dims[0] = A->nrow; dims[1] = A->ncol; |
1046 |
|
Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_pSym, INTSXP, A->ncol + 1)), (int*)A->p, A->ncol + 1); |
1047 |
|
Memcpy(INTEGER(ALLOC_SLOT(ans, Matrix_iSym, INTSXP, nnz)), (int*)A->i, nnz); |
1048 |
|
switch(cls[1]) { |
1049 |
|
case 'd': |
1050 |
|
Memcpy(REAL(ALLOC_SLOT(ans, Matrix_xSym, REALSXP, nnz)), (double*)A->x, nnz); |
1051 |
|
break; |
1052 |
|
case 'l': |
1053 |
|
error(_("code not yet written for cls = \"lgCMatrix\"")); |
1054 |
|
} |
1055 |
|
/* FIXME: dimnames are *NOT* put there yet (if non-NULL) */ |
1056 |
|
cholmod_free_sparse(&A, &c); |
1057 |
|
UNPROTECT(1); |
1058 |
|
return ans; |
1059 |
|
} |