qr.HDF5Matrix: QR decomposition of an HDF5Matrix

View source: R/S3_factorizations.R

qr.HDF5MatrixR Documentation

QR decomposition of an HDF5Matrix

Description

Computes A = Q R block-wise on disk and returns Q and R as HDF5Matrix objects.

Usage

## S3 method for class 'HDF5Matrix'
qr(
  x,
  thin = FALSE,
  block_size = NULL,
  overwrite = FALSE,
  threads = -1L,
  method = "auto",
  compression = NULL,
  ...
)

Arguments

x

An HDF5Matrix.

thin

Logical. Compute thin (economy) QR. Default FALSE. For method = "tsqr" this is strongly recommended; full Q via TSQR requires an expensive basis-completion step (O(m² k)).

block_size

Integer or NULL. Row-block size hint for TSQR; ignored by method = "lapack". NULL = auto.

overwrite

Logical. Overwrite existing results. Default FALSE.

threads

Integer. OpenMP threads (-1 = auto, CRAN-compliant). For method = "lapack" controls Eigen BLAS threads. For method = "tsqr" controls the OpenMP block loop.

method

Character. Algorithm selection:

"auto"

Default. Selects TSQR when m > 5n AND m > 1000 (tall-skinny), LAPACK otherwise.

"lapack"

Eigen HouseholderQR. Reliable for any matrix shape.

"tsqr"

Parallel Tall-Skinny QR. Requires m >= n. Efficient for big-omics tall matrices (e.g. samples × features). The R factor is mathematically equivalent to LAPACK R (up to sign). The Q factor is a valid orthogonal matrix but differs from LAPACK Q.

compression

Integer (0-9) or NULL. gzip compression level for the result datasets. NULL uses the global option set by hdf5matrix_options (default 6). Use 0 to disable compression (faster for benchmarks).

...

Ignored (for S3 compatibility).

Value

Named list: Q (HDF5Matrix), R (HDF5Matrix).

Note

20260304: Added method parameter and TSQR support.

Examples



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)




BigDataStatMeth documentation built on May 15, 2026, 1:07 a.m.