View source: R/S3_factorizations.R
| qr.HDF5Matrix | R Documentation |
Computes A = Q R block-wise on disk and returns Q and R as
HDF5Matrix objects.
## S3 method for class 'HDF5Matrix'
qr(
x,
thin = FALSE,
block_size = NULL,
overwrite = FALSE,
threads = -1L,
method = "auto",
compression = NULL,
...
)
x |
An |
thin |
Logical. Compute thin (economy) QR. Default |
block_size |
Integer or NULL. Row-block size hint for TSQR; ignored by
|
overwrite |
Logical. Overwrite existing results. Default |
threads |
Integer. OpenMP threads (-1 = auto, CRAN-compliant).
For |
method |
Character. Algorithm selection:
|
compression |
Integer (0-9) or NULL. gzip compression level for the
result datasets. NULL uses the global option set by
|
... |
Ignored (for S3 compatibility). |
Named list: Q (HDF5Matrix), R (HDF5Matrix).
20260304: Added method parameter and TSQR support.
tmp <- tempfile(fileext = ".h5")
X <- hdf5_create_matrix(tmp, "data/A", data = matrix(rnorm(10000), 100, 100))
X <- hdf5_matrix(tmp, "data/A")
# Default (auto method)
res <- qr(X)
dim(res$Q) # m x m (or m x min(m,n) if thin = TRUE)
dim(res$R) # m x n (or min(m,n) x n)
# Explicit TSQR for a tall-skinny matrix (recommended: thin = TRUE)
hdf5_close_all()
unlink(tmp)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.