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