Package 'withinr'

Title: High-Performance Fixed Effects Solver
Description: Fast iterative solvers (CG, GMRES) with Schwarz domain-decomposition preconditioners for absorbing high-dimensional fixed effects in panel data regressions. The computational core is written in Rust via the 'within' crate and accessed through 'extendr'.
Authors: Alexander Fischer [aut, cre], Kristof Schroeder [aut]
Maintainer: Alexander Fischer <[email protected]>
License: MIT + file LICENSE
Version: 0.1.0
Built: 2026-05-27 06:09:39 UTC
Source: https://github.com/py-econometrics/within

Help Index


Solve fixed-effects normal equations

Description

Computes fixed-effect coefficients by solving the normal equations DTWDx=DTWyD^T W D x = D^T W y where DD is the dummy-variable design matrix implied by categories and WW is the diagonal weight matrix.

Usage

solve(
  categories,
  y,
  weights = NULL,
  method = c("cg", "gmres"),
  tol = 1e-08,
  maxiter = 1000L,
  restart = 30L,
  preconditioner = c("additive", "multiplicative", "off")
)

Arguments

categories

Integer matrix of shape (n_obs, n_factors). Each column contains 1-based factor level assignments. Values must be positive integers with no NAs.

y

Numeric vector of length n_obs. The response variable.

weights

Numeric vector of length n_obs or NULL (default). Observation weights. NULL means unit weights (unweighted).

method

Character, one of "cg" (default) or "gmres". "cg" uses Conjugate Gradient (requires symmetric preconditioner). "gmres" uses GMRES (supports all preconditioners).

tol

Convergence tolerance on the relative residual norm. Default 1e-8.

maxiter

Maximum number of Krylov iterations. Default 1000L.

restart

GMRES restart parameter (ignored for CG). Default 30L.

preconditioner

Character, one of "additive" (default), "multiplicative", or "off". Schwarz preconditioner variant. "multiplicative" requires method = "gmres".

Value

A named list with components:

coefficients

Numeric vector of fixed-effect coefficient estimates.

demeaned

Numeric vector. Response after subtracting estimated fixed effects.

converged

Logical. Did the solver meet the tolerance?

iterations

Integer. Number of Krylov iterations performed.

residual

Numeric. Final relative residual norm.

time_total

Numeric. Wall-clock seconds (setup + solve).

time_setup

Numeric. Wall-clock seconds for operator/preconditioner construction.

time_solve

Numeric. Wall-clock seconds for the iterative solve.


Solve fixed-effects normal equations for multiple response vectors

Description

Builds the operator and preconditioner once, then solves for each column of Y in parallel. More efficient than calling solve in a loop.

Usage

solve_batch(
  categories,
  Y,
  weights = NULL,
  method = c("cg", "gmres"),
  tol = 1e-08,
  maxiter = 1000L,
  restart = 30L,
  preconditioner = c("additive", "multiplicative", "off")
)

Arguments

categories

Integer matrix of shape (n_obs, n_factors). Each column contains 1-based factor level assignments.

Y

Numeric matrix of shape (n_obs, k). Each column is a separate response vector.

weights

Numeric vector of length n_obs or NULL (default).

method

Character, one of "cg" (default) or "gmres".

tol

Convergence tolerance. Default 1e-8.

maxiter

Maximum Krylov iterations. Default 1000L.

restart

GMRES restart parameter. Default 30L.

preconditioner

Character, one of "additive", "multiplicative", or "off".

Value

A named list with components:

coefficients

Numeric matrix (n_dofs, k).

demeaned

Numeric matrix (n_obs, k).

converged

Logical vector of length k.

iterations

Integer vector of length k.

residual

Numeric vector of length k.

time_total

Numeric scalar. Wall-clock seconds for the entire batch.

time_solve

Numeric vector of length k. Per-RHS solve time.