Matrix r5059
Loading...
Searching...
No Matches
validity.c
Go to the documentation of this file.
1#include "Mdefines.h"
2
3#define MK(_FORMAT_ ) Rf_mkString(_FORMAT_ )
4#define MS(_FORMAT_, ...) Matrix_sprintf(_FORMAT_, __VA_ARGS__)
5
6#define RMK(_FORMAT_ ) \
7 return MK( _FORMAT_ )
8#define RMS(_FORMAT_, ...) \
9 return MS(_FORMAT_, __VA_ARGS__)
10#define RMKMS(_FORMAT_, ...) \
11 return MK(MS(_FORMAT_, __VA_ARGS__))
12
13#define FRMK(_FORMAT_ ) \
14 do { \
15 Matrix_Free(work, lwork); \
16 RMK (_FORMAT_ ); \
17 } while (0)
18#define FRMS(_FORMAT_, ...) \
19 do { \
20 Matrix_Free(work, lwork); \
21 RMS (_FORMAT_, __VA_ARGS__); \
22 } while (0)
23#define FRMKMS(_FORMAT_, ...) \
24 do { \
25 Matrix_Free(work, lwork); \
26 RMKMS(_FORMAT_, __VA_ARGS__); \
27 } while (0)
28
29
30/* Slot validity methods ===============================================
31 Called by various class validity methods (see below).
32*/
33
34static
35char *Dim_validate(SEXP dim)
36{
37 if (TYPEOF(dim) != INTSXP)
38 RMS(_("'%s' slot is not of type \"%s\""), "Dim", "integer");
39 if (XLENGTH(dim) != 2)
40 RMS(_("'%s' slot does not have length %d"), "Dim", 2);
41 int *pdim = INTEGER(dim), m = pdim[0], n = pdim[1];
42 if (m == NA_INTEGER || n == NA_INTEGER)
43 RMS(_("'%s' slot contains NA"), "Dim");
44 if (m < 0 || n < 0)
45 RMS(_("'%s' slot has negative elements"), "Dim");
46
47 return NULL;
48}
49
50SEXP R_Dim_validate(SEXP dim)
51{
52 char *msg = Dim_validate(dim);
53 return (msg) ? Rf_mkString(msg) : Rf_ScalarLogical(1);
54}
55
56static
57char *DimNames_validate(SEXP dimnames, int *pdim)
58{
59 if (TYPEOF(dimnames) != VECSXP)
60 RMS(_("'%s' slot is not a list"), "Dimnames");
61 if (XLENGTH(dimnames) != 2)
62 RMS(_("'%s' slot does not have length %d"), "Dimnames", 2);
63
64 /* Behave as do_matrix() from src/main/array.c:
65 Dimnames[[i]] must be NULL or _coercible to_ character
66 of length Dim[i] or 0 ... see R_Dimnames_fixup() below
67 */
68
69 SEXP s;
70 int i;
71 R_xlen_t ns;
72
73 for (i = 0; i < 2; ++i) {
74 s = VECTOR_ELT(dimnames, i);
75 if (s == R_NilValue)
76 continue;
77 if (!Rf_isVector(s))
78 RMS(_("%s[[%d]] is not NULL or a vector"), "Dimnames", i + 1);
79 ns = XLENGTH(s);
80 if (ns != pdim[i] && ns != 0)
81 RMS(_("length of %s[[%d]] (%lld) is not equal to %s[%d] (%d)"),
82 "Dimnames", i + 1, (long long) ns,
83 "Dim" , i + 1, pdim[i]);
84 }
85
86 return NULL;
87}
88
89SEXP R_DimNames_validate(SEXP dimnames, SEXP dim)
90{
91 char *msg = DimNames_validate(dimnames, INTEGER(dim));
92 return (msg) ? Rf_mkString(msg) : Rf_ScalarLogical(1);
93}
94
95SEXP R_DimNames_fixup(SEXP dimnames)
96{
97 SEXP s;
98 int i, fixup = 0;
99 for (i = 0; i < 2 && !fixup; ++i)
100 fixup =
101 (s = VECTOR_ELT(dimnames, i)) != R_NilValue &&
102 (LENGTH(s) == 0 || TYPEOF(s) != STRSXP);
103 if (!fixup)
104 return dimnames;
105 SEXP dimnames_ = PROTECT(Rf_allocVector(VECSXP, 2));
106 for (i = 0; i < 2; ++i) {
107 if ((s = VECTOR_ELT(dimnames, i)) == R_NilValue || LENGTH(s) == 0)
108 continue;
109 if (TYPEOF(s) == STRSXP)
110 PROTECT(s);
111 else if (TYPEOF(s) == INTSXP && Rf_inherits(s, "factor"))
112 PROTECT(s = Rf_asCharacterFactor(s));
113 else {
114 PROTECT(s = Rf_coerceVector(s, STRSXP));
115 CLEAR_ATTRIB(s);
116 }
117 SET_VECTOR_ELT(dimnames_, i, s);
118 UNPROTECT(1); /* s */
119 }
120 s = Rf_getAttrib(dimnames, R_NamesSymbol);
121 if (s != R_NilValue) {
122 PROTECT(s);
123 Rf_setAttrib(dimnames_, R_NamesSymbol, s);
124 UNPROTECT(1); /* s */
125 }
126 UNPROTECT(1); /* dimnames_ */
127 return dimnames_;
128}
129
130
131/* Class validity methods ==============================================
132 NB: These assume that validity methods for superclasses
133 have already been called via validObject() ...
134*/
135
136SEXP Matrix_validate(SEXP obj)
137{
138 SEXP dim = PROTECT(GET_SLOT(obj, Matrix_DimSym));
139 char *msg = Dim_validate(dim);
140 if (!msg) {
141 SEXP dimnames = PROTECT(GET_SLOT(obj, Matrix_DimNamesSym));
142 msg = DimNames_validate(dimnames, INTEGER(dim));
143 UNPROTECT(1); /* dimnames */
144 }
145 UNPROTECT(1); /* dim */
146 return (msg) ? Rf_mkString(msg) : Rf_ScalarLogical(1);
147}
148
149#define KINDMATRIX_VALIDATE(_PREFIX_, _SEXPTYPE_) \
150SEXP _PREFIX_ ## Matrix_validate(SEXP obj) \
151{ \
152 SEXP x = GET_SLOT(obj, Matrix_xSym); \
153 if (TYPEOF(x) != _SEXPTYPE_) \
154 RMKMS(_("'%s' slot is not of type \"%s\""), \
155 "x", Rf_type2char(_SEXPTYPE_)); \
156 return Rf_ScalarLogical(1); \
157}
158KINDMATRIX_VALIDATE(n, LGLSXP)
159KINDMATRIX_VALIDATE(l, LGLSXP)
160KINDMATRIX_VALIDATE(i, INTSXP)
161KINDMATRIX_VALIDATE(d, REALSXP)
162KINDMATRIX_VALIDATE(z, CPLXSXP)
163#undef KINDMATRIX_VALIDATE
164
166{
167 SEXP factors = GET_SLOT(obj, Matrix_factorsSym);
168 if (TYPEOF(factors) != VECSXP)
169 RMKMS(_("'%s' slot is not of type \"%s\""), "factors", "list");
170 if (XLENGTH(factors) > 0) {
171 PROTECT(factors);
172 SEXP nms = Rf_getAttrib(factors, R_NamesSymbol);
173 UNPROTECT(1); /* factors */
174 if (nms == R_NilValue)
175 RMKMS(_("'%s' slot has no '%s' attribute"), "factors", "names");
176 }
177
178 return Rf_ScalarLogical(1);
179}
180
182{
183 int *pdim = DIM(obj), n = pdim[1];
184 if (pdim[0] != n)
185 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
186
187#ifdef ENFORCE_SYMMETRIC_DIMNAMES
188
189 /* This check can be expensive when both rownames and colnames have
190 nonzero length, and even more so when coercions to character are
191 required ... Users can avoid the expense by setting at least one
192 of rownames and colnames to NULL or by ensuring that they are the
193 same object, as testing for pointer equality is fast ...
194 */
195
196# define ANY_TO_STRING(x) \
197 ((TYPEOF(x) == INTSXP && Rf_inherits(x, "factor")) \
198 ? Rf_asCharacterFactor(x) \
199 : Rf_coerceVector(x, STRSXP))
200
201 SEXP dn = PROTECT(GET_SLOT(obj, Matrix_DimNamesSym)),
202 ndn = Rf_getAttrib(dn, R_NamesSymbol);
203 UNPROTECT(1); /* dn */
204
205 const char *ndn0, *ndn1;
206 if (ndn != R_NilValue &&
207 *(ndn0 = CHAR(STRING_ELT(ndn, 0))) != '\0' &&
208 *(ndn1 = CHAR(STRING_ELT(ndn, 1))) != '\0' &&
209 strcmp(ndn0, ndn1) != 0)
210 RMKMS(_("%s[1] differs from %s[2]"), "Dimnames", "Dimnames");
211 if (n > 0) {
212 /* NB: It is already known that the length of 'dn[[i]]' is 0 or 'n' */
213 SEXP rn, cn;
214 if ((rn = VECTOR_ELT(dn, 0)) != R_NilValue &&
215 (cn = VECTOR_ELT(dn, 1)) != R_NilValue &&
216 LENGTH(rn) == n && LENGTH(cn) == n && rn != cn) {
217 PROTECT(rn);
218 PROTECT(cn);
219 PROTECT(rn = ANY_TO_STRING(rn));
220 PROTECT(cn = ANY_TO_STRING(cn));
221 UNPROTECT(4); /* cn, rn */
222 if (!equalString(rn, cn, n))
223 RMKMS(_("%s[1] differs from %s[2]"), "Dimnames", "Dimnames");
224 }
225 }
226
227# undef ANY_TO_STRING
228
229#endif
230
231 SEXP uplo = GET_SLOT(obj, Matrix_uploSym);
232 if (TYPEOF(uplo) != STRSXP)
233 RMKMS(_("'%s' slot is not of type \"%s\""), "uplo", "character");
234 if (XLENGTH(uplo) != 1)
235 RMKMS(_("'%s' slot does not have length %d"), "uplo", 1);
236 const char *ul = CHAR(STRING_ELT(uplo, 0));
237 if (ul[0] == '\0' || ul[1] != '\0' || (ul[0] != 'U' && ul[0] != 'L'))
238 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "uplo", "U", "L");
239
240 if (HAS_SLOT(obj, Matrix_transSym)) {
241 SEXP trans = GET_SLOT(obj, Matrix_transSym);
242 if (TYPEOF(trans) != STRSXP)
243 RMKMS(_("'%s' slot is not of type \"%s\""), "trans", "character");
244 if (XLENGTH(trans) != 1)
245 RMKMS(_("'%s' slot does not have length %d"), "trans", 1);
246 const char *ct = CHAR(STRING_ELT(trans, 0));
247 if (ct[0] == '\0' || ct[1] != '\0' || (ct[0] != 'C' && ct[0] != 'T'))
248 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "trans", "C", "T");
249 }
250
251 return generalMatrix_validate(obj);
252}
253
255{
256 int *pdim = DIM(obj), n = pdim[1];
257 if (pdim[0] != n)
258 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
259
260 SEXP uplo = GET_SLOT(obj, Matrix_uploSym);
261 if (TYPEOF(uplo) != STRSXP)
262 RMKMS(_("'%s' slot is not of type \"%s\""), "uplo", "character");
263 if (XLENGTH(uplo) != 1)
264 RMKMS(_("'%s' slot does not have length %d"), "uplo", 1);
265 const char *ul = CHAR(STRING_ELT(uplo, 0));
266 if (ul[0] == '\0' || ul[1] != '\0' || (ul[0] != 'U' && ul[0] != 'L'))
267 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "uplo", "U", "L");
268
269 SEXP diag = GET_SLOT(obj, Matrix_diagSym);
270 if (TYPEOF(diag) != STRSXP)
271 RMKMS(_("'%s' slot is not of type \"%s\""), "diag", "character");
272 if (XLENGTH(diag) != 1)
273 RMKMS(_("'%s' slot does not have length %d"), "diag", 1);
274 const char *nu = CHAR(STRING_ELT(diag, 0));
275 if (nu[0] == '\0' || nu[1] != '\0' || (nu[0] != 'N' && nu[0] != 'U'))
276 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "diag", "N", "U");
277
278 return Rf_ScalarLogical(1);
279}
280
282{
283 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
284 SEXP x = GET_SLOT(obj, Matrix_xSym);
285 if (XLENGTH(x) != (int_fast64_t) m * n)
286 RMKMS(_("'%s' slot does not have length %s"), "x", "prod(Dim)");
287 return Rf_ScalarLogical(1);
288}
289
291{
292 int n = DIM(obj)[1];
293 SEXP x = GET_SLOT(obj, Matrix_xSym);
294 if (XLENGTH(x) != n + ((int_fast64_t) n * (n - 1)) / 2)
295 RMKMS(_("'%s' slot does not have length %s"), "x", "Dim[1]*(Dim[1]+1)/2");
296 return Rf_ScalarLogical(1);
297}
298
300{
301 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
302
303 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
304 i = PROTECT(GET_SLOT(obj, Matrix_iSym));
305 UNPROTECT(2); /* i, p */
306
307 if (TYPEOF(p) != INTSXP)
308 RMKMS(_("'%s' slot is not of type \"%s\""), "p", "integer");
309 if (XLENGTH(p) - 1 != n)
310 RMKMS(_("'%s' slot does not have length %s"), "p", "Dim[2]+1");
311 int *pp = INTEGER(p);
312 if (pp[0] != 0)
313 RMKMS(_("first element of '%s' slot is not 0"), "p");
314 int j;
315 for (j = 1; j <= n; ++j) {
316 if (pp[j] == NA_INTEGER)
317 RMKMS(_("'%s' slot contains NA"), "p");
318 if (pp[j] < pp[j - 1])
319 RMKMS(_("'%s' slot is not nondecreasing"), "p");
320 if (pp[j] - pp[j - 1] > m)
321 RMKMS(_("first differences of '%s' slot exceed %s"), "p", "Dim[1]");
322 }
323
324 if (TYPEOF(i) != INTSXP)
325 RMKMS(_("'%s' slot is not of type \"%s\""), "i", "integer");
326 if (XLENGTH(i) < pp[n])
327 RMKMS(_("'%s' slot has length less than %s"), "i", "p[length(p)]");
328 int *pi = INTEGER(i), k, kend, ik, i0;
329 for (j = 1, k = 0; j <= n; ++j) {
330 kend = pp[j];
331 i0 = -1;
332 while (k < kend) {
333 ik = pi[k];
334 if (ik == NA_INTEGER)
335 RMKMS(_("'%s' slot contains NA"), "i");
336 if (ik < 0 || ik >= m)
337 RMKMS(_("'%s' slot has elements not in {%s}"),
338 "i", "0,...,Dim[1]-1");
339 if (ik <= i0)
340 RMKMS(_("'%s' slot is not increasing within columns"), "i");
341 i0 = ik;
342 ++k;
343 }
344 }
345
346 return Rf_ScalarLogical(1);
347}
348
350{
351 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
352
353 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
354 j = PROTECT(GET_SLOT(obj, Matrix_jSym));
355 UNPROTECT(2); /* j, p */
356
357 if (TYPEOF(p) != INTSXP)
358 RMKMS(_("'%s' slot is not of type \"%s\""), "p", "integer");
359 if (XLENGTH(p) - 1 != m)
360 RMKMS(_("'%s' slot does not have length %s"), "p", "Dim[1]+1");
361 int *pp = INTEGER(p);
362 if (pp[0] != 0)
363 RMKMS(_("first element of '%s' slot is not 0"), "p");
364 int i;
365 for (i = 1; i <= m; ++i) {
366 if (pp[i] == NA_INTEGER)
367 RMKMS(_("'%s' slot contains NA"), "p");
368 if (pp[i] < pp[i - 1])
369 RMKMS(_("'%s' slot is not nondecreasing"), "p");
370 if (pp[i] - pp[i - 1] > n)
371 RMKMS(_("first differences of '%s' slot exceed %s"), "p", "Dim[2]");
372 }
373
374 if (TYPEOF(j) != INTSXP)
375 RMKMS(_("'%s' slot is not of type \"%s\""), "j", "integer");
376 if (XLENGTH(j) < pp[m])
377 RMKMS(_("'%s' slot has length less than %s"), "j", "p[length(p)]");
378 int *pj = INTEGER(j), k, kend, jk, j0;
379 for (i = 1, k = 0; i <= m; ++i) {
380 kend = pp[i];
381 j0 = -1;
382 while (k < kend) {
383 jk = pj[k];
384 if (jk == NA_INTEGER)
385 RMKMS(_("'%s' slot contains NA"), "j");
386 if (jk < 0 || jk >= n)
387 RMKMS(_("'%s' slot has elements not in {%s}"),
388 "j", "0,...,Dim[2]-1");
389 if (jk <= j0)
390 RMKMS(_("'%s' slot is not increasing within rows"), "j");
391 j0 = jk;
392 ++k;
393 }
394 }
395
396 return Rf_ScalarLogical(1);
397}
398
400{
401 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
402
403 SEXP i = PROTECT(GET_SLOT(obj, Matrix_iSym)),
404 j = PROTECT(GET_SLOT(obj, Matrix_jSym));
405 UNPROTECT(2); /* j, i */
406
407 if (TYPEOF(i) != INTSXP)
408 RMKMS(_("'%s' slot is not of type \"%s\""), "i", "integer");
409 if (TYPEOF(j) != INTSXP)
410 RMKMS(_("'%s' slot is not of type \"%s\""), "j", "integer");
411 R_xlen_t nnz = XLENGTH(i);
412 if (XLENGTH(j) != nnz)
413 RMKMS(_("'%s' and '%s' slots do not have equal length"), "i", "j");
414 if (nnz > 0) {
415 if (m == 0 || n == 0)
416 RMKMS(_("'%s' slot has nonzero length but %s is 0"), "i", "prod(Dim)");
417 int *pi = INTEGER(i), *pj = INTEGER(j);
418 while (nnz--) {
419 if (*pi == NA_LOGICAL)
420 RMKMS(_("'%s' slot contains NA"), "i");
421 if (*pj == NA_LOGICAL)
422 RMKMS(_("'%s' slot contains NA"), "j");
423 if (*pi < 0 || *pi >= m)
424 RMKMS(_("'%s' slot has elements not in {%s}"),
425 "i", "0,...,Dim[1]-1");
426 if (*pj < 0 || *pj >= n)
427 RMKMS(_("'%s' slot has elements not in {%s}"),
428 "j", "0,...,Dim[2]-1");
429 ++pi;
430 ++pj;
431 }
432 }
433
434 return Rf_ScalarLogical(1);
435}
436
438{
439 int *pdim = DIM(obj), n = pdim[1];
440 if (pdim[0] != n)
441 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
442
443 SEXP diag = GET_SLOT(obj, Matrix_diagSym);
444 if (TYPEOF(diag) != STRSXP)
445 RMKMS(_("'%s' slot is not of type \"%s\""), "diag", "character");
446 if (XLENGTH(diag) != 1)
447 RMKMS(_("'%s' slot does not have length %d"), "diag", 1);
448 const char *nu = CHAR(STRING_ELT(diag, 0));
449 if (nu[0] == '\0' || nu[1] != '\0' || (nu[0] != 'N' && nu[0] != 'U'))
450 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "diag", "N", "U");
451 int nonunit = nu[0] == 'N';
452
453 SEXP x = GET_SLOT(obj, Matrix_xSym);
454 if (XLENGTH(x) != ((nonunit) ? n : 0))
455 RMKMS(_("'%s' slot is \"%s\" but '%s' slot does not have length %s"),
456 "diag", (nonunit) ? "N" : "U", "x", (nonunit) ? "Dim[1]" : "0");
457
458 return Rf_ScalarLogical(1);
459}
460
461SEXP indMatrix_validate(SEXP obj)
462{
463 SEXP margin = GET_SLOT(obj, Matrix_marginSym);
464 if (TYPEOF(margin) != INTSXP)
465 RMKMS(_("'%s' slot is not of type \"%s\""), "margin", "integer");
466 if (XLENGTH(margin) != 1)
467 RMKMS(_("'%s' slot does not have length %d"), "margin", 1);
468 int mg = INTEGER(margin)[0] - 1;
469 if (mg != 0 && mg != 1)
470 RMKMS(_("'%s' slot is not %d or %d"), "margin", 1, 2);
471
472 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
473 if (m > 0 && n == 0)
474 RMKMS(_("%s-by-%s %s invalid for positive '%s' when %s=%d"),
475 (mg == 0) ? "m" : "0", (mg == 0) ? "0" : "n", "indMatrix",
476 (mg == 0) ? "m" : "n", "margin", (mg == 0) ? 1 : 2);
477
478 SEXP perm = GET_SLOT(obj, Matrix_permSym);
479 if (TYPEOF(perm) != INTSXP)
480 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
481 if (XLENGTH(perm) != m)
482 RMKMS(_("'%s' slot does not have length %s"), "perm", "Dim[margin]");
483 int *pperm = INTEGER(perm);
484 while (m--) {
485 if (*pperm == NA_INTEGER)
486 RMKMS(_("'%s' slot contains NA"), "perm");
487 if (*pperm < 1 || *pperm > n)
488 RMKMS(_("'%s' slot has elements not in {%s}"),
489 "perm", "1,...,Dim[1+margin%%2]");
490 ++pperm;
491 }
492
493 return Rf_ScalarLogical(1);
494}
495
496SEXP pMatrix_validate(SEXP obj)
497{
498 int *pdim = DIM(obj), n = pdim[1];
499 if (pdim[0] != n)
500 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
501
502 if (n > 1) {
503 SEXP perm = GET_SLOT(obj, Matrix_permSym);
504 char *work;
505 int lwork = n;
506 Matrix_Calloc(work, lwork, char);
507 int j, *pperm = INTEGER(perm);
508 for (j = 0; j < n; ++j) {
509 if (work[*pperm - 1])
510 FRMKMS(_("'%s' slot contains duplicates"), "perm");
511 work[*(pperm++) - 1] = 1;
512 }
513 Matrix_Free(work, lwork);
514 }
515
516 return Rf_ScalarLogical(1);
517}
518
519SEXP sCMatrix_validate(SEXP obj)
520{
521 SEXP p = GET_SLOT(obj, Matrix_pSym);
522 int *pp = INTEGER(p), n = (int) (XLENGTH(p) - 1);
523 if (pp[n] > 0) {
524 PROTECT(p);
525
526 char ul = UPLO(obj);
527
528 SEXP i = GET_SLOT(obj, Matrix_iSym);
529 int *pi = INTEGER(i), j, k, kend;
530
531 UNPROTECT(1); /* p */
532 if (ul == 'U') {
533 for (j = 0, k = 0; j < n; ++j) {
534 kend = pp[j + 1];
535 while (k < kend) {
536 if (pi[k] > j)
537 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
538 "uplo", "U");
539 ++k;
540 }
541 }
542 } else {
543 for (j = 0, k = 0; j < n; ++j) {
544 kend = pp[j + 1];
545 while (k < kend) {
546 if (pi[k] < j)
547 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
548 "uplo", "L");
549 ++k;
550 }
551 }
552 }
553 }
554
555 return Rf_ScalarLogical(1);
556}
557
558SEXP tCMatrix_validate(SEXP obj)
559{
560 if (DIAG(obj) == 'N')
561 return sCMatrix_validate(obj);
562
563 SEXP p = GET_SLOT(obj, Matrix_pSym);
564 int *pp = INTEGER(p), n = (int) (XLENGTH(p) - 1);
565 if (pp[n] > 0) {
566 PROTECT(p);
567
568 char ul = UPLO(obj);
569
570 SEXP i = GET_SLOT(obj, Matrix_iSym);
571 int *pi = INTEGER(i), j, k, kend;
572
573 UNPROTECT(1); /* p */
574 if (ul == 'U') {
575 for (j = 0, k = 0; j < n; ++j) {
576 kend = pp[j + 1];
577 while (k < kend) {
578 if (pi[k] > j)
579 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
580 "uplo", "U");
581 if (pi[k] == j)
582 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
583 "diag", "U");
584 ++k;
585 }
586 }
587 } else {
588 for (j = 0, k = 0; j < n; ++j) {
589 kend = pp[j + 1];
590 while (k < kend) {
591 if (pi[k] < j)
592 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
593 "uplo", "L");
594 if (pi[k] == j)
595 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
596 "diag", "U");
597 ++k;
598 }
599 }
600 }
601 }
602
603 return Rf_ScalarLogical(1);
604}
605
606SEXP sRMatrix_validate(SEXP obj)
607{
608 SEXP p = GET_SLOT(obj, Matrix_pSym);
609 int *pp = INTEGER(p), m = (int) (XLENGTH(p) - 1);
610 if (pp[m] > 0) {
611 PROTECT(p);
612
613 char ul = UPLO(obj);
614
615 SEXP j = GET_SLOT(obj, Matrix_jSym);
616 int *pj = INTEGER(j), i, k, kend;
617
618 UNPROTECT(1); /* p */
619 if (ul == 'U') {
620 for (i = 0, k = 0; i < m; ++i) {
621 kend = pp[i + 1];
622 while (k < kend) {
623 if (pj[k] < i)
624 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
625 "uplo", "U");
626 ++k;
627 }
628 }
629 } else {
630 for (i = 0, k = 0; i < m; ++i) {
631 kend = pp[i + 1];
632 while (k < kend) {
633 if (pj[k] > i)
634 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
635 "uplo", "L");
636 ++k;
637 }
638 }
639 }
640 }
641
642 return Rf_ScalarLogical(1);
643}
644
645SEXP tRMatrix_validate(SEXP obj)
646{
647 if (DIAG(obj) == 'N')
648 return sRMatrix_validate(obj);
649
650 SEXP p = GET_SLOT(obj, Matrix_pSym);
651 int *pp = INTEGER(p), m = (int) (XLENGTH(p) - 1);
652 if (pp[m] > 0) {
653 PROTECT(p);
654
655 char ul = UPLO(obj);
656
657 SEXP j = GET_SLOT(obj, Matrix_jSym);
658 int *pj = INTEGER(j), i, k, kend;
659
660 UNPROTECT(1); /* p */
661 if (ul == 'U') {
662 for (i = 0, k = 0; i < m; ++i) {
663 kend = pp[i + 1];
664 while (k < kend) {
665 if (pj[k] < i)
666 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
667 "uplo", "U");
668 if (pj[k] == i)
669 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
670 "diag", "U");
671 ++k;
672 }
673 }
674 } else {
675 for (i = 0, k = 0; i < m; ++i) {
676 kend = pp[i + 1];
677 while (k < kend) {
678 if (pj[k] > i)
679 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
680 "uplo", "L");
681 if (pj[k] == i)
682 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
683 "diag", "U");
684 ++k;
685 }
686 }
687 }
688 }
689
690 return Rf_ScalarLogical(1);
691}
692
693SEXP sTMatrix_validate(SEXP obj)
694{
695 SEXP i = GET_SLOT(obj, Matrix_iSym);
696 R_xlen_t nnz = XLENGTH(i);
697 if (nnz > 0) {
698 PROTECT(i);
699
700 char ul = UPLO(obj);
701
702 SEXP j = GET_SLOT(obj, Matrix_jSym);
703 int *pi = INTEGER(i), *pj = INTEGER(j);
704
705 UNPROTECT(1); /* i */
706 if (ul == 'U') {
707 while (nnz--)
708 if (*(pi++) > *(pj++))
709 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
710 "uplo", "U");
711 } else {
712 while (nnz--)
713 if (*(pi++) < *(pj++))
714 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
715 "uplo", "L");
716 }
717 }
718
719 return Rf_ScalarLogical(1);
720}
721
722SEXP tTMatrix_validate(SEXP obj)
723{
724 if (DIAG(obj) == 'N')
725 return sTMatrix_validate(obj);
726
727 SEXP i = GET_SLOT(obj, Matrix_iSym);
728 R_xlen_t nnz = XLENGTH(i);
729 if (nnz > 0) {
730 PROTECT(i);
731
732 char ul = UPLO(obj);
733
734 SEXP j = GET_SLOT(obj, Matrix_jSym);
735 int *pi = INTEGER(i), *pj = INTEGER(j);
736
737 UNPROTECT(1); /* i */
738 if (ul == 'U') {
739 while (nnz--) {
740 if (*pi > *pj)
741 RMKMS(_("%s=\"%s\" but there are entries below the diagonal"),
742 "uplo", "U");
743 if (*pi == *pj)
744 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
745 "diag", "U");
746 ++pi;
747 ++pj;
748 }
749 } else {
750 while (nnz--) {
751 if (*pi < *pj)
752 RMKMS(_("%s=\"%s\" but there are entries above the diagonal"),
753 "uplo", "L");
754 if (*pi == *pj)
755 RMKMS(_("%s=\"%s\" but there are entries on the diagonal"),
756 "diag", "U");
757 ++pi;
758 ++pj;
759 }
760 }
761 }
762
763 return Rf_ScalarLogical(1);
764}
765
766SEXP xgCMatrix_validate(SEXP obj)
767{
768 SEXP x = PROTECT(GET_SLOT(obj, Matrix_xSym)),
769 i = PROTECT(GET_SLOT(obj, Matrix_iSym));
770 UNPROTECT(2); /* i, x */
771 if (XLENGTH(x) != XLENGTH(i))
772 RMKMS(_("'%s' and '%s' slots do not have equal length"), "i", "x");
773 return Rf_ScalarLogical(1);
774}
775
776SEXP xsCMatrix_validate(SEXP obj)
777{
778 SEXP val = xgCMatrix_validate(obj);
779 if (TYPEOF(val) != STRSXP)
780 val = sCMatrix_validate(obj);
781 return val;
782}
783
784SEXP xtCMatrix_validate(SEXP obj)
785{
786 SEXP val = xgCMatrix_validate(obj);
787 if (TYPEOF(val) != STRSXP)
788 val = tCMatrix_validate(obj);
789 return val;
790}
791
792SEXP xgRMatrix_validate(SEXP obj)
793{
794 SEXP x = PROTECT(GET_SLOT(obj, Matrix_xSym)),
795 j = PROTECT(GET_SLOT(obj, Matrix_jSym));
796 UNPROTECT(2); /* j, x */
797 if (XLENGTH(x) != XLENGTH(j))
798 RMKMS(_("'%s' and '%s' slots do not have equal length"), "j", "x");
799 return Rf_ScalarLogical(1);
800}
801
802SEXP xsRMatrix_validate(SEXP obj)
803{
804 SEXP val = xgRMatrix_validate(obj);
805 if (TYPEOF(val) != STRSXP)
806 val = sRMatrix_validate(obj);
807 return val;
808}
809
810SEXP xtRMatrix_validate(SEXP obj)
811{
812 SEXP val = xgRMatrix_validate(obj);
813 if (TYPEOF(val) != STRSXP)
814 val = tRMatrix_validate(obj);
815 return val;
816}
817
818SEXP xgTMatrix_validate(SEXP obj)
819{
820 SEXP x = PROTECT(GET_SLOT(obj, Matrix_xSym)),
821 i = PROTECT(GET_SLOT(obj, Matrix_iSym));
822 UNPROTECT(2); /* i, x */
823 if (XLENGTH(x) != XLENGTH(i))
824 RMKMS(_("'%s' and '%s' slots do not have equal length"), "i", "x");
825 return Rf_ScalarLogical(1);
826}
827
828SEXP xsTMatrix_validate(SEXP obj)
829{
830 SEXP val = xgTMatrix_validate(obj);
831 if (TYPEOF(val) != STRSXP)
832 val = sTMatrix_validate(obj);
833 return val;
834}
835
836SEXP xtTMatrix_validate(SEXP obj)
837{
838 SEXP val = xgTMatrix_validate(obj);
839 if (TYPEOF(val) != STRSXP)
840 val = tTMatrix_validate(obj);
841 return val;
842}
843
844/* NB: Non-finite entries are "valid" because we consider
845 crossprod(x) and tcrossprod(x) to be positive semidefinite
846 even if 'x' contains non-finite entries (for speed) ...
847*/
848
849SEXP xpoMatrix_validate(SEXP obj)
850{
851 int n = DIM(obj)[1], j;
852 R_xlen_t n1a = (R_xlen_t) n + 1;
853
854 SEXP x = GET_SLOT(obj, Matrix_xSym);
855 if (TYPEOF(x) == REALSXP) {
856 double *px = REAL(x);
857 for (j = 0; j < n; ++j, px += n1a)
858 if (!ISNAN(*px) && *px < 0.0)
859 RMK(_("matrix has negative diagonal elements"));
860 } else {
861 if (TRANS(obj) != 'C')
862 RMKMS(_("'%s' slot is not \"%s\""), "trans", "C");
863 Rcomplex *px = COMPLEX(x);
864 for (j = 0; j < n; ++j, px += n1a)
865 if (!ISNAN((*px).r) && (*px).r < 0.0)
866 RMK(_("matrix has diagonal elements with negative real part"));
867 }
868
869 return Rf_ScalarLogical(1);
870}
871
872SEXP xppMatrix_validate(SEXP obj)
873{
874 int n = DIM(obj)[1], j;
875 char ul = UPLO(obj);
876
877 SEXP x = GET_SLOT(obj, Matrix_xSym);
878 if (TYPEOF(x) == REALSXP) {
879 double *px = REAL(x);
880 if (ul == 'U') {
881 for (j = 0; j < n; px += (++j)+1)
882 if (!ISNAN(*px) && *px < 0.0)
883 RMK(_("matrix has negative diagonal elements"));
884 } else {
885 for (j = 0; j < n; px += n-(j++))
886 if (!ISNAN(*px) && *px < 0.0)
887 RMK(_("matrix has negative diagonal elements"));
888 }
889 } else {
890 if (TRANS(obj) != 'C')
891 RMKMS(_("'%s' slot is not \"%s\""), "trans", "C");
892 Rcomplex *px = COMPLEX(x);
893 if (ul == 'U') {
894 for (j = 0; j < n; px += (++j)+1)
895 if (!ISNAN((*px).r) && (*px).r < 0.0)
896 RMK(_("matrix has diagonal elements with negative real part"));
897 } else {
898 for (j = 0; j < n; px += n-(j++))
899 if (!ISNAN((*px).r) && (*px).r < 0.0)
900 RMK(_("matrix has diagonal elements with negative real part"));
901 }
902 }
903
904 return Rf_ScalarLogical(1);
905}
906
907SEXP xpCMatrix_validate(SEXP obj)
908{
909 int n = DIM(obj)[1], j;
910 char ul = UPLO(obj);
911
912 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
913 i = PROTECT(GET_SLOT(obj, Matrix_iSym)),
914 x = PROTECT(GET_SLOT(obj, Matrix_xSym));
915 UNPROTECT(3); /* x, i, p */
916 int *pp = INTEGER(p), *pi = INTEGER(i);
917 if (TYPEOF(x) == REALSXP) {
918 double *px = REAL(x);
919 if (ul == 'U') {
920 for (j = 0; j < n; ++j)
921 if (pp[j + 1] - pp[j] > 0 && pi[pp[j + 1] - 1] == j &&
922 !ISNAN(px[pp[j + 1] - 1]) && px[pp[j + 1] - 1] < 0.0)
923 RMK(_("matrix has negative diagonal elements"));
924 } else {
925 for (j = 0; j < n; ++j)
926 if (pp[j + 1] - pp[j] > 0 && pi[pp[j]] == j &&
927 !ISNAN(px[pp[j]]) && px[pp[j]] < 0.0)
928 RMK(_("matrix has negative diagonal elements"));
929 }
930 } else {
931 if (TRANS(obj) != 'C')
932 RMKMS(_("'%s' slot is not \"%s\""), "trans", "C");
933 Rcomplex *px = COMPLEX(x);
934 if (ul == 'U') {
935 for (j = 0; j < n; ++j)
936 if (pp[j + 1] - pp[j] > 0 && pi[pp[j + 1] - 1] == j &&
937 !ISNAN(px[pp[j + 1] - 1].r) && px[pp[j + 1] - 1].r < 0.0)
938 RMK(_("matrix has diagonal elements with negative real part"));
939 } else {
940 for (j = 0; j < n; ++j)
941 if (pp[j + 1] - pp[j] > 0 && pi[pp[j]] == j &&
942 !ISNAN(px[pp[j]].r) && px[pp[j]].r < 0.0)
943 RMK(_("matrix has diagonal elements with negative real part"));
944 }
945 }
946
947 return Rf_ScalarLogical(1);
948}
949
950SEXP xpRMatrix_validate(SEXP obj)
951{
952 int m = DIM(obj)[0], i;
953 char ul = UPLO(obj);
954
955 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
956 j = PROTECT(GET_SLOT(obj, Matrix_jSym)),
957 x = PROTECT(GET_SLOT(obj, Matrix_xSym));
958 UNPROTECT(3); /* x, j, p */
959 int *pp = INTEGER(p), *pj = INTEGER(j);
960 if (TYPEOF(x) == REALSXP) {
961 double *px = REAL(x);
962 if (ul == 'U') {
963 for (i = 0; i < m; ++i)
964 if (pp[i + 1] - pp[i] > 0 && pj[pp[i]] == i &&
965 !ISNAN(px[pp[i]]) && px[pp[i]] < 0.0)
966 RMK(_("matrix has negative diagonal elements"));
967 } else {
968 for (i = 0; i < m; ++i)
969 if (pp[i + 1] - pp[i] > 0 && pj[pp[i + 1] - 1] == i &&
970 !ISNAN(px[pp[i + 1] - 1]) && px[pp[i + 1] - 1] < 0.0)
971 RMK(_("matrix has negative diagonal elements"));
972 }
973 } else {
974 if (TRANS(obj) != 'C')
975 RMKMS(_("'%s' slot is not \"%s\""), "trans", "C");
976 Rcomplex *px = COMPLEX(x);
977 if (ul == 'U') {
978 for (i = 0; i < m; ++i)
979 if (pp[i + 1] - pp[i] > 0 && pj[pp[i]] == i &&
980 !ISNAN(px[pp[i]].r) && px[pp[i]].r < 0.0)
981 RMK(_("matrix has diagonal elements with negative real part"));
982 } else {
983 for (i = 0; i < m; ++i)
984 if (pp[i + 1] - pp[i] > 0 && pj[pp[i + 1] - 1] == i &&
985 !ISNAN(px[pp[i + 1] - 1].r) && px[pp[i + 1] - 1].r < 0.0)
986 RMK(_("matrix has diagonal elements with negative real part"));
987 }
988 }
989
990 return Rf_ScalarLogical(1);
991}
992
993SEXP xpTMatrix_validate(SEXP obj)
994{
995 int n = DIM(obj)[1];
996
997 SEXP i = PROTECT(GET_SLOT(obj, Matrix_iSym)),
998 j = PROTECT(GET_SLOT(obj, Matrix_jSym)),
999 x = PROTECT(GET_SLOT(obj, Matrix_xSym));
1000 UNPROTECT(3); /* x, j, i */
1001 int *pi = INTEGER(i), *pj = INTEGER(j);
1002 R_xlen_t nnz = XLENGTH(x);
1003
1004 double *work;
1005 int lwork = n;
1006 Matrix_Calloc(work, lwork, double);
1007 if (TYPEOF(x) == REALSXP) {
1008 double *px = REAL(x);
1009 while (nnz--) {
1010 if (*pi == *pj)
1011 work[*pi] += *px;
1012 ++pi; ++pj; ++px;
1013 }
1014 while (n--)
1015 if (!ISNAN(work[n]) && work[n] < 0.0)
1016 FRMK(_("matrix has negative diagonal elements"));
1017 } else {
1018 if (TRANS(obj) != 'C')
1019 RMKMS(_("'%s' slot is not \"%s\""), "trans", "C");
1020 Rcomplex *px = COMPLEX(x);
1021 while (nnz--) {
1022 if (*pi == *pj)
1023 work[*pi] += (*px).r;
1024 ++pi; ++pj; ++px;
1025 }
1026 while (n--)
1027 if (!ISNAN(work[n]) && work[n] < 0.0)
1028 FRMK(_("matrix has diagonal elements with negative real part"));
1029 }
1030 Matrix_Free(work, lwork);
1031
1032 return Rf_ScalarLogical(1);
1033}
1034
1036{
1037 int n = DIM(obj)[1], j;
1038 R_xlen_t n1a = (R_xlen_t) n + 1;
1039
1040 SEXP x = GET_SLOT(obj, Matrix_xSym);
1041 double *px = REAL(x);
1042 for (j = 0; j < n; ++j, px += n1a)
1043 if (ISNAN(*px) || *px != 1.0)
1044 RMK(_("matrix has nonunit diagonal elements"));
1045
1046 SEXP sd = GET_SLOT(obj, Matrix_sdSym);
1047 if (TYPEOF(sd) != REALSXP)
1048 RMKMS(_("'%s' slot is not of type \"%s\""), "sd", "double");
1049 if (XLENGTH(sd) != n)
1050 RMKMS(_("'%s' slot does not have length %s"), "sd", "Dim[1]");
1051 double *psd = REAL(sd);
1052 for (j = 0; j < n; ++j)
1053 if (!ISNAN(psd[j]) && psd[j] < 0.0)
1054 RMKMS(_("'%s' slot has negative elements"), "sd");
1055
1056 return Rf_ScalarLogical(1);
1057}
1058
1060{
1061 int n = DIM(obj)[1], j;
1062 char ul = UPLO(obj);
1063
1064 SEXP x = GET_SLOT(obj, Matrix_xSym);
1065 double *px = REAL(x);
1066 if (ul == 'U') {
1067 for (j = 0; j < n; px += (++j)+1)
1068 if (ISNAN(*px) || *px != 1.0)
1069 RMK(_("matrix has nonunit diagonal elements"));
1070 } else {
1071 for (j = 0; j < n; px += n-(j++))
1072 if (ISNAN(*px) || *px != 1.0)
1073 RMK(_("matrix has nonunit diagonal elements"));
1074 }
1075
1076 SEXP sd = GET_SLOT(obj, Matrix_sdSym);
1077 if (TYPEOF(sd) != REALSXP)
1078 RMKMS(_("'%s' slot is not of type \"%s\""), "sd", "double");
1079 if (XLENGTH(sd) != n)
1080 RMKMS(_("'%s' slot does not have length %s"), "sd", "Dim[1]");
1081 double *psd = REAL(sd);
1082 for (j = 0; j < n; ++j)
1083 if (!ISNAN(psd[j]) && psd[j] < 0.0)
1084 RMKMS(_("'%s' slot has negative elements"), "sd");
1085
1086 return Rf_ScalarLogical(1);
1087}
1088
1090{
1091 SEXP length = GET_SLOT(obj, Matrix_lengthSym);
1092 if (TYPEOF(length) != INTSXP && TYPEOF(length) != REALSXP)
1093 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1094 "length", "integer", "double");
1095 if (XLENGTH(length) != 1)
1096 RMKMS(_("'%s' slot does not have length %d"), "length", 1);
1097 int_fast64_t n;
1098 if (TYPEOF(length) == INTSXP) {
1099 int n_ = INTEGER(length)[0];
1100 if (n_ == NA_INTEGER)
1101 RMKMS(_("'%s' slot is NA"), "length");
1102 if (n_ < 0)
1103 RMKMS(_("'%s' slot is negative"), "length");
1104 n = (int_fast64_t) n_;
1105 } else {
1106 double n_ = REAL(length)[0];
1107 if (ISNAN(n_))
1108 RMKMS(_("'%s' slot is NA"), "length");
1109 if (n_ < 0.0)
1110 RMKMS(_("'%s' slot is negative"), "length");
1111 if (n_ > 0x1.0p+53)
1112 RMKMS(_("'%s' slot exceeds %s"), "length", "2^53");
1113 n = (int_fast64_t) n_;
1114 }
1115
1116 SEXP i = GET_SLOT(obj, Matrix_iSym);
1117 if (TYPEOF(i) != INTSXP && TYPEOF(i) != REALSXP)
1118 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1119 "i", "integer", "double");
1120 R_xlen_t nnz = XLENGTH(i);
1121 if (nnz > n)
1122 RMKMS(_("'%s' slot has length greater than '%s' slot"), "i", "length");
1123 if (TYPEOF(i) == INTSXP) {
1124 int *pi = INTEGER(i), max = (n > INT_MAX) ? INT_MAX : (int) n, last = 0;
1125 while (nnz--) {
1126 if (*pi == NA_INTEGER)
1127 RMKMS(_("'%s' slot contains NA"), "i");
1128 if (*pi < 1 || *pi > max)
1129 RMKMS(_("'%s' slot has elements not in {%s}"),
1130 "i", "1,...,length");
1131 if (*pi <= last)
1132 RMKMS(_("'%s' slot is not increasing"), "i");
1133 last = *(pi++);
1134 }
1135 } else {
1136 double *pi = REAL(i), max = (double) n, last = 0.0, tmp;
1137 while (nnz--) {
1138 if (ISNAN(*pi))
1139 RMKMS(_("'%s' slot contains NA"), "i");
1140 tmp = trunc(*(pi++));
1141 if (tmp < 1.0 || tmp > max)
1142 RMKMS(_("'%s' slot has elements not in {%s} after truncation towards zero"),
1143 "i", "1,...,length");
1144 if (tmp <= last)
1145 RMKMS(_("'%s' slot is not increasing after truncation towards zero"), "i");
1146 last = tmp;
1147 }
1148 }
1149
1150 return Rf_ScalarLogical(1);
1151}
1152
1153#define KINDVECTOR_VALIDATE(_PREFIX_, _SEXPTYPE_) \
1154SEXP _PREFIX_ ## sparseVector_validate(SEXP obj) \
1155{ \
1156 SEXP x = PROTECT(GET_SLOT(obj, Matrix_xSym)), \
1157 i = PROTECT(GET_SLOT(obj, Matrix_iSym)); \
1158 UNPROTECT(2); /* i, x */ \
1159 if (TYPEOF(x) != _SEXPTYPE_) \
1160 RMKMS(_("'%s' slot is not of type \"%s\""), \
1161 "x", Rf_type2char(_SEXPTYPE_)); \
1162 if (XLENGTH(x) != XLENGTH(i)) \
1163 RMKMS(_("'%s' and '%s' slots do not have equal length"), \
1164 "i", "x"); \
1165 return Rf_ScalarLogical(1); \
1166}
1167KINDVECTOR_VALIDATE(l, LGLSXP)
1168KINDVECTOR_VALIDATE(i, INTSXP)
1169KINDVECTOR_VALIDATE(d, REALSXP)
1170KINDVECTOR_VALIDATE(z, CPLXSXP)
1171#undef KINDVECTOR_VALIDATE
1172
1174{
1175 return Matrix_validate(obj);
1176}
1177
1179{
1180 int *pdim = DIM(obj), n = pdim[1];
1181 if (pdim[0] != n)
1182 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
1183 int_fast64_t nn = (int_fast64_t) n * n;
1184
1185 SEXP x = GET_SLOT(obj, Matrix_xSym);
1186 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1187 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1188 "x", "double", "complex");
1189 if (XLENGTH(x) != nn && XLENGTH(x) != 0)
1190 RMKMS(_("'%s' slot does not have length %s or length %s"),
1191 "x", "prod(Dim)", "0");
1192 int normal = n > 0 && XLENGTH(x) == 0;
1193
1194 SEXP vectors = GET_SLOT(obj, Matrix_vectorsSym);
1195 if (TYPEOF(vectors) != TYPEOF(x))
1196 RMKMS(_("'%s' and '%s' slots do not have the same type"),
1197 "x", "vectors");
1198 if (XLENGTH(vectors) != nn && XLENGTH(vectors) != 0)
1199 RMKMS(_("'%s' slot does not have length %s or length %s"),
1200 "vectors", "prod(Dim)", "0");
1201
1202 SEXP values = GET_SLOT(obj, Matrix_valuesSym);
1203 if (TYPEOF(values) != REALSXP) {
1204 if (normal)
1205 RMKMS(_("'%s' slot is not of type \"%s\""),
1206 "values", "double");
1207 else if (TYPEOF(x) != CPLXSXP)
1208 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1209 "values", "double", "complex");
1210 }
1211 if (XLENGTH(values) != n)
1212 RMKMS(_("'%s' slot does not have length %s"), "values", "Dim[1]");
1213
1214 return Rf_ScalarLogical(1);
1215}
1216
1217SEXP denseQR_validate(SEXP obj)
1218{
1219 int *pdim = DIM(obj), m = pdim[0], n = pdim[1], r = (m < n) ? m : n;
1220
1221 SEXP x = GET_SLOT(obj, Matrix_xSym);
1222 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1223 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1224 "x", "double", "complex");
1225 if (XLENGTH(x) != (int_fast64_t) m * n)
1226 RMKMS(_("'%s' slot does not have length %s"), "x", "prod(Dim)");
1227
1228 SEXP beta = GET_SLOT(obj, Matrix_betaSym);
1229 if (TYPEOF(beta) != REALSXP && TYPEOF(beta) != CPLXSXP)
1230 RMKMS(_("'%s' and '%s' slots do not have the same type"),
1231 "x", "beta");
1232 if (XLENGTH(beta) != r)
1233 RMKMS(_("'%s' slot does not have length %s"), "beta", "min(Dim)");
1234
1235 SEXP perm = GET_SLOT(obj, Matrix_permSym);
1236 if (TYPEOF(perm) != INTSXP)
1237 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
1238 if (XLENGTH(perm) != n)
1239 RMKMS(_("'%s' slot does not have length %s"), "perm", "Dim[2]");
1240 char *work;
1241 int lwork = n;
1242 Matrix_Calloc(work, lwork, char);
1243 int j, *pperm = INTEGER(perm);
1244 for (j = 0; j < n; ++j) {
1245 if (*pperm == NA_INTEGER)
1246 FRMKMS(_("'%s' slot contains NA"), "perm");
1247 if (*pperm < 1 || *pperm > n)
1248 FRMKMS(_("'%s' slot has elements not in {%s}"),
1249 "perm", "1,...,Dim[2]");
1250 if (work[*pperm - 1])
1251 FRMKMS(_("'%s' slot contains duplicates"), "perm");
1252 work[*(pperm++)] = 1;
1253 }
1254 Matrix_Free(work, lwork);
1255
1256 return Rf_ScalarLogical(1);
1257}
1258
1259SEXP denseLU_validate(SEXP obj)
1260{
1261 int *pdim = DIM(obj), m = pdim[0], n = pdim[1], r = (m < n) ? m : n;
1262
1263 SEXP x = GET_SLOT(obj, Matrix_xSym);
1264 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1265 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1266 "x", "double", "complex");
1267 if (XLENGTH(x) != (int_fast64_t) m * n)
1268 RMKMS(_("'%s' slot does not have length %s"), "x", "prod(Dim)");
1269
1270 SEXP perm = GET_SLOT(obj, Matrix_permSym);
1271 if (TYPEOF(perm) != INTSXP)
1272 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
1273 if (XLENGTH(perm) != r)
1274 RMKMS(_("'%s' slot does not have length %s"), "perm", "min(Dim)");
1275 int *pperm = INTEGER(perm);
1276 while (r--) {
1277 if (*pperm == NA_INTEGER)
1278 RMKMS(_("'%s' slot contains NA"), "perm");
1279 if (*pperm < 1 || *pperm > m)
1280 RMKMS(_("'%s' slot has elements not in {%s}"),
1281 "perm", "1,...,Dim[1]");
1282 ++pperm;
1283 }
1284
1285 return Rf_ScalarLogical(1);
1286}
1287
1289{
1290 int *pdim = DIM(obj), n = pdim[1];
1291 if (pdim[0] != n)
1292 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
1293
1294 SEXP uplo = GET_SLOT(obj, Matrix_uploSym);
1295 if (TYPEOF(uplo) != STRSXP)
1296 RMKMS(_("'%s' slot is not of type \"%s\""), "uplo", "character");
1297 if (XLENGTH(uplo) != 1)
1298 RMKMS(_("'%s' slot does not have length %d"), "uplo", 1);
1299 const char *ul = CHAR(STRING_ELT(uplo, 0));
1300 if (ul[0] == '\0' || ul[1] != '\0' || (ul[0] != 'U' && ul[0] != 'L'))
1301 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "uplo", "U", "L");
1302
1303 if (HAS_SLOT(obj, Matrix_transSym)) {
1304 SEXP trans = GET_SLOT(obj, Matrix_transSym);
1305 if (TYPEOF(trans) != STRSXP)
1306 RMKMS(_("'%s' slot is not of type \"%s\""), "trans", "character");
1307 if (XLENGTH(trans) != 1)
1308 RMKMS(_("'%s' slot does not have length %d"), "trans", 1);
1309 const char *ct = CHAR(STRING_ELT(trans, 0));
1310 if (ct[0] == '\0' || ct[1] != '\0' || (ct[0] != 'C' && ct[0] != 'T'))
1311 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "trans", "C", "T");
1312 }
1313
1314 SEXP x = GET_SLOT(obj, Matrix_xSym);
1315 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1316 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1317 "x", "double", "complex");
1318 int packed = XLENGTH(x) != (int_fast64_t) n * n;
1319 if (packed && XLENGTH(x) != n + ((int_fast64_t) n * (n - 1)) / 2)
1320 RMKMS(_("'%s' slot does not have length %s or length %s"),
1321 "x", "prod(Dim)", "Dim[1]*(Dim[1]+1)/2");
1322
1323 SEXP perm = GET_SLOT(obj, Matrix_permSym);
1324 if (TYPEOF(perm) != INTSXP)
1325 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
1326 if (XLENGTH(perm) != n)
1327 RMKMS(_("'%s' slot does not have length %s"), "perm", "Dim[1]");
1328 int n_ = n, *pperm = INTEGER(perm);
1329 while (n_) {
1330 if (*pperm == NA_INTEGER)
1331 RMKMS(_("'%s' slot contains NA"), "perm");
1332 if (*pperm < -n || *pperm == 0 || *pperm > n)
1333 RMKMS(_("'%s' slot has elements not in {%s}\\{%s}"),
1334 "perm", "-Dim[1],...,Dim[1]", "0");
1335 if (*pperm > 0) {
1336 pperm += 1;
1337 n_ -= 1;
1338 } else if (n_ > 1 && *(pperm + 1) == *pperm) {
1339 pperm += 2;
1340 n_ -= 2;
1341 } else
1342 RMKMS(_("'%s' slot has unpaired negative elements"), "perm");
1343 }
1344
1345 return Rf_ScalarLogical(1);
1346}
1347
1349{
1350 int *pdim = DIM(obj), n = pdim[1];
1351 if (pdim[0] != n)
1352 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
1353
1354 SEXP uplo = GET_SLOT(obj, Matrix_uploSym);
1355 if (TYPEOF(uplo) != STRSXP)
1356 RMKMS(_("'%s' slot is not of type \"%s\""), "uplo", "character");
1357 if (XLENGTH(uplo) != 1)
1358 RMKMS(_("'%s' slot does not have length %d"), "uplo", 1);
1359 const char *ul = CHAR(STRING_ELT(uplo, 0));
1360 if (ul[0] == '\0' || ul[1] != '\0' || (ul[0] != 'U' && ul[0] != 'L'))
1361 RMKMS(_("'%s' slot is not \"%s\" or \"%s\""), "uplo", "U", "L");
1362
1363 SEXP x = GET_SLOT(obj, Matrix_xSym);
1364 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1365 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1366 "x", "double", "complex");
1367 int j, packed = XLENGTH(x) != (int_fast64_t) n * n;
1368 if (packed && XLENGTH(x) != n + ((int_fast64_t) n * (n - 1)) / 2)
1369 RMKMS(_("'%s' slot does not have length %s or length %s"),
1370 "x", "prod(Dim)", "Dim[1]*(Dim[1]+1)/2");
1371
1372 /* Real, non-negative diagonal elements are necessary and sufficient */
1373 if (TYPEOF(x) == REALSXP) {
1374 double *px = REAL(x);
1375 if (!packed) {
1376 R_xlen_t n1a = (R_xlen_t) n + 1;
1377 for (j = 0; j < n; ++j, px += n1a)
1378 if (!ISNAN(*px) && *px < 0.0)
1379 RMK(_("Cholesky factor has negative diagonal elements"));
1380 } else if (*ul == 'U') {
1381 for (j = 0; j < n; ++j, px += (++j)+1)
1382 if (!ISNAN(*px) && *px < 0.0)
1383 RMK(_("Cholesky factor has negative diagonal elements"));
1384 } else {
1385 for (j = 0; j < n; ++j, px += n-(j++))
1386 if (!ISNAN(*px) && *px < 0.0)
1387 RMK(_("Cholesky factor has negative diagonal elements"));
1388 }
1389 } else {
1390 Rcomplex *px = COMPLEX(x);
1391 if (!packed) {
1392 R_xlen_t n1a = (R_xlen_t) n + 1;
1393 for (j = 0; j < n; ++j, px += n1a)
1394 if (!ISNAN((*px).r) && (*px).r < 0.0)
1395 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1396 } else if (*ul == 'U') {
1397 for (j = 0; j < n; ++j, px += (++j)+1)
1398 if (!ISNAN((*px).r) && (*px).r < 0.0)
1399 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1400 } else {
1401 for (j = 0; j < n; ++j, px += n-(j++))
1402 if (!ISNAN((*px).r) && (*px).r < 0.0)
1403 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1404 }
1405 }
1406
1407 SEXP perm = GET_SLOT(obj, Matrix_permSym);
1408 if (TYPEOF(perm) != INTSXP)
1409 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
1410 if (XLENGTH(perm) != n && XLENGTH(perm) != 0)
1411 RMKMS(_("'%s' slot does not have length %s or length %s"),
1412 "perm", "Dim[1]", "0");
1413 if (LENGTH(perm) == n) {
1414 char *work;
1415 int lwork = n;
1416 Matrix_Calloc(work, lwork, char);
1417 int *pperm = INTEGER(perm);
1418 for (j = 0; j < n; ++j) {
1419 if (*pperm == NA_INTEGER)
1420 FRMKMS(_("'%s' slot contains NA"), "perm");
1421 if (*pperm < 0 || *pperm >= n)
1422 FRMKMS(_("'%s' slot has elements not in {%s}"),
1423 "perm", "0,...,Dim[1]-1");
1424 if (work[*pperm])
1425 FRMKMS(_("'%s' slot contains duplicates"), "perm");
1426 work[*(pperm++)] = 1;
1427 }
1428 Matrix_Free(work, lwork);
1429 }
1430
1431 return Rf_ScalarLogical(1);
1432}
1433
1434SEXP sparseQR_validate(SEXP obj)
1435{
1436 /* MJ: assuming for simplicity that 'V' and 'R' slots are formally valid */
1437
1438 int *pdim = DIM(obj), m = pdim[0], n = pdim[1];
1439 if (m < n)
1440 RMK(_("matrix has more columns than rows"));
1441
1442 SEXP beta = GET_SLOT(obj, Matrix_betaSym);
1443 if (TYPEOF(beta) != REALSXP)
1444 RMKMS(_("'%s' slot is not of type \"%s\""), "beta", "double");
1445 if (XLENGTH(beta) != n)
1446 RMKMS(_("'%s' slot does not have length %s"), "beta", "Dim[2]");
1447
1448 SEXP p, i, q;
1449 int *pp, *pi, *pq, j, k, kend, m0;
1450
1451 SEXP V = PROTECT(GET_SLOT(obj, Matrix_VSym));
1452 PROTECT(p = GET_SLOT(V, Matrix_pSym));
1453 PROTECT(i = GET_SLOT(V, Matrix_iSym));
1454 pdim = DIM(V); m0 = pdim[0];
1455 UNPROTECT(3); /* i, p, V */
1456 if (m0 < m)
1457 RMKMS(_("'%s' slot has fewer than %s rows"), "V", "Dim[1]");
1458 if (m0 > m + n)
1459 RMKMS(_("'%s' slot has more than %s rows"), "V", "Dim[1]+Dim[2]");
1460 if (pdim[1] != n)
1461 RMKMS(_("'%s' slot does not have %s columns"), "V", "Dim[2]");
1462 pp = INTEGER(p);
1463 pi = INTEGER(i);
1464 for (j = 0, k = 0; j < n; ++j) {
1465 kend = pp[j + 1];
1466 if (k < kend) {
1467 if (pi[k] < j)
1468 RMKMS(_("'%s' slot must be lower trapezoidal but has entries above the diagonal"), "V");
1469 }
1470 k = kend;
1471 }
1472
1473 SEXP R = PROTECT(GET_SLOT(obj, Matrix_RSym));
1474 PROTECT(p = GET_SLOT(R, Matrix_pSym));
1475 PROTECT(i = GET_SLOT(R, Matrix_iSym));
1476 pdim = DIM(R);
1477 UNPROTECT(3); /* i, p, R */
1478 if (pdim[0] != m0)
1479 RMKMS(_("'%s' slot does not have %s row"), "R", "nrow(V)");
1480 if (pdim[1] != n)
1481 RMKMS(_("'%s' slot does not have %s columns"), "R", "Dim[2]");
1482 pp = INTEGER(p);
1483 pi = INTEGER(i);
1484 for (j = 0, k = 0; j < n; ++j) {
1485 kend = pp[j + 1];
1486 if (k < kend && pi[kend - 1] > j)
1487 RMKMS(_("'%s' slot must be upper trapezoidal but has entries below the diagonal"), "R");
1488 k = kend;
1489 }
1490
1491 PROTECT(p = GET_SLOT(obj, Matrix_pSym));
1492 PROTECT(q = GET_SLOT(obj, Matrix_qSym));
1493 UNPROTECT(2); /* q, p */
1494 if (TYPEOF(p) != INTSXP)
1495 RMKMS(_("'%s' slot is not of type \"%s\""), "p", "integer");
1496 if (TYPEOF(q) != INTSXP)
1497 RMKMS(_("'%s' slot is not of type \"%s\""), "q", "integer");
1498 if (XLENGTH(p) != m0)
1499 RMKMS(_("'%s' slot does not have length %s"), "p", "nrow(V)");
1500 if (XLENGTH(q) != n && XLENGTH(q) != 0)
1501 RMKMS(_("'%s' slot does not have length %s or length %s"),
1502 "q", "Dim[2]", "0");
1503 char *work;
1504 int lwork = m0; /* n <= m <= m0 */
1505 Matrix_Calloc(work, lwork, char);
1506 pp = INTEGER(p);
1507 for (j = 0; j < m0; ++j) {
1508 if (*pp == NA_INTEGER)
1509 FRMKMS(_("'%s' slot contains NA"), "p");
1510 if (*pp < 0 || *pp >= m0)
1511 FRMKMS(_("'%s' slot has elements not in {%s}"),
1512 "p", "0,...,nrow(V)-1");
1513 if (work[*pp])
1514 FRMKMS(_("'%s' slot contains duplicates"), "p");
1515 work[*(pp++)] = 1;
1516 }
1517 if (LENGTH(q) == n) {
1518 pq = INTEGER(q);
1519 for (j = 0; j < n; ++j) {
1520 if (*pq == NA_INTEGER)
1521 FRMKMS(_("'%s' slot contains NA"), "q");
1522 if (*pq < 0 || *pq >= n)
1523 FRMKMS(_("'%s' slot has elements not in {%s}"),
1524 "q", "0,...,Dim[2]-1");
1525 if (!work[*pq])
1526 FRMKMS(_("'%s' slot contains duplicates"), "q");
1527 work[*(pq++)] = 0;
1528 }
1529 }
1530 Matrix_Free(work, lwork);
1531
1532 return Rf_ScalarLogical(1);
1533}
1534
1535SEXP sparseLU_validate(SEXP obj)
1536{
1537 /* MJ: assuming for simplicity that 'L' and 'U' slots are formally valid */
1538
1539 int *pdim = DIM(obj), n = pdim[1];
1540 if (pdim[0] != n)
1541 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
1542 char ul, nu;
1543
1544 SEXP L = PROTECT(GET_SLOT(obj, Matrix_LSym));
1545 ul = UPLO(L);
1546 nu = DIAG(L);
1547 pdim = DIM(L);
1548 UNPROTECT(1); /* L */
1549 if (pdim[0] != n || pdim[1] != n)
1550 RMKMS(_("dimensions of '%s' slot are not identical to '%s'"), "L", "Dim");
1551 if (ul == 'U')
1552 RMKMS(_("'%s' slot is upper (not lower) triangular"), "L");
1553 if (nu == 'N') {
1554 PROTECT(L);
1555 SEXP p = PROTECT(GET_SLOT(L, Matrix_pSym)),
1556 i = PROTECT(GET_SLOT(L, Matrix_iSym)),
1557 x = PROTECT(GET_SLOT(L, Matrix_xSym));
1558 UNPROTECT(4); /* x, i, p, L */
1559
1560 int *pp = INTEGER(p), *pi = INTEGER(i), j, k, kend;
1561 if (TYPEOF(x) == REALSXP) {
1562 double *px = REAL(x);
1563 for (j = 0, k = 0; j < n; ++j) {
1564 kend = pp[j + 1];
1565 if (kend == k || pi[k] != j || px[k] != 1.0)
1566 RMKMS(_("'%s' slot has nonunit diagonal elements"), "L");
1567 k = kend;
1568 }
1569 } else {
1570 Rcomplex *px = COMPLEX(x);
1571 for (j = 0, k = 0; j < n; ++j) {
1572 kend = pp[j + 1];
1573 if (kend == k || pi[k] != j || px[k].r != 1.0 || px[k].i != 0.0)
1574 RMKMS(_("'%s' slot has nonunit diagonal elements"), "L");
1575 k = kend;
1576 }
1577 }
1578 }
1579
1580 SEXP U = PROTECT(GET_SLOT(obj, Matrix_USym));
1581 ul = UPLO(U);
1582 pdim = DIM(U);
1583 UNPROTECT(1); /* U */
1584 if (pdim[0] != n || pdim[1] != n)
1585 RMKMS(_("dimensions of '%s' slot are not identical to '%s'"), "U", "Dim");
1586 if (ul != 'U')
1587 RMKMS(_("'%s' slot is lower (not upper) triangular"), "U");
1588
1589 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
1590 q = PROTECT(GET_SLOT(obj, Matrix_qSym));
1591 UNPROTECT(2); /* q, p */
1592 if (TYPEOF(p) != INTSXP)
1593 RMKMS(_("'%s' slot is not of type \"%s\""), "p", "integer");
1594 if (TYPEOF(q) != INTSXP)
1595 RMKMS(_("'%s' slot is not of type \"%s\""), "q", "integer");
1596 if (XLENGTH(p) != n)
1597 RMKMS(_("'%s' slot does not have length %s"), "p", "Dim[1]");
1598 if (XLENGTH(q) != n && XLENGTH(q) != 0)
1599 RMKMS(_("'%s' slot does not have length %s or length %s"),
1600 "q", "Dim[2]", "0");
1601 char *work;
1602 int lwork = n;
1603 Matrix_Calloc(work, lwork, char);
1604 int j, *pp = INTEGER(p);
1605 for (j = 0; j < n; ++j) {
1606 if (*pp == NA_INTEGER)
1607 FRMKMS(_("'%s' slot contains NA"), "p");
1608 if (*pp < 0 || *pp >= n)
1609 FRMKMS(_("'%s' slot has elements not in {%s}"),
1610 "p", "0,...,Dim[1]-1");
1611 if (work[*pp])
1612 FRMKMS(_("'%s' slot contains duplicates"), "p");
1613 work[*(pp++)] = 1;
1614 }
1615 if (LENGTH(q) == n) {
1616 int *pq = INTEGER(q);
1617 for (j = 0; j < n; ++j) {
1618 if (*pq == NA_INTEGER)
1619 FRMKMS(_("'%s' slot contains NA"), "q");
1620 if (*pq < 0 || *pq >= n)
1621 FRMKMS(_("'%s' slot has elements not in {%s}"),
1622 "q", "0,...,Dim[2]-1");
1623 if (!work[*pq])
1624 FRMKMS(_("'%s' slot contains duplicates"), "q");
1625 work[*(pq++)] = 0;
1626 }
1627 }
1628 Matrix_Free(work, lwork);
1629
1630 return Rf_ScalarLogical(1);
1631}
1632
1634{
1635 int *pdim = DIM(obj), n = pdim[1];
1636 if (pdim[0] != n)
1637 RMKMS(_("%s[1] != %s[2] (matrix is not square)"), "Dim", "Dim");
1638
1639 SEXP ordering = GET_SLOT(obj, Matrix_orderingSym);
1640 if (TYPEOF(ordering) != INTSXP)
1641 RMKMS(_("'%s' slot is not of type \"%s\""), "ordering", "integer");
1642 if (XLENGTH(ordering) != 1)
1643 RMKMS(_("'%s' slot does not have length %d"), "ordering", 1);
1644 int or = INTEGER(ordering)[0];
1645 if (or < 0 || or > 6)
1646 RMKMS(_("'%s' is not in %d:%d"), "ordering", 0, 6);
1647
1648 SEXP perm = GET_SLOT(obj, Matrix_permSym);
1649 if (TYPEOF(perm) != INTSXP)
1650 RMKMS(_("'%s' slot is not of type \"%s\""), "perm", "integer");
1651 if (or == 0) {
1652 if (XLENGTH(perm) != 0)
1653 RMKMS(_("'%s' slot does not have length %d"), "perm", 0);
1654 } else {
1655 if (XLENGTH(perm) != n)
1656 RMKMS(_("'%s' slot does not have length %s"), "perm", "Dim[1]");
1657 char *work;
1658 int lwork = n;
1659 Matrix_Calloc(work, lwork, char);
1660 int j, *pperm = INTEGER(perm);
1661 for (j = 0; j < n; ++j) {
1662 if (*pperm == NA_INTEGER)
1663 FRMKMS(_("'%s' slot contains NA"), "perm");
1664 if (*pperm < 0 || *pperm >= n)
1665 FRMKMS(_("'%s' slot has elements not in {%s}"),
1666 "perm", "0,...,Dim[1]-1");
1667 if (work[*pperm])
1668 FRMKMS(_("'%s' slot contains duplicates"), "perm");
1669 work[*(pperm++)] = 1;
1670 }
1671 Matrix_Free(work, lwork);
1672 }
1673
1674 SEXP colcount = GET_SLOT(obj, Matrix_colcountSym);
1675 if (TYPEOF(colcount) != INTSXP)
1676 RMKMS(_("'%s' slot is not of type \"%s\""), "colcount", "integer");
1677 if (XLENGTH(colcount) != n)
1678 RMKMS(_("'%s' slot does not have length %s"), "colcount", "Dim[2]");
1679 int j, *pcolcount = INTEGER(colcount);
1680 for (j = 0; j < n; ++j) {
1681 if (pcolcount[j] == NA_INTEGER)
1682 RMKMS(_("'%s' slot contains NA"), "colcount");
1683 if (pcolcount[j] < 0 || pcolcount[j] > n - j)
1684 RMKMS(_("%s is not in {%s}"), "colcount[j]", "0,...,Dim[2]-j+1");
1685 }
1686
1687 return Rf_ScalarLogical(1);
1688}
1689
1691{
1692 int pattern = !HAS_SLOT(obj, Matrix_minorSym);
1693 if (pattern)
1694 return Rf_ScalarLogical(1);
1695
1696 int n = DIM(obj)[1];
1697 if (n == INT_MAX)
1698 RMKMS(_("%s is not representable as \"%s\""), "Dim[2]+1", "integer");
1699
1700 SEXP minor = GET_SLOT(obj, Matrix_minorSym);
1701 int mr = INTEGER(minor)[0];
1702 if (mr < 0 || mr > n)
1703 RMKMS(_("'%s' slot is not in {%s}"), "minor", "{0,...,Dim[2]}");
1704
1705 SEXP is_ll = GET_SLOT(obj, Matrix_isllSym);
1706 if (TYPEOF(is_ll) != LGLSXP)
1707 RMKMS(_("'%s' slot is not of type \"%s\""), "is_ll", "logical");
1708 if (XLENGTH(is_ll) != 1)
1709 RMKMS(_("'%s' slot does not have length %d"), "is_ll", 1);
1710 int ll = LOGICAL(is_ll)[0];
1711 if (ll == NA_LOGICAL)
1712 RMKMS(_("'%s' slot is NA"), "is_ll");
1713
1714 SEXP is_monotonic = GET_SLOT(obj, Matrix_ismtSym);
1715 if (TYPEOF(is_monotonic) != LGLSXP)
1716 RMKMS(_("'%s' slot is not of type \"%s\""), "is_monotonic", "logical");
1717 if (XLENGTH(is_monotonic) != 1)
1718 RMKMS(_("'%s' slot does not have length %d"), "is_monotonic", 1);
1719 int mt = LOGICAL(is_monotonic)[0];
1720 if (mt == NA_LOGICAL)
1721 RMKMS(_("'%s' slot is NA"), "is_monotonic");
1722
1723 SEXP p = PROTECT(GET_SLOT(obj, Matrix_pSym)),
1724 i = PROTECT(GET_SLOT(obj, Matrix_iSym)),
1725 x = PROTECT(GET_SLOT(obj, Matrix_xSym)),
1726 nz = PROTECT(GET_SLOT(obj, Matrix_nzSym)),
1727 next = PROTECT(GET_SLOT(obj, Matrix_nextSym)),
1728 prev = PROTECT(GET_SLOT(obj, Matrix_prevSym));
1729 UNPROTECT(6); /* prev, next, nz, x, i, p */
1730
1731 if (TYPEOF(next) != INTSXP)
1732 RMKMS(_("'%s' slot is not of type \"%s\""), "next", "integer");
1733 if (TYPEOF(prev) != INTSXP)
1734 RMKMS(_("'%s' slot is not of type \"%s\""), "prev", "integer");
1735 if (XLENGTH(next) - 2 != n)
1736 RMKMS(_("'%s' slot does not have length %s"), "next", "Dim[2]+2");
1737 if (XLENGTH(prev) - 2 != n)
1738 RMKMS(_("'%s' slot does not have length %s"), "prev", "Dim[2]+2");
1739 int *pnext = INTEGER(next), *pprev = INTEGER(prev),
1740 j1 = pnext[n + 1], j2 = pprev[n], count = n + 1;
1741 while (count--) {
1742 if (j1 < 0 || j1 > n)
1743 RMKMS(_("%s has elements not in {%s}"),
1744 "`next`[-(Dim[2]+1)]", "0,...,Dim[2]");
1745 if (j2 < 0 || j2 > n + 1 || j2 == n)
1746 RMKMS(_("%s has elements not in {%s}\\{%s}"),
1747 "`prev`[-(Dim[2]+2)]", "0,...,Dim[2]+1", "Dim[2]");
1748 if ((count > 1) && mt && (pnext[j1] != j1 + 1 || pprev[j2] != j2 - 1))
1749 RMKMS(_("'%s' slot is %s but columns are not stored in increasing order"),
1750 "is_monotonic", "TRUE");
1751 if ((count >= 1) ? j1 == n : j1 != n)
1752 RMKMS(_("traversal of '%s' slot does not complete in exactly %s steps"),
1753 "next", "length(`next`)");
1754 if ((count >= 1) ? j2 == n + 1 : j2 != n + 1)
1755 RMKMS(_("traversal of '%s' slot does not complete in exactly %s steps"),
1756 "prev", "length(`prev`)");
1757 j1 = pnext[j1];
1758 j2 = pprev[j2];
1759 }
1760 if (j1 != -1)
1761 RMKMS(_("%s is not %d"), "`next`[Dim[2]+1]", -1);
1762 if (j2 != -1)
1763 RMKMS(_("%s is not %d"), "`prev`[Dim[2]+2]", -1);
1764
1765 if (TYPEOF(nz) != INTSXP)
1766 RMKMS(_("'%s' slot is not of type \"%s\""), "nz", "integer");
1767 if (XLENGTH(nz) != n)
1768 RMKMS(_("'%s' slot does not have length %s"), "nz", "Dim[2]");
1769 int j, *pnz = INTEGER(nz);
1770 for (j = 0; j < n; ++j) {
1771 if (pnz[j] == NA_INTEGER)
1772 RMKMS(_("'%s' slot contains NA"), "nz");
1773 if (pnz[j] < 1 || pnz[j] > n - j)
1774 RMKMS(_("%s is not in {%s}"), "nz[j]", "1,...,Dim[1]-j+1");
1775 }
1776
1777 if (TYPEOF(p) != INTSXP)
1778 RMKMS(_("'%s' slot is not of type \"%s\""), "p", "integer");
1779 if (XLENGTH(p) - 1 != n)
1780 RMKMS(_("'%s' slot does not have length %s"), "p", "Dim[2]+1");
1781 j1 = pnext[n + 1];
1782 int *pp = INTEGER(p);
1783 if (pp[j1] != 0)
1784 RMKMS(_("column '%s' is stored first but %s is not 0"), "j", "p[j]");
1785 for (j = 0; j < n; ++j) {
1786 j2 = pnext[j1];
1787 if (pp[j2] == NA_INTEGER)
1788 RMKMS(_("'%s' slot contains NA"), "p");
1789 if (pp[j2] < pp[j1])
1790 RMKMS(_("'%s' slot is not increasing when traversed in stored column order"), "p");
1791 if (pp[j2] - pp[j1] < pnz[j1])
1792 RMKMS(_("'%s' slot allocates fewer than %s elements for column '%s'"),
1793 "i", "nz[j]", "j");
1794 if (pp[j2] - pp[j1] > n - j1)
1795 RMKMS(_("'%s' slot allocates more than %s elements for column '%s'"),
1796 "i", "Dim[2]-j+1", "j");
1797 j1 = j2;
1798 }
1799
1800 if (TYPEOF(i) != INTSXP)
1801 RMKMS(_("'%s' slot is not of type \"%s\""), "i", "integer");
1802 if (XLENGTH(i) != pp[n])
1803 RMKMS(_("'%s' slot does not have length %s"), "i", "p[length(p)]");
1804 int *pi = INTEGER(i), *pi_, k;
1805 j1 = pnext[n + 1];
1806 for (j = 0; j < n; ++j) {
1807 pi_ = pi + pp[j1];
1808 if (pi_[0] != j1)
1809 RMKMS(_("first entry in column '%s' does not have row index '%s'"),
1810 "j", "j");
1811 for (k = 1; k < pnz[j1]; ++k) {
1812 if (pi_[k] == NA_INTEGER)
1813 RMKMS(_("'%s' slot contains NA"), "i");
1814 if (pi_[k] < 0 || pi_[k] >= n)
1815 RMKMS(_("'%s' slot has elements not in {%s}"),
1816 "i", "0,...,Dim[1]-1");
1817 if (pi_[k] <= pi_[k - 1])
1818 RMKMS(_("'%s' slot is not increasing within columns"), "i");
1819 }
1820 j1 = pnext[j1];
1821 }
1822
1823 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1824 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1825 "x", "double", "complex");
1826 if (XLENGTH(x) != pp[n])
1827 RMKMS(_("'%s' slot does not have length %s"), "x", "p[length(p)]");
1828
1829 if (ll) {
1830 /* Real, non-negative diagonal elements are necessary and sufficient */
1831 if (TYPEOF(x) == REALSXP) {
1832 double *px = REAL(x);
1833 for (j = 0; j < n; ++j)
1834 if (!ISNAN(px[pp[j]]) && px[pp[j]] < 0.0)
1835 RMK(_("Cholesky factor has negative diagonal elements"));
1836 } else {
1837 Rcomplex *px = COMPLEX(x);
1838 for (j = 0; j < n; ++j)
1839 if (!ISNAN(px[pp[j]].r) && px[pp[j]].r < 0.0)
1840 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1841 }
1842 }
1843
1844 return Rf_ScalarLogical(1);
1845}
1846
1848{
1849 int pattern = !HAS_SLOT(obj, Matrix_minorSym), n = DIM(obj)[1];
1850
1851 if (!pattern) {
1852 SEXP minor = GET_SLOT(obj, Matrix_minorSym);
1853 int mr = INTEGER(minor)[0];
1854 if (mr < 0 || mr > n)
1855 RMKMS(_("'%s' slot is not in {%s}"), "minor", "{0,...,Dim[2]}");
1856 }
1857
1858 SEXP maxcsize = PROTECT(GET_SLOT(obj, Matrix_maxcsizeSym)),
1859 maxesize = PROTECT(GET_SLOT(obj, Matrix_maxesizeSym)),
1860 super = PROTECT(GET_SLOT(obj, Matrix_superSym)),
1861 pi = PROTECT(GET_SLOT(obj, Matrix_piSym)),
1862 px = PROTECT(GET_SLOT(obj, Matrix_pxSym)),
1863 s = PROTECT(GET_SLOT(obj, Matrix_sSym)),
1864 x = (HAS_SLOT(obj, Matrix_xSym)) ? GET_SLOT(obj, Matrix_xSym) : NULL;
1865 UNPROTECT(6); /* s, px, pi, super, maxesize, maxcsize */
1866
1867 if (TYPEOF(super) != INTSXP)
1868 RMKMS(_("'%s' slot is not of type \"%s\""), "super", "integer");
1869 R_xlen_t nsuper1a = XLENGTH(super);
1870 if (nsuper1a - 1 < ((n > 0) ? 1 : 0))
1871 RMKMS(_("'%s' slot has length less than %d"), "super", 2);
1872 if (nsuper1a - 1 > n)
1873 RMKMS(_("'%s' slot has length greater than %s"), "super", "Dim[2]+1");
1874 int k, nsuper = (int) (nsuper1a - 1), *psuper = INTEGER(super);
1875 if (psuper[0] != 0)
1876 RMKMS(_("first element of '%s' slot is not 0"), "super");
1877 if (psuper[nsuper] != n)
1878 RMKMS(_("last element of '%s' slot is not %s"), "super", "Dim[2]");
1879 for (k = 1; k <= nsuper; ++k) {
1880 if (psuper[k] == NA_INTEGER)
1881 RMKMS(_("'%s' slot contains NA"), "super");
1882 if (psuper[k] <= psuper[k - 1])
1883 RMKMS(_("'%s' slot is not increasing"), "super");
1884 }
1885
1886 if (TYPEOF(pi) != INTSXP)
1887 RMKMS(_("'%s' slot is not of type \"%s\""), "pi", "integer");
1888 if (TYPEOF(px) != INTSXP)
1889 RMKMS(_("'%s' slot is not of type \"%s\""), "px", "integer");
1890 if (XLENGTH(pi) != nsuper1a)
1891 RMKMS(_("'%s' and '%s' slots do not have equal length"), "pi", "super");
1892 if (XLENGTH(px) != nsuper1a)
1893 RMKMS(_("'%s' and '%s' slots do not have equal length"), "px", "super");
1894 int *ppi = INTEGER(pi), *ppx = INTEGER(px), nr, nc, l, ml = 0;
1895 if (ppi[0] != 0)
1896 RMKMS(_("first element of '%s' slot is not 0"), "pi");
1897 if (ppx[0] != 0)
1898 RMKMS(_("first element of '%s' slot is not 0"), "px");
1899 for (k = 1; k <= nsuper; ++k) {
1900 if (ppi[k] == NA_INTEGER)
1901 RMKMS(_("'%s' slot contains NA"), "pi");
1902 if (ppx[k] == NA_INTEGER)
1903 RMKMS(_("'%s' slot contains NA"), "px");
1904 if (ppi[k] <= ppi[k - 1])
1905 RMKMS(_("'%s' slot is not increasing"), "pi");
1906 if (ppx[k] <= ppx[k - 1])
1907 RMKMS(_("'%s' slot is not increasing"), "px");
1908 nr = ppi[k] - ppi[k - 1];
1909 nc = psuper[k] - psuper[k - 1];
1910 if (nr < nc)
1911 RMKMS(_("first differences of '%s' slot are less than those of '%s' slot"),
1912 "pi", "super");
1913 if ((int_fast64_t) nr * nc > INT_MAX)
1914 RMKMS(_("supernode lengths exceed %s"), "2^31-1");
1915 l = nr * nc;
1916 if (ppx[k] - ppx[k - 1] != l)
1917 RMKMS(_("first differences of '%s' slot are not equal to supernode lengths"),
1918 "px");
1919 if (l > ml)
1920 ml = l;
1921 }
1922
1923 if (TYPEOF(s) != INTSXP)
1924 RMKMS(_("'%s' slot is not of type \"%s\""), "s", "integer");
1925 if (XLENGTH(s) != ppi[nsuper])
1926 RMKMS(_("'%s' slot does not have length %s"), "s", "pi[length(pi)]");
1927 int i, j, *ps = INTEGER(s);
1928 for (k = 1; k <= nsuper; ++k) {
1929 nr = ppi[k] - ppi[k-1];
1930 nc = psuper[k] - (j = psuper[k-1]);
1931 for (i = 0; i < nr; ++i) {
1932 if (ps[i] == NA_INTEGER)
1933 RMKMS(_("'%s' slot contains NA"), "s");
1934 if (ps[i] < 0 || ps[i] >= n)
1935 RMKMS(_("'%s' slot has elements not in {%s}"),
1936 "s", "0,...,Dim[1]-1");
1937 if (i < nc) {
1938 if (ps[i] != j + i)
1939 RMKMS(_("'%s' slot is wrong within diagonal blocks (row and column indices do not coincide)"), "s");
1940 } else {
1941 if (ps[i] <= ps[i-1])
1942 RMKMS(_("'%s' slot is not increasing within supernodes"), "s");
1943 }
1944 }
1945 ps += nr;
1946 }
1947
1948 if (TYPEOF(maxcsize) != INTSXP)
1949 RMKMS(_("'%s' slot is not of type \"%s\""), "maxesize", "integer");
1950 if (XLENGTH(maxcsize) != 1)
1951 RMKMS(_("'%s' slot does not have length %d"), "maxcsize", 1);
1952 int mc = INTEGER(maxcsize)[0];
1953 if (mc < 0 || mc > ml)
1954 RMKMS(_("'%s' slot is negative or exceeds maximum supernode length"), "maxcsize");
1955
1956 if (TYPEOF(maxesize) != INTSXP)
1957 RMKMS(_("'%s' slot is not of type \"%s\""), "maxesize", "integer");
1958 if (XLENGTH(maxesize) != 1)
1959 RMKMS(_("'%s' slot does not have length %d"), "maxesize", 1);
1960 int me = INTEGER(maxesize)[0];
1961 if (me < 0 || me > n)
1962 RMKMS(_("'%s' slot is negative or exceeds %s"), "maxesize", "Dim[1]");
1963
1964 if (pattern)
1965 return Rf_ScalarLogical(1);
1966
1967 if (TYPEOF(x) != REALSXP && TYPEOF(x) != CPLXSXP)
1968 RMKMS(_("'%s' slot is not of type \"%s\" or \"%s\""),
1969 "x", "double", "complex");
1970 if (XLENGTH(x) != ppx[nsuper])
1971 RMKMS(_("'%s' slot does not have length %s"), "x", "px[length(px)]");
1972
1973 /* Real, non-negative diagonal elements are necessary and sufficient */
1974 if (TYPEOF(x) == REALSXP) {
1975 double *pu = REAL(x), *pv;
1976 for (k = 0; k < nsuper; ++k) {
1977 nr = ppi[k+1] - ppi[k];
1978 nc = psuper[k+1] - psuper[k];
1979 pv = pu + ppx[k];
1980 for (j = 0; j < nc; ++j) {
1981 if (!ISNAN(*pv) && *pv < 0.0)
1982 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1983 pv += (R_xlen_t) nr + 1;
1984 }
1985 }
1986 } else {
1987 Rcomplex *pu = COMPLEX(x), *pv;
1988 for (k = 0; k < nsuper; ++k) {
1989 nr = ppi[k+1] - ppi[k];
1990 nc = psuper[k+1] - psuper[k];
1991 pv = pu + ppx[k];
1992 for (j = 0; j < nc; ++j) {
1993 if (!ISNAN((*pv).r) && (*pv).r < 0.0)
1994 RMK(_("Cholesky factor has diagonal elements with negative real part"));
1995 pv += (R_xlen_t) nr + 1;
1996 }
1997 }
1998 }
1999
2000 return Rf_ScalarLogical(1);
2001}
2002
2003/* where 'cl' must be an element of 'valid_matrix' */
2004void validObject(SEXP obj, const char *cl)
2005{
2006#ifndef MATRIX_DISABLE_VALIDITY
2007
2008 SEXP status;
2009
2010# define IS_VALID(_CLASS_) \
2011 do { \
2012 status = _CLASS_ ## _validate(obj); \
2013 if (TYPEOF(status) == STRSXP) \
2014 Rf_error(_("invalid class \"%s\" object: %s"), \
2015 cl, CHAR(STRING_ELT(status, 0))); \
2016 } while (0)
2017
2018#define IS_VALID_SPARSE(_C_) \
2019 do { \
2020 IS_VALID(_C_ ## sparseMatrix); \
2021 if (cl[0] == 'n') { \
2022 if (cl[1] == 's') \
2023 IS_VALID(s ## _C_ ## Matrix); \
2024 else if (cl[1] == 't') \
2025 IS_VALID(t ## _C_ ## Matrix); \
2026 } else { \
2027 if (cl[1] == 'g') \
2028 IS_VALID(xg ## _C_ ## Matrix); \
2029 else if (cl[1] == 's' || cl[1] == 'p') \
2030 IS_VALID(xs ## _C_ ## Matrix); \
2031 else if (cl[1] == 't') \
2032 IS_VALID(xt ## _C_ ## Matrix); \
2033 if (cl[1] == 'p') \
2034 IS_VALID(xp ## _C_ ## Matrix); \
2035 } \
2036 } while (0)
2037
2038 IS_VALID(Matrix);
2039
2040 const char *cl_ = cl;
2041 if (cl[0] == 'c')
2042 cl = (cl[2] != 'p') ? "dpoMatrix" : "dppMatrix";
2043 else if (cl[0] == 'p')
2044 cl = "indMatrix";
2045
2046 if (cl[0] == 'i' && cl[1] == 'n' && cl[2] == 'd') {
2047 IS_VALID(indMatrix);
2048 if (cl_[0] == 'p')
2049 IS_VALID(pMatrix);
2050 return;
2051 }
2052
2053 if (cl[0] == 'n' && cl[2] != 'C' && cl[2] != 'R' && cl[2] != 'T')
2054 IS_VALID(nMatrix);
2055 else if (cl[0] == 'l')
2056 IS_VALID(lMatrix);
2057 else if (cl[0] == 'i')
2058 IS_VALID(iMatrix);
2059 else if (cl[0] == 'd')
2060 IS_VALID(dMatrix);
2061 else if (cl[0] == 'z')
2062 IS_VALID(zMatrix);
2063
2064 if (cl[1] == 'g')
2065 IS_VALID(generalMatrix);
2066 else if (cl[1] == 's' || cl[1] == 'p')
2067 IS_VALID(symmetricMatrix);
2068 else if (cl[1] == 't')
2069 IS_VALID(triangularMatrix);
2070 else if (cl[1] == 'd') {
2071 IS_VALID(diagonalMatrix);
2072 return;
2073 }
2074
2075 if (cl[2] == 'C')
2076 IS_VALID_SPARSE(C);
2077 else if (cl[2] == 'R')
2078 IS_VALID_SPARSE(R);
2079 else if (cl[2] == 'T')
2080 IS_VALID_SPARSE(T);
2081 else if (cl[2] != 'p') {
2082 IS_VALID(unpackedMatrix);
2083 if (cl[1] == 'p') {
2084 IS_VALID(xpoMatrix);
2085 if (cl_[0] == 'c')
2086 IS_VALID(corMatrix);
2087 }
2088 } else {
2089 IS_VALID(packedMatrix);
2090 if (cl[1] == 'p') {
2091 IS_VALID(xppMatrix);
2092 if (cl_[0] == 'c')
2093 IS_VALID(copMatrix);
2094 }
2095 }
2096
2097# undef IS_VALID_SPARSE
2098# undef IS_VALID
2099
2100#endif
2101
2102 return;
2103}
#define Matrix_Calloc(p, n, t)
Definition Mdefines.h:45
#define _(String)
Definition Mdefines.h:66
#define DIAG(x)
Definition Mdefines.h:111
#define UPLO(x)
Definition Mdefines.h:101
#define Matrix_Free(p, n)
Definition Mdefines.h:56
#define TRANS(x)
Definition Mdefines.h:106
#define DIM(x)
Definition Mdefines.h:85
#define GET_SLOT(x, name)
Definition Mdefines.h:72
int equalString(SEXP, SEXP, R_xlen_t)
Definition utils.c:28
#define HAS_SLOT(x, name)
Definition Mdefines.h:71
#define TYPEOF(s)
Definition Mdefines.h:123
cholmod_common cl
Definition cholmod-etc.c:6
SEXP Matrix_nextSym
Definition init.c:618
SEXP Matrix_nzSym
Definition init.c:619
SEXP Matrix_transSym
Definition init.c:631
SEXP Matrix_sdSym
Definition init.c:629
SEXP Matrix_permSym
Definition init.c:623
SEXP Matrix_DimSym
Definition init.c:598
SEXP Matrix_factorsSym
Definition init.c:606
SEXP Matrix_betaSym
Definition init.c:603
SEXP Matrix_valuesSym
Definition init.c:633
SEXP Matrix_pxSym
Definition init.c:626
SEXP Matrix_marginSym
Definition init.c:614
SEXP Matrix_prevSym
Definition init.c:625
SEXP Matrix_orderingSym
Definition init.c:621
SEXP Matrix_RSym
Definition init.c:600
SEXP Matrix_xSym
Definition init.c:635
SEXP Matrix_lengthSym
Definition init.c:612
SEXP Matrix_iSym
Definition init.c:607
SEXP Matrix_sSym
Definition init.c:628
SEXP Matrix_jSym
Definition init.c:610
SEXP Matrix_DimNamesSym
Definition init.c:597
SEXP Matrix_VSym
Definition init.c:602
SEXP Matrix_LSym
Definition init.c:599
SEXP Matrix_maxesizeSym
Definition init.c:616
SEXP Matrix_diagSym
Definition init.c:605
SEXP Matrix_superSym
Definition init.c:630
SEXP Matrix_uploSym
Definition init.c:632
SEXP Matrix_minorSym
Definition init.c:617
SEXP Matrix_isllSym
Definition init.c:608
SEXP Matrix_colcountSym
Definition init.c:604
SEXP Matrix_qSym
Definition init.c:627
SEXP Matrix_ismtSym
Definition init.c:609
SEXP Matrix_USym
Definition init.c:601
SEXP Matrix_pSym
Definition init.c:622
SEXP Matrix_vectorsSym
Definition init.c:634
SEXP Matrix_maxcsizeSym
Definition init.c:615
SEXP Matrix_piSym
Definition init.c:624
SEXP triangularMatrix_validate(SEXP obj)
Definition validity.c:254
SEXP MatrixFactorization_validate(SEXP obj)
Definition validity.c:1173
SEXP symmetricMatrix_validate(SEXP obj)
Definition validity.c:181
static char * Dim_validate(SEXP dim)
Definition validity.c:35
static char * DimNames_validate(SEXP dimnames, int *pdim)
Definition validity.c:57
SEXP tCMatrix_validate(SEXP obj)
Definition validity.c:558
SEXP denseBunchKaufman_validate(SEXP obj)
Definition validity.c:1288
#define RMK(_FORMAT_)
Definition validity.c:6
SEXP TsparseMatrix_validate(SEXP obj)
Definition validity.c:399
SEXP corMatrix_validate(SEXP obj)
Definition validity.c:1035
SEXP pMatrix_validate(SEXP obj)
Definition validity.c:496
SEXP RsparseMatrix_validate(SEXP obj)
Definition validity.c:349
SEXP xgCMatrix_validate(SEXP obj)
Definition validity.c:766
SEXP diagonalMatrix_validate(SEXP obj)
Definition validity.c:437
SEXP copMatrix_validate(SEXP obj)
Definition validity.c:1059
#define FRMKMS(_FORMAT_,...)
Definition validity.c:23
SEXP tRMatrix_validate(SEXP obj)
Definition validity.c:645
void validObject(SEXP obj, const char *cl)
Definition validity.c:2004
SEXP sCMatrix_validate(SEXP obj)
Definition validity.c:519
SEXP simplicialCholesky_validate(SEXP obj)
Definition validity.c:1690
#define RMS(_FORMAT_,...)
Definition validity.c:8
SEXP xpCMatrix_validate(SEXP obj)
Definition validity.c:907
SEXP sparseQR_validate(SEXP obj)
Definition validity.c:1434
SEXP xppMatrix_validate(SEXP obj)
Definition validity.c:872
SEXP sparseVector_validate(SEXP obj)
Definition validity.c:1089
SEXP packedMatrix_validate(SEXP obj)
Definition validity.c:290
SEXP sRMatrix_validate(SEXP obj)
Definition validity.c:606
SEXP xgRMatrix_validate(SEXP obj)
Definition validity.c:792
SEXP sparseCholesky_validate(SEXP obj)
Definition validity.c:1633
SEXP R_Dim_validate(SEXP dim)
Definition validity.c:50
#define FRMK(_FORMAT_)
Definition validity.c:13
SEXP sparseLU_validate(SEXP obj)
Definition validity.c:1535
SEXP xsTMatrix_validate(SEXP obj)
Definition validity.c:828
SEXP Matrix_validate(SEXP obj)
Definition validity.c:136
SEXP R_DimNames_fixup(SEXP dimnames)
Definition validity.c:95
SEXP xtCMatrix_validate(SEXP obj)
Definition validity.c:784
SEXP supernodalCholesky_validate(SEXP obj)
Definition validity.c:1847
#define RMKMS(_FORMAT_,...)
Definition validity.c:10
#define KINDMATRIX_VALIDATE(_PREFIX_, _SEXPTYPE_)
Definition validity.c:149
SEXP R_DimNames_validate(SEXP dimnames, SEXP dim)
Definition validity.c:89
#define KINDVECTOR_VALIDATE(_PREFIX_, _SEXPTYPE_)
Definition validity.c:1153
#define IS_VALID_SPARSE(_C_)
SEXP xpoMatrix_validate(SEXP obj)
Definition validity.c:849
SEXP xgTMatrix_validate(SEXP obj)
Definition validity.c:818
SEXP denseQR_validate(SEXP obj)
Definition validity.c:1217
SEXP unpackedMatrix_validate(SEXP obj)
Definition validity.c:281
SEXP generalMatrix_validate(SEXP obj)
Definition validity.c:165
SEXP denseSchur_validate(SEXP obj)
Definition validity.c:1178
SEXP xtTMatrix_validate(SEXP obj)
Definition validity.c:836
SEXP xsRMatrix_validate(SEXP obj)
Definition validity.c:802
SEXP CsparseMatrix_validate(SEXP obj)
Definition validity.c:299
SEXP xpRMatrix_validate(SEXP obj)
Definition validity.c:950
SEXP tTMatrix_validate(SEXP obj)
Definition validity.c:722
SEXP xtRMatrix_validate(SEXP obj)
Definition validity.c:810
SEXP indMatrix_validate(SEXP obj)
Definition validity.c:461
SEXP denseLU_validate(SEXP obj)
Definition validity.c:1259
SEXP xpTMatrix_validate(SEXP obj)
Definition validity.c:993
SEXP denseCholesky_validate(SEXP obj)
Definition validity.c:1348
SEXP sTMatrix_validate(SEXP obj)
Definition validity.c:693
#define IS_VALID(_CLASS_)
SEXP xsCMatrix_validate(SEXP obj)
Definition validity.c:776