pedmat: Genetic Relationship Matrices and Inbreeding Coefficients

View source: R/pedmatrix.R

pedmatR Documentation

Genetic Relationship Matrices and Inbreeding Coefficients

Description

Optimized calculation of additive (A), dominance (D), epistatic (AA) relationship matrices, their inverses, and inbreeding coefficients (f). Uses Rcpp with Meuwissen & Luo (1992) algorithm for efficient computation.

Usage

pedmat(
  ped,
  method = "A",
  sparse = TRUE,
  invert_method = "auto",
  threads = 0,
  compact = FALSE
)

Arguments

ped

A tidied pedigree from tidyped. Must be a single pedigree, not a splitped object. For splitped results, use pedmat(ped_split$GP1, ...) to process individual groups.

method

Character, one of:

  • "A": Additive (numerator) relationship matrix (default)

  • "f": Inbreeding coefficients (returns named vector)

  • "Ainv": Inverse of A using Henderson's rules (O(n) complexity)

  • "D": Dominance relationship matrix

  • "Dinv": Inverse of D (requires matrix inversion)

  • "AA": Additive-by-additive epistatic matrix (A # A)

  • "AAinv": Inverse of AA

sparse

Logical, if TRUE returns sparse Matrix (recommended for large pedigrees). Default is TRUE.

invert_method

Character, method for matrix inversion (Dinv/AAinv only):

  • "auto": Auto-detect and use optimal method (default)

  • "sympd": Force Cholesky decomposition (faster for SPD matrices)

  • "general": Force general LU decomposition

threads

Integer. Number of OpenMP threads to use. Use 0 to keep the system/default setting. Currently, multi-threading is explicitly implemented for:

  • "D": Dominance relationship matrix (significant speedup).

  • "Ainv": Inverse of A (only for large pedigrees, n >= 5000).

For "Dinv", "AA", and "AAinv", parallelism depends on the linked BLAS/LAPACK library (e.g., OpenBLAS, MKL, Accelerate) and is not controlled by this parameter. Methods "A" and "f" are single-threaded.

compact

Logical, if TRUE compacts full-sibling families by selecting one representative per family. This dramatically reduces matrix dimensions for pedigrees with large full-sib groups. See Details.

Details

API Design:

Only a single method may be requested per call. This design prevents accidental heavy computations. If multiple matrices are needed, call pedmat() separately for each method.

Compact Mode (compact = TRUE):

Full-siblings share identical relationships with all other individuals. Compact mode exploits this by selecting one representative per full-sib family, dramatically reducing matrix size. For example, a pedigree with 170,000 individuals might compact to 1,800 unique relationship patterns.

Key features:

  • query_relationship(x, id1, id2): Query any individual pair, including merged siblings (automatic lookup)

  • expand_pedmat(x): Restore full matrix dimensions

  • vismat(x): Visualize directly (auto-expands compact)

Performance Notes:

  • Ainv: O(n) complexity using Henderson's rules. Fast for any size.

  • Dinv/AAinv: O(n³) matrix inversion. Practical limits:

    • n < 500: ~10-20 ms

    • n = 1,000: ~40-60 ms

    • n = 2,000: ~130-150 ms

    • n > 2,000: Consider using compact = TRUE

  • Memory: Sparse matrices use ~O(nnz) memory; dense use O(n²)

Value

Returns a matrix or vector with S3 class "pedmat".

Object type by method:

  • method="f": Named numeric vector of inbreeding coefficients

  • All other methods: Sparse or dense matrix (depending on sparse)

S3 Methods:

  • print(x): Display matrix with metadata header

  • summary_pedmat(x): Detailed statistics (size, compression, mean, density)

  • dim(x), length(x), Matrix::diag(x), t(x): Standard operations

  • x[i, j]: Subsetting (behaves like underlying matrix)

  • as.matrix(x): Convert to base matrix

Accessing Metadata (use attr(), not $):

  • attr(x, "ped"): The pedigree used (or compact pedigree if compact=TRUE)

  • attr(x, "ped_compact"): Compact pedigree (when compact=TRUE)

  • attr(x, "method"): Calculation method used

  • attr(x, "call_info"): Full calculation metadata (timing, sizes, etc.)

  • names(attributes(x)): List all available attributes

Additional attributes when compact = TRUE:

  • attr(x, "compact_map"): data.table mapping individuals to representatives

  • attr(x, "family_summary"): data.table summarizing merged families

  • attr(x, "compact_stats"): Compression statistics (ratio, n_families, etc.)

References

Meuwissen, T. H. E., & Luo, Z. (1992). Computing inbreeding coefficients in large populations. Genetics Selection Evolution, 24(4), 305-313.

Henderson, C. R. (1976). A simple method for computing the inverse of a numerator relationship matrix used in prediction of breeding values. Biometrics, 32(1), 69-83.

See Also

tidyped for preparing pedigree data, query_relationship for querying individual pairs, expand_pedmat for restoring full dimensions, vismat for visualization, inbreed for simple inbreeding calculation

Examples

# Basic usage with small pedigree
library(visPedigree)
tped <- tidyped(small_ped)

# --- Additive Relationship Matrix (default) ---
A <- pedmat(tped)
A["A", "B"]      # Relationship between A and B
Matrix::diag(A)  # Diagonal = 1 + F (inbreeding)

# --- Inbreeding Coefficients ---
f <- pedmat(tped, method = "f")
f["Z1"]  # Inbreeding of individual Z1

# --- Using summary_pedmat() ---
summary_pedmat(A)   # Detailed matrix statistics

# --- Accessing Metadata ---
attr(A, "ped")              # Original pedigree
attr(A, "method")           # "A"
names(attributes(A))        # All available attributes

# --- Compact Mode (for large full-sib families) ---
A_compact <- pedmat(tped, method = "A", compact = TRUE)

# Query relationships (works for any individual, including merged sibs)
query_relationship(A_compact, "Z1", "Z2")  # Full-sibs Z1 and Z2

# View compression statistics
attr(A_compact, "compact_stats")
attr(A_compact, "family_summary")

# Expand back to full size
A_full <- expand_pedmat(A_compact)
dim(A_full)  # Original dimensions restored

# --- Inverse Matrices ---
Ainv <- pedmat(tped, method = "Ainv")  # Henderson's rules (fast)

# --- Dominance and Epistatic ---
D <- pedmat(tped, method = "D")
AA <- pedmat(tped, method = "AA")

# --- Visualization (requires display device) ---
## Not run: 
vismat(A)                       # Heatmap of relationship matrix
vismat(A_compact)               # Works with compact matrices
vismat(A, by = "Gen")     # Group by generation

## End(Not run)


visPedigree documentation built on March 30, 2026, 9:07 a.m.