genio
The genio
(GenIO = Genetics I/O) package aims to facilitate reading and writing genetics data. The focus of this vignette is processing plink BED/BIM/FAM files.
There are some limited alternatives for reading and/or writing BED files in R, which are slower and harder to use, which motivated me to write this package. Here we make direct comparisons to those packages, to illustrate the advantages of genio
.
Let’s begin by creating a large genotype matrix with completely random data, and with missing values so this becomes non-trivial to read and write.
# Data dimensions.
# Choose non-multiples of 4 to test edge cases of BED parsers.
# Number of loci.
m_loci <- 10001
# Number of individuals
n_ind <- 1001
# Overall allele frequency
# (we don't do any actual genetics in this vignette,
# so just set to a reasonable value)
p <- 0.5
# Missingness rate
miss <- 0.1
# Total number of genotypes
n_data <- n_ind * m_loci
# Draw random genotypes from Binomial
X <- rbinom( n_data, 2, p)
# Add missing values
X[ sample(n_data, n_data * miss) ] <- NA
# Turn into matrix
X <- matrix(X, nrow = m_loci, ncol = n_ind)
# Inspect the first 10 individuals at the first 10 loci
X[1:10, 1:10]
#> [,1] [,2] [,3] [,4] [,5] [,6] [,7] [,8] [,9] [,10]
#> [1,] 2 1 0 2 2 NA 1 1 NA 1
#> [2,] NA 2 0 1 NA 0 1 1 1 2
#> [3,] 1 NA 2 1 1 0 2 NA 1 1
#> [4,] 1 2 1 1 0 2 2 2 0 NA
#> [5,] 0 1 1 NA 2 0 NA 2 1 0
#> [6,] 2 1 NA 0 1 1 0 0 0 2
#> [7,] 1 2 2 0 1 1 0 NA 2 2
#> [8,] 1 1 1 1 NA NA 2 2 0 1
#> [9,] 0 1 0 NA NA 0 2 1 1 1
#> [10,] 2 NA 1 1 1 1 NA 1 2 NA
To create annotation tables that look a bit more interesting than the defaults, let us create some slightly more realistic values. Here we use our first genio
function!
First we create and edit the locus annotations table.
# We have to specify the number of loci
bim <- make_bim( n = m_loci )
# Inspect the default values
bim
#> # A tibble: 10,001 x 6
#> chr id posg pos ref alt
#> <dbl> <int> <dbl> <int> <dbl> <dbl>
#> 1 1 1 0 1 1 2
#> 2 1 2 0 2 1 2
#> 3 1 3 0 3 1 2
#> 4 1 4 0 4 1 2
#> 5 1 5 0 5 1 2
#> 6 1 6 0 6 1 2
#> 7 1 7 0 7 1 2
#> 8 1 8 0 8 1 2
#> 9 1 9 0 9 1 2
#> 10 1 10 0 10 1 2
#> # … with 9,991 more rows
# Let's add the "chr" prefix to the chromosome values,
# so we recognize them when we see them later.
bim$chr <- paste0('chr', bim$chr)
# Make SNP IDs look like "rs" IDs
bim$id <- paste0('rs', bim$id)
# Make positions 1000 bigger
bim$pos <- bim$pos * 1000
# Select randomly between Cs and Gs for the reference alleles
bim$ref <- sample(c('C', 'G'), m_loci, replace = TRUE)
# Select randomly between As and Ts for the alternative alleles
bim$alt <- sample(c('A', 'T'), m_loci, replace = TRUE)
# Inspect the table with our changes
bim
#> # A tibble: 10,001 x 6
#> chr id posg pos ref alt
#> <chr> <chr> <dbl> <dbl> <chr> <chr>
#> 1 chr1 rs1 0 1000 G T
#> 2 chr1 rs2 0 2000 C A
#> 3 chr1 rs3 0 3000 G A
#> 4 chr1 rs4 0 4000 C T
#> 5 chr1 rs5 0 5000 C T
#> 6 chr1 rs6 0 6000 G T
#> 7 chr1 rs7 0 7000 G T
#> 8 chr1 rs8 0 8000 G A
#> 9 chr1 rs9 0 9000 G T
#> 10 chr1 rs10 0 10000 C T
#> # … with 9,991 more rows
Now we similarly create and edit the individual annotations table.
# Specify the number of individuals
fam <- make_fam( n = n_ind )
# Inspect the default values
fam
#> # A tibble: 1,001 x 6
#> fam id pat mat sex pheno
#> <int> <int> <dbl> <dbl> <dbl> <dbl>
#> 1 1 1 0 0 0 0
#> 2 2 2 0 0 0 0
#> 3 3 3 0 0 0 0
#> 4 4 4 0 0 0 0
#> 5 5 5 0 0 0 0
#> 6 6 6 0 0 0 0
#> 7 7 7 0 0 0 0
#> 8 8 8 0 0 0 0
#> 9 9 9 0 0 0 0
#> 10 10 10 0 0 0 0
#> # … with 991 more rows
# Add prefixes to families and IDs to recognize them later
fam$fam <- paste0('fam', fam$fam)
fam$id <- paste0('id', fam$id)
# Sex values are usually 1 and 2
fam$sex <- sample(1:2, n_ind, replace = TRUE)
# Let's make phenotypes continuous.
# Draw independently from Standard Normal.
fam$pheno <- rnorm(n_ind)
# Let's leave maternal and paternal IDs as missing
# Inspect again
fam
#> # A tibble: 1,001 x 6
#> fam id pat mat sex pheno
#> <chr> <chr> <dbl> <dbl> <int> <dbl>
#> 1 fam1 id1 0 0 2 -0.899
#> 2 fam2 id2 0 0 2 1.49
#> 3 fam3 id3 0 0 2 -1.28
#> 4 fam4 id4 0 0 2 -0.0518
#> 5 fam5 id5 0 0 1 0.455
#> 6 fam6 id6 0 0 2 -0.431
#> 7 fam7 id7 0 0 2 0.342
#> 8 fam8 id8 0 0 2 1.11
#> 9 fam9 id9 0 0 2 -0.235
#> 10 fam10 id10 0 0 2 1.68
#> # … with 991 more rows
Lastly, let’s copy the locus and individual IDs as row and column names of the genotype matrix, respectively. Although this step is not required, it is encouraged for consistency (if present, values are checked when writing files).
# Add column and row names from bim and fam tables we just created.
rownames(X) <- bim$id
colnames(X) <- fam$id
# Inspect again the first 10 individuals and loci
X[1:10, 1:10]
#> id1 id2 id3 id4 id5 id6 id7 id8 id9 id10
#> rs1 2 1 0 2 2 NA 1 1 NA 1
#> rs2 NA 2 0 1 NA 0 1 1 1 2
#> rs3 1 NA 2 1 1 0 2 NA 1 1
#> rs4 1 2 1 1 0 2 2 2 0 NA
#> rs5 0 1 1 NA 2 0 NA 2 1 0
#> rs6 2 1 NA 0 1 1 0 0 0 2
#> rs7 1 2 2 0 1 1 0 NA 2 2
#> rs8 1 1 1 1 NA NA 2 2 0 1
#> rs9 0 1 0 NA NA 0 2 1 1 1
#> rs10 2 NA 1 1 1 1 NA 1 2 NA
Let’s write this random data to a file. This mode is intended for simulated data, generating dummy BIM and FAM files in the process and writing out all three.
# Will delete at the end of the vignette
file_plink <- tempfile('vignette-random-data')
# Write genotypes, along with the BIM and FAM files we created.
# Omiting them would result in writing the original dummy version of these tables,
# before we edited them.
time_write_genio <- system.time(
write_plink(file_plink, X, bim, fam)
)
#> Writing: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.bed
#> Writing: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.fam
#> Writing: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.bim
time_write_genio
#> user system elapsed
#> 0.151 0.014 0.164
Here we demonstrate how easy it is to read the data back. To compare to other packages, we shall time loading all this data.
# Read the data back in memory.
# Time this step
time_read_genio <- system.time(
data_genio <- read_plink(file_plink)
)
#> Reading: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.bim
#> Reading: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.fam
#> Reading: /tmp/RtmpMsCEXK/vignette-random-data778e5ef14ce6.bed
time_read_genio
#> user system elapsed
#> 0.078 0.013 0.091
# Inspect the data we just read back
# The same random genotypes (first 10 individuals and loci, now with row and column names):
data_genio$X[1:10, 1:10]
#> id1 id2 id3 id4 id5 id6 id7 id8 id9 id10
#> rs1 2 1 0 2 2 NA 1 1 NA 1
#> rs2 NA 2 0 1 NA 0 1 1 1 2
#> rs3 1 NA 2 1 1 0 2 NA 1 1
#> rs4 1 2 1 1 0 2 2 2 0 NA
#> rs5 0 1 1 NA 2 0 NA 2 1 0
#> rs6 2 1 NA 0 1 1 0 0 0 2
#> rs7 1 2 2 0 1 1 0 NA 2 2
#> rs8 1 1 1 1 NA NA 2 2 0 1
#> rs9 0 1 0 NA NA 0 2 1 1 1
#> rs10 2 NA 1 1 1 1 NA 1 2 NA
# The locus annotations
data_genio$bim
#> # A tibble: 10,001 x 6
#> chr id posg pos ref alt
#> <chr> <chr> <dbl> <int> <chr> <chr>
#> 1 chr1 rs1 0 1000 G T
#> 2 chr1 rs2 0 2000 C A
#> 3 chr1 rs3 0 3000 G A
#> 4 chr1 rs4 0 4000 C T
#> 5 chr1 rs5 0 5000 C T
#> 6 chr1 rs6 0 6000 G T
#> 7 chr1 rs7 0 7000 G T
#> 8 chr1 rs8 0 8000 G A
#> 9 chr1 rs9 0 9000 G T
#> 10 chr1 rs10 0 10000 C T
#> # … with 9,991 more rows
# The individual annotations
data_genio$fam
#> # A tibble: 1,001 x 6
#> fam id pat mat sex pheno
#> <chr> <chr> <chr> <chr> <int> <dbl>
#> 1 fam1 id1 0 0 2 -0.899
#> 2 fam2 id2 0 0 2 1.49
#> 3 fam3 id3 0 0 2 -1.28
#> 4 fam4 id4 0 0 2 -0.0518
#> 5 fam5 id5 0 0 1 0.455
#> 6 fam6 id6 0 0 2 -0.431
#> 7 fam7 id7 0 0 2 0.342
#> 8 fam8 id8 0 0 2 1.11
#> 9 fam9 id9 0 0 2 -0.235
#> 10 fam10 id10 0 0 2 1.68
#> # … with 991 more rows
# Quickly test that the inputs and outputs are identical.
# Genotypes have NAs, so we have to compare this way.
stopifnot( all( X == data_genio$X, na.rm = TRUE) )
stopifnot( bim == data_genio$bim )
# FAM has mixed content (chars, ints, and doubles).
# First 5 columns should be exact match:
stopifnot( fam[,1:5] == data_genio$fam[,1:5] )
# Exact equality may fail for pheno due to precisions, so test this way instead:
stopifnot( max(abs(fam$pheno - data_genio$fam$pheno)) < 1e-4 )
One drawback of genio
and other related approaches is high memory consumption (see more on that at the end of this vignette). It is entirely possible that an attempt to parse a genotype matrix into memory will fail with an “out of memory” error message. Let’s estimate memory usage here.
Genotypes are stored as integers by genio
, which in R on a modern machine (64 bit architecture) consumes 4 bytes! So an \(n \times m\) matrix takes up at least \(4 n m\) bytes (an R matrix contains an additional constant overhead, which does not depend on \(n,m\) and is relatively small for large matrices).
To get an idea of how much this is, let’s assume a standard genotyping array, where \(m\) is about half a million. In this scenario, we get a simple rule of thumb that every 1000 individuals consume a little less than 2G:
# Constants
bytes_per_genotype <- 4
bytes_per_gb <- 1024 ^ 3
# Example data dimensions
num_ind <- 1000
num_loci <- 500000
# Gigabytes per 1000 individuals for a typical genotyping array
bytes_per_genotype * num_ind * num_loci / bytes_per_gb
#> [1] 1.862645
This is the amount of free memory required just to load the matrix. Several common matrix operations for genotypes consume even more memory, so more free memory will be needed to accomplish anything.
There are several R packages that provide BED readers, and to a lesser extent some of the other functionality of genio
. Here we compare directly to BEDMatrix
, snpStats
, and lfa
. Each of these packages has different goals so they are optimized for their use cases. The genio
parser is optimized to read the entire genotype data into a native R matrix, so that it is easy to inspect and manipulate for R beginners. Here we demonstrate that genio
is not only the fastest at this task, but also the easiest to obtain in terms of requiring the least amount of coding.
The BEDMatrix
package allows access to genotypes from a BED file without loading the entire matrix in memory. This is very convenient for large datasets, as long as only small portions of the data are pulled at once. However, there are some disadvantages:
BEDMatrix
return object is not a regular R matrix, which confuses some users.BEDMatrix
package provides no way to access the full annotation tables (BIM and FAM files).
BEDMatrix
does not have any write functions.Here is an example of that usage and how the data compares. Note in particular that the BEDMatrix
representation of the genotype matrix is transposed compared to the genio
matrix:
library(BEDMatrix)
# Time it too.
# Although the BIM and FAM tables are not returned,
# they are partially parsed and kept in memory,
# which can take time for extremely large files
time_read_bedmatrix_1 <- system.time(
X_BEDMatrix <- BEDMatrix(file_plink)
)
#> Extracting number of samples and rownames from vignette-random-data778e5ef14ce6.fam...
#> Extracting number of variants and colnames from vignette-random-data778e5ef14ce6.bim...
# Inspect the first 10 loci and individuals as usual.
# Note it is transposed compared to our X.
# Also note the column and row names are different from genio's.
X_BEDMatrix[1:10, 1:10]
#> rs1_G rs2_C rs3_G rs4_C rs5_C rs6_G rs7_G rs8_G rs9_G rs10_C
#> fam1_id1 2 NA 1 1 0 2 1 1 0 2
#> fam2_id2 1 2 NA 2 1 1 2 1 1 NA
#> fam3_id3 0 0 2 1 1 NA 2 1 0 1
#> fam4_id4 2 1 1 1 NA 0 0 1 NA 1
#> fam5_id5 2 NA 1 0 2 1 1 NA NA 1
#> fam6_id6 NA 0 0 2 0 1 1 NA 0 1
#> fam7_id7 1 1 2 2 NA 0 0 2 2 NA
#> fam8_id8 1 1 NA 2 2 0 NA 2 1 1
#> fam9_id9 NA 1 1 0 1 0 2 0 1 2
#> fam10_id10 1 2 1 NA 0 2 2 1 1 NA
Therefore, for very large datasets, if it suffices to access the genotype data in small portions and the user is willing to deal with the limitations of the object, and no detailed annotation table information is required, then BEDMatrix
is a better solution than the read_bed
function provided in the genio
package.
To compare the two genotype reading functions on an equal footing, let us assume that the entire genotype matrix is required in memory and stored in a regular R matrix.
# This turns it into a regular R matrix.
# Since most of the reading is actually happening now,
# we time this step now.
time_read_bedmatrix_2 <- system.time(
X_BEDMatrix_Rmat <- as.matrix(X_BEDMatrix)
)
time_read_bedmatrix_2
#> user system elapsed
#> 0.070 0.007 0.078
# Now we can test that the BEDMatrix agrees with the original matrix we simulated.
# Need to transpose first.
stopifnot( all( X == t(X_BEDMatrix_Rmat), na.rm = TRUE) )
The snpStats
package is the most fully-featured alternative to genio
. One of its advantages is its memory efficiency, encoding the entire data in less memory than a native R matrix. However, this memory efficiency comes at the cost of making the data harder to access, especially for inexperienced R users. Like genio
, and in contrast to the other packages, snpStats
also reads the BIM and FAM tables fully, and provides a BED writer, but it is considerably harder to use.
First we illustrate data parsing. The annotation tables are similar but column names are different, and certain missing values (zero in text files) are converted to NA
instead.
library(snpStats)
#> Loading required package: survival
#> Loading required package: Matrix
# Read data, time it.
time_read_snpStats_1 <- system.time(
data_snpStats <- read.plink(file_plink)
)
time_read_snpStats_1
#> user system elapsed
#> 0.030 0.003 0.032
# Inspect the data
# Genotypes.
# Note data is not visualized this way.
# This matrix is also transposed compared to the genio matrix.
data_snpStats$genotypes
#> A SnpMatrix with 1001 rows and 10001 columns
#> Row names: fam1 ... fam1001
#> Col names: rs1 ... rs10001
# Locus annotations
head( data_snpStats$map )
#> chromosome snp.name cM position allele.1 allele.2
#> rs1 chr1 rs1 NA 1000 G T
#> rs2 chr1 rs2 NA 2000 C A
#> rs3 chr1 rs3 NA 3000 G A
#> rs4 chr1 rs4 NA 4000 C T
#> rs5 chr1 rs5 NA 5000 C T
#> rs6 chr1 rs6 NA 6000 G T
# Individual annotations
head (data_snpStats$fam )
#> pedigree member father mother sex affected
#> fam1 fam1 id1 NA NA 2 -0.89945429
#> fam2 fam2 id2 NA NA 2 1.48550249
#> fam3 fam3 id3 NA NA 2 -1.27863368
#> fam4 fam4 id4 NA NA 2 -0.05182039
#> fam5 fam5 id5 NA NA 1 0.45474930
#> fam6 fam6 id6 NA NA 2 -0.43133813
As for BEDMatrix
, assuming we ultimately desire to convert the entire data into a regular R matrix, an extra step is required:
# Transpose, then convert to a regular R matrix.
# Let's time this step too.
time_read_snpStats_2 <- system.time(
X_snpStats <- as( t(data_snpStats$genotypes), 'numeric')
)
time_read_snpStats_2
#> user system elapsed
#> 0.071 0.015 0.086
# Now we can visualize the matrix.
# First 10 loci and individuals, as before.
# Note that, compared to (genio, BEDMatrix, lfa), alleles are encoded in reverse,
# so 0s and 2s are flipped in this matrix.
X_snpStats[1:10, 1:10]
#> fam1 fam2 fam3 fam4 fam5 fam6 fam7 fam8 fam9 fam10
#> rs1 0 1 2 0 0 NA 1 1 NA 1
#> rs2 NA 0 2 1 NA 2 1 1 1 0
#> rs3 1 NA 0 1 1 2 0 NA 1 1
#> rs4 1 0 1 1 2 0 0 0 2 NA
#> rs5 2 1 1 NA 0 2 NA 0 1 2
#> rs6 0 1 NA 2 1 1 2 2 2 0
#> rs7 1 0 0 2 1 1 2 NA 0 0
#> rs8 1 1 1 1 NA NA 0 0 2 1
#> rs9 2 1 2 NA NA 2 0 1 1 1
#> rs10 0 NA 1 1 1 1 NA 1 0 NA
# Again verify that the matrices are identical.
# (Here 2-X flips back 0s and 2s)
stopifnot( all( X == 2 - X_snpStats, na.rm = TRUE) )
snpStats
is the only package other than genio
to provide a BED writer. Here we demonstrate how hard it is to use it to write our data.
# Let's write this to another file
file_plink_copy <- tempfile('vignette-random-data-copy')
# Copy objects to not change originals
X_snpStats <- X
bim_snpStats <- as.data.frame(bim) # to use rownames
fam_snpStats <- as.data.frame(fam) # ditto
# All data requires matching row and/or column names.
# These first two were already done above.
#rownames(X_snpStats) <- bim$id
#colnames(X_snpStats) <- fam$id
# Row names here are redundant but required.
rownames(bim_snpStats) <- bim$id
rownames(fam_snpStats) <- fam$id
# We shall time several required steps in order to write genotypes in a standard R matrix,
# and the related annotation tables, to BED.
time_write_snpStats <- system.time({
# Transpose and convert our genotypes to SnpMatrix object.
# We flip 0s and 2s before converting
X_snpStats <- as(2 - t(X_snpStats), 'SnpMatrix')
# This complicated command is required to write the data.
# Although X, fam, and bim are passed as for genio's write_plink,
# here the name of every column must be specified (there are no reasonable defaults).
# Interestingly, the parameter names of snpStats' write.plink
# do not agree with read.plink from the same package.
write.plink(
file_plink_copy,
snps = X_snpStats,
subject.data = fam_snpStats,
pedigree = fam,
id = id,
father = pat,
mother = mat,
sex = sex,
phenotype = pheno,
snp.data = bim_snpStats,
chromosome = chr,
genetic.distance = posg,
position = pos,
allele.1 = ref,
allele.2 = alt
)
})
#> Writing FAM file to /tmp/RtmpMsCEXK/vignette-random-data-copy778e12960a6c.fam
#> Writing extended MAP file to /tmp/RtmpMsCEXK/vignette-random-data-copy778e12960a6c.bim
#> Writing BED file to /tmp/RtmpMsCEXK/vignette-random-data-copy778e12960a6c.bed (SNP-major mode)
# remove the new file, no longer need it
delete_files_plink(file_plink_copy)
The lfa
package also provides a minimal BED parser, ignoring the BIM and FAM files except to count loci and individuals. Its genotype matrix is exactly like genio
’s: same class (R matrix) and orientation (loci along rows). The package does not provide file writers or other such utilities.
library(lfa)
# Parse the genotype matrix
time_read_lfa <- system.time(
X_lfa <- read.bed(file_plink)
)
#> [1] "reading in 1001 individuals"
#> [1] "reading in 10001 snps"
#> [1] "snp major mode"
time_read_lfa
#> user system elapsed
#> 0.206 0.059 0.267
# This genotype matrix does not have column or row names:
X_lfa[1:10, 1:10]
#> [,1] [,2] [,3] [,4] [,5] [,6] [,7] [,8] [,9] [,10]
#> [1,] 2 1 0 2 2 NA 1 1 NA 1
#> [2,] NA 2 0 1 NA 0 1 1 1 2
#> [3,] 1 NA 2 1 1 0 2 NA 1 1
#> [4,] 1 2 1 1 0 2 2 2 0 NA
#> [5,] 0 1 1 NA 2 0 NA 2 1 0
#> [6,] 2 1 NA 0 1 1 0 0 0 2
#> [7,] 1 2 2 0 1 1 0 NA 2 2
#> [8,] 1 1 1 1 NA NA 2 2 0 1
#> [9,] 0 1 0 NA NA 0 2 1 1 1
#> [10,] 2 NA 1 1 1 1 NA 1 2 NA
# Again verify that the matrices are identical
stopifnot( all( X == X_lfa, na.rm = TRUE) )
Note that reading performance varies on different machines depending on the balance between hard drive access and processor speeds (where the relative bottleneck is). That being said, the genio
reader is consistently the fastest, if not a close second (see below for tests on your machine), for several reasons. Genotypes are read very efficiently using Rcpp
code, and stored directly into an R matrix, which is most efficient if that is the end goal. The annotation tables are also read most efficiently, using the readr
package internally.
The BEDMatrix
reader is also written in Rcpp
, but its main goal is to provide efficient random access to genotypes stored in a file, which makes obtaining an R
matrix more awkward, though surprisingly without paying a time penalty for it. The annotation tables are read with data.table
, which is actually the fastest, though the difference in performance is small compared to readr
for reasonably-sized files. However, BEDMatrix
does not process or return full annotation tables, which gives it an unfair advantage compared to genio
, as BEDMatrix
takes shortcuts to only read the columns it needs. If the annotation tables are needed, reading them will incur an additional (and in some cases considerable) time penalty.
The snpStats
reader is written in C
, so it is also very fast for its initial step, but converting the genotypes from the snpStats
format into a native R matrix proves too expensive. The annotation tables are read with the base function read.table
, which is also the slowest. In terms of completeness of output (full genotypes and annotations tables), only this package matches genio
, which makes it the fairest comparison.
The lfa
reader is written in pure R
, which explains why it is the slowest in this comparison. The annotation tables are read with read.table
, as in snpStats
, but are not returned at all.
# Extract component 3 of each time object,
# which is is total time elapsed.
# Sum the two steps it takes for each of BEDMatrix and snpStats to obtain a native R matrix.
times_read <- c(
time_read_genio[3],
time_read_bedmatrix_1[3] + time_read_bedmatrix_2[3],
time_read_snpStats_1[3] + time_read_snpStats_2[3],
time_read_lfa[3]
)
names_read <- c(
'genio',
'BEDMatrix',
'snpStats',
'lfa'
)
# Create barplot
barplot(
times_read,
names.arg = names_read,
main = 'BED reader runtimes',
xlab = 'packages',
ylab = 'runtime (s)'
)
Now we compare the writers. Only snpStats
and genio
have a BED writer. Not only was the genio
writer much easier to use, it is also considerably faster.
times_write <- c(
time_write_genio[3],
time_write_snpStats[3]
)
names_write <- c(
'genio',
'snpStats'
)
# Create barplot
barplot(
times_write,
names.arg = names_write,
main = 'BED writer runtimes',
xlab = 'packages',
ylab = 'runtime (s)'
)
High memory consumption is the main sacrifice for the ease of use of native R matrices, and in this respect genio
consumes the most memory (second only to lfa
).
We use the pryr
package to compare the memory usage of each native genotype object.
library(pryr)
# Store directly into a vector
sizes <- c(
object_size( X ),
object_size( data_genio$X ),
object_size( X_BEDMatrix ),
object_size( data_snpStats$genotypes ),
object_size( X_lfa )
)
names_sizes <- c(
'original',
'genio',
'BEDMatrix',
'snpStats',
'lfa'
)
# Create barplot
barplot(
sizes,
names.arg = names_sizes,
main = 'Native genotype object sizes',
xlab = 'packages',
ylab = 'memory (bytes)'
)
The BEDMatrix
package is ideal for very large files, since data is loaded into memory only as needed. So the main object does not actually hold any genotypes in memory. It is up to the user to decide how much data to access at once to trade-off speed and memory consumption.
snpStats
manages memory better than genio
, but strictly by a factor of 4. Each genotype is encoded as an integer in R
(in general) and genio
(in particular), which uses 4 bytes (this is actually platform dependent). On the other hand, snpStats
uses strictly one byte per genotype.
lfa
uses 8 bytes per genotype because its genotype matrix treats is elements as doubles rather than integers, so it uses twice as much memory as genio
.
Overall, there is a narrow window in data sizes in which a genotype matrix is too large for a native R matrix but small enough for a snpStats
object. That, and the much more complex interface of snpStats
, makes it not very worthwhile in my opinion, unless there is need for functions provided only by that package. I prefer to use BEDMatrix
for large datasets.
This handy function removes the three BED/BIM/FAM files we generated at the beginning of this vignette.
Let’s close by showing the package versions used when this vignette was last built, as the implementations compared here could change in the future.
sessionInfo()
#> R version 3.6.1 (2019-07-05)
#> Platform: x86_64-redhat-linux-gnu (64-bit)
#> Running under: Fedora 30 (Workstation Edition)
#>
#> Matrix products: default
#> BLAS/LAPACK: /usr/lib64/R/lib/libRblas.so
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_US.UTF-8 LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> attached base packages:
#> [1] stats graphics grDevices utils datasets methods base
#>
#> other attached packages:
#> [1] pryr_0.1.4 lfa_1.14.0 snpStats_1.34.0 Matrix_1.2-17
#> [5] survival_2.44-1.1 BEDMatrix_2.0.1 genio_1.0.12
#>
#> loaded via a namespace (and not attached):
#> [1] Rcpp_1.0.3 compiler_3.6.1 pillar_1.4.2
#> [4] tools_3.6.1 zlibbioc_1.30.0 crochet_2.2.0
#> [7] zeallot_0.1.0 digest_0.6.23 evaluate_0.14
#> [10] tibble_2.1.3 lattice_0.20-38 pkgconfig_2.0.3
#> [13] rlang_0.4.2 cli_2.0.0 yaml_2.2.0
#> [16] parallel_3.6.1 xfun_0.11 stringr_1.4.0
#> [19] knitr_1.26 hms_0.5.2 vctrs_0.2.1
#> [22] grid_3.6.1 glue_1.3.1 R6_2.4.1
#> [25] fansi_0.4.0 rmarkdown_2.0 readr_1.3.1
#> [28] corpcor_1.6.9 magrittr_1.5 backports_1.1.5
#> [31] codetools_0.2-16 htmltools_0.4.0 splines_3.6.1
#> [34] BiocGenerics_0.30.0 assertthat_0.2.1 utf8_1.1.4
#> [37] stringi_1.4.3 crayon_1.3.4