Population Structure and Relatedness Inference using the GENESIS Package

Overview

GENESIS provides statistical methodology for analyzing genetic data from samples with population structure and/or familial relatedness. This vignette provides a description of how to use GENESIS for inferring population structure, as well as estimating relatedness measures such as kinship coefficients, identity by descent (IBD) sharing probabilities, and inbreeding coefficients. GENESIS uses PC-AiR for population structure inference that is robust to known or cryptic relatedness, and it uses PC-Relate for accurate relatedness estimation in the presence of population structure, admixutre, and departures from Hardy-Weinberg equilibrium.

Relatedness Estimation Adjusted for Principal Components (PC-Relate)

Many estimators exist that use genome-wide SNP genotype data from samples in genetic studies to estimate measures of recent genetic relatedness such as pairwise kinship coefficients, pairwise IBD sharing probabilities, and individual inbreeding coefficients. However, many of these estimators either make simplifying assumptions that do not hold in the presence of population structure and/or ancestry admixture, or they require reference population panels of known ancestry from pre-specified populations. When assumptions are violated, these estimators can provide highly biased estimates.

The PC-Relate method is used to accurately estimate measures of recent genetic relatedness in samples with unknown or unspecified population structure without requiring reference population panels. PC-Relate uses ancestry representative principal components to account for sample ancestry differences and provide estimates that are robust to population structure, ancestry admixture, and departures from Hardy-Weinberg equilibirum.

Data

Reading in Genotype Data

The functions in the GENESIS package can read genotype data from a GenotypeData class object as created by the GWASTools package. Through the use of GWASTools, a GenotypeData class object can easily be created from:

  • an R matrix of SNP genotype data
  • a GDS file
  • PLINK files

Example R code for creating a GenotypeData object is presented below. Much more detail can be found in the GWASTools package reference manual.

GENESIS can also work with genotype data from sequencing, starting with a VCF file. For examples using this format, see the vignette “Analyzing Sequence Data using the GENESIS Package”.

R Matrix

geno <- MatrixGenotypeReader(genotype = genotype, snpID = snpID, 
                             chromosome = chromosome, position = position, 
                             scanID = scanID)
genoData <- GenotypeData(geno)
  • genotype is a matrix of genotype values coded as 0 / 1 / 2, where rows index SNPs and columns index samples
  • snpID is an integer vector of unique SNP IDs
  • chromosome is an integer vector specifying the chromosome of each SNP
  • position is an integer vector specifying the position of each SNP
  • scanID is a vector of unique individual IDs

GDS files

geno <- GdsGenotypeReader(filename = "genotype.gds")
genoData <- GenotypeData(geno)
  • filename is the file path to the GDS object

HapMap Data

To demonstrate PC-AiR and PC-Relate analyses with the GENESIS package, we analyze SNP data from the Mexican Americans in Los Angeles, California (MXL) and African American individuals in the southwestern USA (ASW) population samples of HapMap 3. Mexican Americans and African Americans have a diverse ancestral background, and familial relatives are present in these data. Genotype data at a subset of 20K autosomal SNPs for 173 individuals are provided as a GDS file.

gdsfile <- system.file("extdata", "HapMap_ASW_MXL_geno.gds", package="GENESIS")

Relatedness Estimation Adjusted for Principal Components (PC-Relate)

Running PC-Relate

PC-Relate uses the ancestry representative principal components (PCs) calculated from PC-AiR to adjust for the population structure and ancestry of individuals in the sample and provide accurate estimates of recent genetic relatedness due to family structure. The pcrelate function performs the PC-Relate analysis.

The training.set input of the pcrelate function allows for the specification of which samples are used to estimate the ancestry adjustment at each SNP. The adjustment tends to perform best when close relatives are excluded from training.set, so the individuals in the “unrelated subset” from the PC-AiR analysis are typically a good choice. However, when an “unrelated subset” is not available, the adjustment still works well when estimated using all samples (training.set = NULL).

In order to run PC-Relate, we first need to create an iterator object to read SNPs in blocks. We create the iterator such that only pruned SNPs are returned in each block.

# run PC-Relate
HapMap_genoData <- GenotypeBlockIterator(HapMap_genoData, snpInclude=pruned)
mypcrelate <- pcrelate(HapMap_genoData, pcs = mypcair$vectors[,1:2], 
                       training.set = mypcair$unrels,
                       BPPARAM = BiocParallel::SerialParam())
  • genoData is a GenotypeIterator class object
  • pcs is a matrix whose columns are the PCs used for the ancestry adjustment
  • training.set is a vector of individual IDs specifying which samples are used to esimate the ancestry adjustment at each SNP

If estimates of IBD sharing probabilities are not desired, then setting the input ibd.probs = FALSE will speed up the computation.

Output from PC-Relate

The pcrelate function will return an object of class pcrelate, which is a list of two data.frames: kinBtwn with pairwise kinship values, and kinSelf with inbreeding coefficients.

plot(mypcrelate$kinBtwn$k0, mypcrelate$kinBtwn$kin, xlab="k0", ylab="kinship")

A function is provided for making a genetic relationship matrix (GRM). Using a threshold for kinship will create a sparse matrix by setting kinship for pairs less than the threshold to 0. This can be useful to reduce memory usage for very large sample sizes.

iids <- as.character(getScanID(HapMap_genoData))
pcrelateToMatrix(mypcrelate, sample.include = iids[1:5], thresh = 2^(-11/2), scaleKin = 2)
## 5 x 5 sparse Matrix of class "dsCMatrix"
##           NA19919  NA19916   NA19835   NA20282   NA19703
## NA19919 0.9960379 .        .         .         .        
## NA19916 .         1.009738 .         .         .        
## NA19835 .         .        0.9794893 .         .        
## NA20282 .         .        .         0.9283949 .        
## NA19703 .         .        .         .         0.9753479
  • pcrelobj is the output from pcrelate; either a class pcrelate object or a GDS file
  • sample.include is a vector of individual IDs specifying which individuals to include in the table or matrix
  • thresh specifies a minimum kinship coefficient value to include in the GRM
  • scaleKin specifies a factor to multiply the kinship coefficients by in the GRM

References

  • Conomos M.P., Reiner A.P., Weir B.S., & Thornton T.A. (2016). Model-free Estimation of Recent Genetic Relatedness. American Journal of Human Genetics, 98(1), 127-148.

  • Conomos M.P., Miller M.B., & Thornton T.A. (2015). Robust Inference of Population Structure for Ancestry Prediction and Correction of Stratification in the Presence of Relatedness. Genetic Epidemiology, 39(4), 276-293.

  • Gogarten, S.M., Bhangale, T., Conomos, M.P., Laurie, C.A., McHugh, C.P., Painter, I., … & Laurie, C.C. (2012). GWASTools: an R/Bioconductor package for quality control and analysis of Genome-Wide Association Studies. Bioinformatics, 28(24), 3329-3331.

  • Manichaikul, A., Mychaleckyj, J.C., Rich, S.S., Daly, K., Sale, M., & Chen, W.M. (2010). Robust relationship inference in genome-wide association studies. Bioinformatics, 26(22), 2867-2873.