SCM

[#6610] Matrix::rankMatrix(, method = "qr") segfaults

Date:
2019-01-16 16:48
Priority:
3
State:
Closed
Submitted by:
Benjamin Tyner (btyner)
Assigned to:
Martin Maechler (mmaechler)
Hardware:
None
Product:
None
Operating System:
Linux
Component:
None
Version:
None
Severity:
None
Resolution:
Fixed
URL:
Summary:
Matrix::rankMatrix(, method = "qr") segfaults

Detailed description
Apologies if this is a duplicate of one of the existing segfault reports, but for version 1.2-15 of Matrix, if I run:

set.seed(6860)
n <- 2669041L

ni <- 296236L
nj <- 54L

dat <- data.frame(i = structure(sample(seq_len(ni), size = n, replace = TRUE),
.Label = as.character(seq_len(ni)),
class = "factor"
),
j = structure(sample(seq_len(nj), size = n, replace = TRUE),
.Label = as.character(seq_len(nj)),
class = "factor"
)
)

X <- MatrixModels::model.Matrix(~0+i+j, data = dat, sparse = TRUE)

library(Matrix)
XtX <- crossprod(X)

# this segfaults!!!
rnk <- Matrix::rankMatrix(XtX, method = "qr")

then it gives:

*** caught segfault ***
address 0x7f597f74f26c, cause 'memory not mapped'

Traceback:
1: .Call(dgCMatrix_QR, x, if (verbose) -1L else TRUE, keep.dimnames)
2: .local(x, ...)
3: qr(as(as(as(x, "CsparseMatrix"), "dsparseMatrix"), "dgCMatrix"), ...)
4: qr(as(as(as(x, "CsparseMatrix"), "dsparseMatrix"), "dgCMatrix"), ...)
5: qr(if (do.t) t(x) else x, tol = tol, LAPACK = method != "qrLINPACK")
6: qr(if (do.t) t(x) else x, tol = tol, LAPACK = method != "qrLINPACK")
7: structure(if (useGrad) which.min(diff1) else if (Meth == "qr") { if ((do.t <- (d[1L] < d[2L])) && warn.t) warning(gettextf("rankMatrix(x, method='qr'): computing t(x) as nrow(x) < ncol(x)")) q.r <- qr(if (do.t) t(x) else x, tol = tol, LAPACK = method != "qrLINPACK") if (x.dense && (method == "qrLINPACK")) q.r$rank else { diagR <- if (x.dense) diag(q.r$qr) else diag(q.r@R) d.i <- abs(diagR) if ((mdi <- max(d.i)) > 0) sum(d.i >= tol * mdi) else 0L }} else if (sval[1] > 0) sum(sval >= tol * sval[1]) else 0L, method = method, useGrad = useGrad, tol = if (useGrad) NA else tol)
8: Matrix::rankMatrix(XtX, method = "qr")
An irrecoverable exception occurred. R is aborting now ...
Segmentation fault

Comments:

Message  ↓
Date: 2022-08-31 10:38
Sender: Martin Maechler

Thank you. We are now giving a warning "Out of memory" and an error in this situation.

If anybody finds a reproducible example that does not involve the 'MatrixModels' package, I'd be glad, so could add it to the -- only_if_enough_GB_RAM regression checks.

Date: 2022-08-21 14:39
Sender: Mikael Jagan

I can reproduce with the latest Matrix (r3597). The cause seems to be an integer overflow inside of cs_multiply(), at line 1272 of src/cs.c (more precisely, when j == 6330):

if (nz + m > C->nzmax && !cs_sprealloc (C, 2*(C->nzmax)+m))

which does not check whether 2*(C->nzmax) or 2*(C->nzmax)+m would exceed INT_MAX. Interestingly, cs_sprealloc() does not return 0 when its second argument 'nzmax' is negative; it silently uses the value of p[n] instead, which may very well be less than what is needed by the caller.

I've just checked a version our package with cs_sprealloc() modified to return 0 instead of "coping". No regressions, but optimally we would contact the author of the library to find out if such a change is really safe.

Date: 2021-01-05 16:13
Sender: Martin Maechler

Thank you, Benjamin.
Somehow, this got completely forgotten (and I'm not even sure if I ever read the report in full ...), I'm really sorry.

The bug still triggers the segfault, also for me with current
R-devel and "next version" Matrix.

Actually, I'd be grateful for help here.... as the new year, and upcoming new semester means a lot of other duties to be taken care of ...

Date: 2019-01-17 03:43
Sender: Benjamin Tyner

Another way to trigger is with:

QR <- qr(XtX)

valgrind has this to say:

==25984== Warning: set address range perms: large range [0x6decd55c, 0x82078ffc) (undefined)
==25984== Warning: set address range perms: large range [0x41d7d028, 0x55e07574) (noaccess)
==25984== Warning: set address range perms: large range [0xaa2afffc, 0xd260753c) (undefined)
==25984== Warning: set address range perms: large range [0x59e43028, 0x82079014) (noaccess)
==25984== Warning: set address range perms: large range [0x122b9553c, 0x173243fbc) (undefined)
==25984== Warning: set address range perms: large range [0x8207a028, 0xd2607554) (noaccess)
==25984== Warning: set address range perms: large range [0x213e7ffbc, 0x2b4bdd4bc) (undefined)
==25984== Warning: set address range perms: large range [0xd2608028, 0x173243fd4) (noaccess)
==25984== Conditional jump or move depends on uninitialised value(s)
==25984== at 0x4C31CE5: realloc (in /usr/lib/valgrind/vgpreload_memcheck-amd64-linux.so)
==25984== by 0x12AABB91: cs_realloc (cs.c:1158)
==25984== by 0x12AACA9C: cs_sprealloc (cs.c:1933)
==25984== by 0x12AAD8A2: cs_multiply (cs.c:1271)
==25984== by 0x12AAEC9E: cs_amd (cs.c:82)
==25984== by 0x12AB04B3: cs_sqr (cs.c:1748)
==25984== by 0x12AB5EAC: dgCMatrix_QR (dgCMatrix.c:233)
==25984== by 0x1F5735: do_dotcall (dotcode.c:1252)
==25984== by 0x239F0C: Rf_eval (eval.c:719)
==25984== by 0x23B37E: R_execClosure (eval.c:1765)
==25984== by 0x230686: bcEval (eval.c:6743)
==25984== by 0x23994F: Rf_eval (eval.c:616)
==25984==
==25984== Warning: set address range perms: large range [0x173244028, 0x2b4bdd4d4) (noaccess)
==25984== Invalid write of size 4
==25984== at 0x12AAC25C: cs_scatter (cs.c:1583)
==25984== by 0x12AAD836: cs_multiply (cs.c:1279)
==25984== by 0x12AAEC9E: cs_amd (cs.c:82)
==25984== by 0x12AB04B3: cs_sqr (cs.c:1748)
==25984== by 0x12AB5EAC: dgCMatrix_QR (dgCMatrix.c:233)
==25984== by 0x1F5735: do_dotcall (dotcode.c:1252)
==25984== by 0x239F0C: Rf_eval (eval.c:719)
==25984== by 0x23B37E: R_execClosure (eval.c:1765)
==25984== by 0x230686: bcEval (eval.c:6743)
==25984== by 0x23994F: Rf_eval (eval.c:616)
==25984== by 0x23B37E: R_execClosure (eval.c:1765)
==25984== by 0x23C2D1: R_execMethod (eval.c:1935)
==25984== Address 0x14cad953c is not stack'd, malloc'd or (recently) free'd
==25984==

*** caught segfault ***
address 0x14cad953c, cause 'memory not mapped'

Traceback:
1: .Call(dgCMatrix_QR, x, if (verbose) -1L else TRUE, keep.dimnames)
2: .local(x, ...)
3: qr(as(as(as(x, "CsparseMatrix"), "dsparseMatrix"), "dgCMatrix"), ...)
4: qr(as(as(as(x, "CsparseMatrix"), "dsparseMatrix"), "dgCMatrix"), ...)
5: qr(XtX)
6: qr(XtX)

Attached Files:

Changes

Field Old Value Date By
status_idOpen2022-08-31 10:38mmaechler
close_dateNone2022-08-31 10:38mmaechler
ResolutionAccepted As Bug2022-08-31 10:38mmaechler
assigned_tonone2021-01-05 16:13mmaechler
ResolutionNone2021-01-05 16:13mmaechler
Thanks to:
Vienna University of Economics and Business Powered By FusionForge