Skip to contents

`aligned_rrr()` fits a shared reduced-rank response model for multiblock predictors `X` paired with blockwise multivariate responses `Y`. It is the clean supervised baseline against `response_aligned_mfa()`: each block learns a block-specific regression weight matrix `W_k`, all blocks share a common response basis `B`, and latent scores are defined by `Z_k = X_k W_k`.

Usage

aligned_rrr(
  Y,
  X,
  preproc = multivarious::center(),
  response_preproc = multivarious::center(),
  ncomp = 2,
  block_weight = 1,
  response_weights = NULL,
  max_iter = 50,
  tol = 1e-06,
  ridge = 1e-08,
  verbose = FALSE,
  use_future = FALSE,
  ...
)

Arguments

Y

A list of numeric matrices/data.frames. Each element `Y[[k]]` is `n_k x q` and all blocks must share the same response column dimension.

X

A list of numeric matrices/data.frames. Each element `X[[k]]` is `n_k x p_k`.

preproc

A `multivarious` preprocessing pipeline (a `pre_processor` or `prepper`) or a list of them for the `X` blocks.

response_preproc

A single shared `multivarious` preprocessing pipeline used to define the common response space.

ncomp

Integer latent rank. In v1, `aligned_rrr()` requires `ncomp <= ncol(Y[[1]])`.

block_weight

Optional numeric scalar or vector of per-block weights.

response_weights

Optional rowwise response weights. May be `NULL`, a single non-negative scalar, a length-`length(X)` vector of blockwise constants, or a list mirroring `X` with one non-negative weight vector per block row.

max_iter

Maximum ALS iterations.

tol

Relative convergence tolerance on the objective.

ridge

Non-negative ridge stabilization applied to the block regression weight updates.

verbose

Logical; if `TRUE`, prints iteration diagnostics.

use_future

Logical; if `TRUE`, block-wise computations are performed via `furrr::future_map()` when available.

...

Unused (reserved for future extensions).

Value

An object inheriting from `multivarious::multiblock_biprojector` with additional class `"aligned_rrr"`. The object stores block-specific regression weights in `W_list`, shared response basis `B`, block-specific latent scores in `Z_list`, rowwise response weights in `response_weights`, and pooled score normalization weights in `score_weights`. For compatibility with generic score accessors, the `s` slot stores a concatenated stack of `Z_list`, `score_index` maps those rows back to blocks, and `score_representation = "stacked_block_scores"` records that this is not a shared observation-level score table.

Details

The fitted model minimizes $$ \sum_k \beta_k \|R_k^{1/2}(Y_k - X_k W_k B^\top)\|_F^2 $$ plus ridge stabilization on the block regression weights `W_k`. There is no additional predictor-reconstruction term and no anchor machinery.

Identifiability is imposed on the score side: after each iteration, the block score matrices `Z_k = X_k W_k` are orthonormalized in the pooled weighted score space and the shared response basis `B` is rotated accordingly. A deterministic sign convention is then applied componentwise.

Predictors are preprocessed blockwise, while responses are preprocessed in a shared pooled response space. Out-of-sample prediction is pure: `x_new -> z_hat -> y_hat`, with no conditional completion path.

Examples

# \donttest{
set.seed(1)
X1 <- matrix(rnorm(40 * 6), 40, 6)
X2 <- matrix(rnorm(36 * 5), 36, 5)
B <- qr.Q(qr(matrix(rnorm(3 * 2), 3, 2)))
W1 <- matrix(rnorm(6 * 2), 6, 2)
W2 <- matrix(rnorm(5 * 2), 5, 2)
Y1 <- X1 %*% W1 %*% t(B) + matrix(rnorm(40 * 3, sd = 0.05), 40, 3)
Y2 <- X2 %*% W2 %*% t(B) + matrix(rnorm(36 * 3, sd = 0.05), 36, 3)

fit <- aligned_rrr(
  Y = list(X1 = Y1, X2 = Y2),
  X = list(X1 = X1, X2 = X2),
  ncomp = 2
)
#> Applying the same preprocessor definition independently to each block.

pred <- predict(fit, X1[1:5, , drop = FALSE], block = "X1", type = "response")
stopifnot(nrow(pred) == 5)
# }