MungeSumstats: Getting startedIf you use the MungeSumstats package, please cite
The MungeSumstats package is designed to facilitate the standardisation of GWAS summary statistics as utilised in our Nature Genetics paper.1
The package is designed to handle the lack of standardisation of output files by the GWAS community. There is a group who have now manually standardised many GWAS: R interface to the IEU GWAS database API • ieugwasr and gwasvcf but because a lot of GWAS remain closed access, these repositories are not all encompassing.
The GWAS-Download project has collated summary statistics from 200+ GWAS. This repository has been utilsed to identify the most common formats, all of which can be standardised with MungeSumstats.
Moreover, there is an emerging standard of VCF format for summary statistics files with multiple, useful, associated R packages such as vcfR. However, there is currently no method to convert VCF formats to a standardised format that matches older approaches.
The MungeSumstats package standardises both VCF and the most common summary statistic file formats to enable downstream integration and analysis.
MungeSumstats also offers comprehensive Quality Control (QC) steps which are important prerequisites for downstream analysis like Linkage disequilibrium score regression (LDSC) and MAGMA.
Moreover, MungeSumstats is efficiently written resulting in all
reformatting and quality control checks completing in minutes for GWAS
summary statistics with 500k SNPs on a standard desktop machine. This
speed can be increased further by increasing the number of threads
(nThread) for data.table to use.
Currently MungeSumstats only works on data from humans, as it uses human-based genome references.
MungeSumstats will ensure that the all essential columns for analysis are present and syntactically correct. Generally, summary statistic files include (but are not limited to) the columns:
MungeSumstats uses a mapping file to infer the inputted column names
(run data("sumstatsColHeaders") to view these). This mapping file is
far more comprehensive than any other publicly available munging tool
containing more than 200 unique mappings at the time of writing this
vignette. However, if your column headers are missing or if you want to
change the mapping, you can do so by passing your own mapping file (see
format_sumstats(mapping_file)).
MungeSumstats offers unmatched levels of quality control to ensure, for example, consistency of allele assignment and direction of effects. Tests run by MungeSumstats include:
Users can specify which checks to run on their data. A note on the allele flipping check: MungeSumstats infers the effect allele will always be the A2 allele, this is the approach done for IEU GWAS VCF and has such also been adopted here. This inference is first from the inputted file’s column headers however, the allele flipping check ensures this by comparing A1, what should be the reference allele, to the reference genome. If a SNP’s A1 DNA base doesn’t match the reference genome but it’s A2 (what should be the alternative allele) does, the alleles will be flipped along with the effect information (e.g. Beta, Odds Ratio, signed summary statistics, FRQ, Z-score*).
*-by default the Z-score is assumed to be calculated off the effect size not the P-value and so will be flipped. This can be changed by a user.
If a test is failed, the user will be notified and if possible, the
input will be corrected. The QC steps from the checks above can also be
adjusted to suit the user’s analysis, see
MungeSumstats::format_sumstats.
MungeSumstats can handle VCF, txt, tsv, csv file types or .gz/.bgz versions of these file types. The package also gives the user the flexibility to export the reformatted file as tab-delimited, VCF or R native objects such as data.table, GRanges or VRanges objects. The output can also be outputted in an LDSC ready format which means the file can be fed directly into LDSC without the need for additional munging. NOTE - If LDSC format is used, the naming convention of A1 as the reference (genome build) allele and A2 as the effect allele will be reversed to match LDSC (A1 will now be the effect allele). See more info on this here. Note that any effect columns (e.g. Z) will be inrelation to A1 now instead of A2.
Please read carefully through our FAQ Website to gain insight on how best to run MungeSumstats on your data.
The MungeSumstats package contains small subsets of GWAS summary statistics files. Firstly, on Educational Attainment by Okbay et al 2016: PMID: 27898078 PMCID: PMC5509058 DOI: 10.1038/ng1216-1587b.
Secondly, a VCF file (VCFv4.2) relating to the GWAS Amyotrophic lateral sclerosis from ieu open GWAS project. Dataset: ebi-a-GCST005647: https://gwas.mrcieu.ac.uk/datasets/ebi-a-GCST005647/
These datasets will be used to showcase MungeSumstats functionality.
MungeSumstats is available on Bioconductor. To install the package on Bioconductor run the following lines of code:
if (!require("BiocManager")) install.packages("BiocManager")
BiocManager::install("MungeSumstats")Once installed, load the package:
library(MungeSumstats)To standardise the summary statistics’ file format, simply call
format_sumstats() passing in the path to your summary statistics file
or directly pass the summary statistics as a dataframe or datatable. You
can specify which genome build was used in the GWAS(GRCh37 or GRCh38)
or, as default, infer the genome build from the data.The reference
genome is used for multiple checks like deriving missing data such
SNP/BP/CHR/A1/A2 and for QC steps like removing non-biallelic SNPs,
strand-ambiguous SNPs or ensuring correct allele and direction of SNP
effects. The path to the reformatted summary statistics file can be
returned by the function call, the user can specify a location to save
the file or the user can return an R native object for the data:
data.table, VRanges or GRanges object.
Note that for a number of the checks implored by MungeSumstats a
reference genome is used. If your GWAS summary statistics file of
interest relates to GRCh38, you will need to install
SNPlocs.Hsapiens.dbSNP155.GRCh38 and BSgenome.Hsapiens.NCBI.GRCh38
from Bioconductor as follows:
#increase permissible time to download data, in case of slow internet access
options(timeout=2000)
BiocManager::install("SNPlocs.Hsapiens.dbSNP155.GRCh38")
BiocManager::install("BSgenome.Hsapiens.NCBI.GRCh38")If your GWAS summary statistics file of interest relates to GRCh37,
you will need to install SNPlocs.Hsapiens.dbSNP155.GRCh37 and
BSgenome.Hsapiens.1000genomes.hs37d5 from Bioconductor as follows:
BiocManager::install("SNPlocs.Hsapiens.dbSNP155.GRCh37")
BiocManager::install("BSgenome.Hsapiens.1000genomes.hs37d5")These may take some time to install and are not included in the package as some users may only need one of GRCh37/GRCh38.
The Educational Attainment by Okbay GWAS summary statistics file is saved as a text document in the package’s external data folder so we can just pass the file path to it straight to MungeSumstats.
NOTE - By default, Formatted results will be saved to tempdir().
This means all formatted summary stats will be deleted upon ending the R
session if not copied to a local file path. Otherwise, to keep formatted
summary stats, change save_path (
e.g.file.path('./formatted',basename(path))), or make sure to copy
files elsewhere after processing (
e.g.file.copy(save_path, './formatted/' ).
eduAttainOkbayPth <- system.file("extdata","eduAttainOkbay.txt",
                                  package="MungeSumstats")
reformatted <- 
  MungeSumstats::format_sumstats(path=eduAttainOkbayPth,
                                 ref_genome="GRCh37")## 
## 
## ******::NOTE::******
##  - Formatted results will be saved to `tempdir()` by default.
##  - This means all formatted summary stats will be deleted upon ending the R session.
##  - To keep formatted summary stats, change `save_path`  ( e.g. `save_path=file.path('./formatted',basename(path))` ),   or make sure to copy files elsewhere after processing  ( e.g. `file.copy(save_path, './formatted/' )`.
##  ********************## Formatted summary statistics will be saved to ==>  /tmp/RtmpCpT9Tr/file221b343a459308.tsv.gz## Warning: replacing previous import 'utils::findMatches' by
## 'S4Vectors::findMatches' when loading 'SNPlocs.Hsapiens.dbSNP155.GRCh37'## Importing tabular file: /tmp/RtmpsSO0dx/Rinst220cdc78d9feb1/MungeSumstats/extdata/eduAttainOkbay.txt## Checking for empty columns.## Infer Effect Column## First line of summary statistics file:## MarkerName   CHR POS A1  A2  EAF Beta    SE  Pval    ## Allele columns are ambiguous, attempting to infer direction## Can't infer allele columns from sumstats## Standardising column headers.## First line of summary statistics file:## MarkerName   CHR POS A1  A2  EAF Beta    SE  Pval    ## Summary statistics report:
##    - 93 rows
##    - 93 unique variants
##    - 70 genome-wide significant variants (P<5e-8)
##    - 20 chromosomes## Checking for multi-GWAS.## Checking for multiple RSIDs on one row.## Checking SNP RSIDs.## Checking for merged allele column.## Checking A1 is uppercase## Checking A2 is uppercase## Checking for incorrect base-pair positions## Checking for missing data.## Checking for duplicate columns.## Checking for duplicated rows.## INFO column not available. Skipping INFO score filtering step.## Filtering SNPs, ensuring SE>0.## Ensuring all SNPs have N<5 std dev above mean.## 47 SNPs (50.5%) have FRQ values > 0.5. Conventionally the FRQ column is intended to show the minor/effect allele frequency.
## The FRQ column was mapped from one of the following from the inputted  summary statistics file:
## FRQ, EAF, FREQUENCY, FRQ_U, F_U, MAF, FREQ, FREQ_TESTED_ALLELE, FRQ_TESTED_ALLELE, FREQ_EFFECT_ALLELE, FRQ_EFFECT_ALLELE, EFFECT_ALLELE_FREQUENCY, EFFECT_ALLELE_FREQ, EFFECT_ALLELE_FRQ, A2FREQ, A2FRQ, ALLELE_FREQUENCY, ALLELE_FREQ, ALLELE_FRQ, AF, MINOR_AF, EFFECT_AF, A2_AF, EFF_AF, ALT_AF, ALTERNATIVE_AF, INC_AF, A_2_AF, TESTED_AF, ALLELEFREQ, ALT_FREQ, EAF_HRC, EFFECTALLELEFREQ, FREQ.B, FREQ_EUROPEAN_1000GENOMES, FREQ_HAPMAP, FREQ_TESTED_ALLELE_IN_HRS, FRQ_U_113154, FRQ_U_31358, FRQ_U_344901, FRQ_U_43456, POOLED_ALT_AF, AF_ALT, AF.ALT, AF-ALT, ALT.AF, ALT-AF, A2.AF, A2-AF, AF.EFF, AF_EFF, ALL_AF## As frq_is_maf=TRUE, the FRQ column will not be renamed. If the FRQ values were intended to represent major allele frequency,
## set frq_is_maf=FALSE to rename the column as MAJOR_ALLELE_FRQ and differentiate it from minor/effect allele frequency.## Sorting coordinates with 'data.table'.## Writing in tabular format ==> /tmp/RtmpCpT9Tr/file221b343a459308.tsv.gz## Summary statistics report:
##    - 93 rows (100% of original 93 rows)
##    - 93 unique variants
##    - 70 genome-wide significant variants (P<5e-8)
##    - 20 chromosomes## Done munging in 0.072 minutes.## Successfully finished preparing sumstats file, preview:## Reading header.##           SNP   CHR       BP     A1     A2     FRQ   BETA    SE         P
##        <char> <int>    <int> <char> <char>   <num>  <num> <num>     <num>
## 1:   rs301800     1  8490603      T      C 0.17910  0.019 0.003 1.794e-08
## 2: rs11210860     1 43982527      A      G 0.36940  0.017 0.003 2.359e-10
## 3: rs34305371     1 72733610      A      G 0.08769  0.035 0.005 3.762e-14
## 4:  rs2568955     1 72762169      T      C 0.23690 -0.017 0.003 1.797e-08## Returning path to saved data.Here we know the summary statistics are based on the reference genome GRCh37,
GRCh38 can also be inputted. Moreover, if you are unsure of the genome build,
leave it as NULL and Mungesumstats will infer it from the data.
Also note that the default dbSNP version used along with the reference genome is
the latest version available on Bioconductor (currently dbSNP 155) but older
versions are also availble. Use the dbSNP input parameter to control this.
The arguments format_sumstats in that control the level of QC
conducted by MungeSumstats are:
p-values < 5e-324 be converted
to 0? Small p-values pass the R limit and can cause errors with
LDSC/MAGMA and should be converted. Default is TRUE.1-22, X, Y, MT; the UCSC
style is chr1-chr22, chrX, chrY, chrM; and the dbSNP style is
ch1-ch22, chX, chY, chMT. Default is Ensembl.c("X", "Y", "MT") which removes all non-autosomal SNPs.NULL, all columns will be checked for
missing values. Default columns are SNP, chromosome, position, allele 1,
allele 2, effect columns (frequency, beta, Z-score, standard error,
log odds, signed sumstats, odds ratio), p value and N columns.data.table, GRanges or VRangesdirectly
to user. Otherwise, return the path to the save data. Default is
FALSE.data(sumstatsColHeaders) for
default mapping and necessary format.See ?MungeSumstats::format_sumstats() for the full list of parameters
to control MungeSumstats QC and standardisation steps.
VCF files can also be standardised to the same format as other summary statistic files. A subset of the Amyotrophic lateral sclerosis GWAS from the ieu open GWAS project (a .vcf file) has been added to MungeSumstats to demonstrate this functionality.Simply pass the path to the file in the same manner you would for other summary statistic files:
#save ALS GWAS from the ieu open GWAS project to a temp directory
ALSvcfPth <- system.file("extdata","ALSvcf.vcf", package="MungeSumstats")reformatted_vcf <- 
  MungeSumstats::format_sumstats(path=ALSvcfPth, 
                                 ref_genome="GRCh37")You can also get more information on the SNPs which have had data
imputed or have been filtered out by MungeSumstats by using the
imputation_ind and log_folder_ind parameters respectively. For
example:
#set
reformatted_vcf_2 <- 
  MungeSumstats::format_sumstats(path=ALSvcfPth,
                                 ref_genome="GRCh37",
                                 log_folder_ind=TRUE,
                                 imputation_ind=TRUE,
                                 log_mungesumstats_msgs=TRUE)## Time difference of 0.1 secs## Time difference of 0.5 secsCheck the file snp_bi_allelic.tsv.gz in the log_folder directory you supply
(by default a temp directory), for a list of SNPs removed as they are
non-bi-allelic. The text files containing the console output and messages are
also stored in the same directory.
Note you can also control the dbSNP version used as a reference dataset by
MungeSumstats using the dbSNP parameter. By default this will be set to the
most recent dbSNP version available (155).
Note that using log_folder_ind returns a list from format_sumstats
which includes the file locations of the differing classes of removed
SNPs. Using log_mungesumstats_msgs saves the messages to the console
to a file which is returned in the same list. Note that not all the
messages will also print to screen anymore when you set
log_mungesumstats_msgs:
names(reformatted_vcf_2)## [1] "sumstats"  "log_files"A user can load a file to view the excluded SNPs.
In this case, SNPs were filtered based on non-bi-allelic criterion:
print(reformatted_vcf_2$log_files$snp_bi_allelic)## NULLThe different types of exclusion which lead to the names are explained below:
Note to export to another type such as R native objects including
data.table, GRanges, VRanges or save as a VCF file, set return_data=TRUE and
choose your return_format:
#set
reformatted_vcf_2 <- 
  MungeSumstats::format_sumstats(path=ALSvcfPth,
                                 ref_genome="GRCh37", 
                                 log_folder_ind=TRUE,
                                 imputation_ind=TRUE,
                                 log_mungesumstats_msgs=TRUE,
                                 return_data=TRUE,
                                 return_format="GRanges")Also you can now output a VCF compatible with IEU OpenGWAS format (Note that currently all IEU OpenGWAS sumstats are GRCh37, MungeSumstats will throw a warning if your data isn’t GRCh37 when saving):
#set
reformatted_vcf_2 <- 
  MungeSumstats::format_sumstats(path=ALSvcfPth,
                                 ref_genome="GRCh37", 
                                 write_vcf=TRUE,
                                 save_format ="openGWAS")See our publication for further discussion of these checks and options:
MungeSumstats also contains a function to quickly infer the genome
build of multiple summary statistic files. This can be called separately
to format_sumstats() which is useful if you want to just quickly check
the genome build:
# Pass path to Educational Attainment Okbay sumstat file to a temp directory
eduAttainOkbayPth <- system.file("extdata", "eduAttainOkbay.txt",
                                  package = "MungeSumstats")
ALSvcfPth <- system.file("extdata","ALSvcf.vcf", package="MungeSumstats")
sumstats_list <- list(ss1 = eduAttainOkbayPth, ss2 = ALSvcfPth)
ref_genomes <- MungeSumstats::get_genome_builds(sumstats_list = sumstats_list)MungeSumstats exposes the liftover() function as a general utility
for users.
Useful features include: - Automatic standardisation of genome build
names (i.e. “hg19”, “hg37”, and “GRCh37” will all be recognized as the
same genome build.) - Ability to specify chrom_col as well as both
start_col and end_col (for variants that span >1bp). - Ability to
return in data.table or GRanges format. - Ability to specify which
chromosome format (e.g. “chr1” vs. 1) to return GRanges as.
sumstats_dt <- MungeSumstats::formatted_example()## Standardising column headers.## First line of summary statistics file:## MarkerName   CHR POS A1  A2  EAF Beta    SE  Pval    ## Sorting coordinates with 'data.table'.sumstats_dt_hg38 <- MungeSumstats::liftover(sumstats_dt = sumstats_dt, 
                                            ref_genome = "hg19",
                                            convert_ref_genome = "hg38")## Performing data liftover from hg19 to hg38.## Converting summary statistics to GenomicRanges.## Downloading chain file...## Downloading chain file from Ensembl.## /tmp/RtmpCpT9Tr/GRCh37_to_GRCh38.chain.gz## Reordering so first three column headers are SNP, CHR and BP in this order.## Reordering so the fourth and fifth columns are A1 and A2.knitr::kable(head(sumstats_dt_hg38))| SNP | CHR | BP | A1 | A2 | FRQ | BETA | SE | P | IMPUTATION_gen_build | 
|---|---|---|---|---|---|---|---|---|---|
| rs301800 | 1 | 8430543 | T | C | 0.17910 | 0.019 | 0.003 | 0e+00 | TRUE | 
| rs11210860 | 1 | 43516856 | A | G | 0.36940 | 0.017 | 0.003 | 0e+00 | TRUE | 
| rs34305371 | 1 | 72267927 | A | G | 0.08769 | 0.035 | 0.005 | 0e+00 | TRUE | 
| rs2568955 | 1 | 72296486 | T | C | 0.23690 | -0.017 | 0.003 | 0e+00 | TRUE | 
| rs1008078 | 1 | 90724174 | T | C | 0.37310 | -0.016 | 0.003 | 0e+00 | TRUE | 
| rs61787263 | 1 | 98153158 | T | C | 0.76120 | 0.016 | 0.003 | 1e-07 | TRUE | 
In some cases, users may not want to run the full munging pipeline
provided by
MungeSumstats::format_sumstats, but still would like to take advantage
of the file type conversion and column header standardisation features.
This will not be nearly as robust as the full pipeline, but can still be
helpful.
To do this, simply run the following:
eduAttainOkbayPth <- system.file("extdata", "eduAttainOkbay.txt",
                                  package = "MungeSumstats")
formatted_path <- tempfile(fileext = "_eduAttainOkbay_standardised.tsv.gz")
#### 1. Read in the data and standardise header names ####
dat <- MungeSumstats::read_sumstats(path = eduAttainOkbayPth, 
                                    standardise_headers = TRUE)## Importing tabular file: /tmp/RtmpsSO0dx/Rinst220cdc78d9feb1/MungeSumstats/extdata/eduAttainOkbay.txt## Checking for empty columns.## Standardising column headers.## First line of summary statistics file:## MarkerName   CHR POS A1  A2  EAF Beta    SE  Pval    knitr::kable(head(dat))| SNP | CHR | BP | A1 | A2 | FRQ | BETA | SE | P | 
|---|---|---|---|---|---|---|---|---|
| rs10061788 | 5 | 87934707 | A | G | 0.2164 | 0.021 | 0.004 | 0e+00 | 
| rs1007883 | 16 | 51163406 | T | C | 0.3713 | -0.015 | 0.003 | 1e-07 | 
| rs1008078 | 1 | 91189731 | T | C | 0.3731 | -0.016 | 0.003 | 0e+00 | 
| rs1043209 | 14 | 23373986 | A | G | 0.6026 | 0.018 | 0.003 | 0e+00 | 
| rs10496091 | 2 | 61482261 | A | G | 0.2705 | -0.018 | 0.003 | 0e+00 | 
| rs10930008 | 2 | 161854736 | A | G | 0.7183 | -0.016 | 0.003 | 1e-07 | 
#### 2. Write to disk as a compressed, tab-delimited, tabix-indexed file ####
formatted_path <- MungeSumstats::write_sumstats(sumstats_dt = dat,
                                                save_path = formatted_path,
                                                tabix_index = TRUE,
                                                write_vcf = FALSE,
                                                return_path = TRUE)   ## Sorting coordinates with 'data.table'.## Writing in tabular format ==> /tmp/RtmpCpT9Tr/file221b344212def2_eduAttainOkbay_standardised.tsv## Writing uncompressed instead of gzipped to enable tabix indexing.## Converting full summary stats file to tabix format for fast querying...## Reading header.## Ensuring file is bgzipped.## Tabix-indexing file.## Removing temporary .tsv file.data.tableIf you already have your data imported as an data.table, you can also
standardise its headers like so:
#### Mess up some column names ####
dat_raw <- data.table::copy(dat)
data.table::setnames(dat_raw, c("SNP","CHR"), c("rsID","Seqnames"))
#### Add a non-standard column that I want to keep the casing for ####
dat_raw$Support <- runif(nrow(dat_raw))
dat2 <- MungeSumstats::standardise_header(sumstats_dt = dat_raw,
                                          uppercase_unmapped = FALSE, 
                                          return_list = FALSE )## Standardising column headers.## First line of summary statistics file:## rsID Seqnames    BP  A1  A2  FRQ BETA    SE  P   Support ## Returning unmapped column names without making them uppercase.knitr::kable(head(dat2))| SNP | CHR | BP | A1 | A2 | FRQ | BETA | SE | P | Support | 
|---|---|---|---|---|---|---|---|---|---|
| rs301800 | 1 | 8490603 | T | C | 0.17910 | 0.019 | 0.003 | 0e+00 | 0.8291361 | 
| rs11210860 | 1 | 43982527 | A | G | 0.36940 | 0.017 | 0.003 | 0e+00 | 0.4046099 | 
| rs34305371 | 1 | 72733610 | A | G | 0.08769 | 0.035 | 0.005 | 0e+00 | 0.7170452 | 
| rs2568955 | 1 | 72762169 | T | C | 0.23690 | -0.017 | 0.003 | 0e+00 | 0.6517453 | 
| rs1008078 | 1 | 91189731 | T | C | 0.37310 | -0.016 | 0.003 | 0e+00 | 0.4898946 | 
| rs61787263 | 1 | 98618714 | T | C | 0.76120 | 0.016 | 0.003 | 1e-07 | 0.8882005 | 
The MungeSumstats package aims to be able to handle the most common summary statistic file formats including VCF. If your file can not be formatted by MungeSumstats feel free to report the bug on github: https://github.com/neurogenomics/MungeSumstats along with your summary statistic file header.
We also encourage people to edit the code to resolve their particular issues too and are happy to incorporate these through pull requests on github. If your summary statistic file headers are not recognised by MungeSumstats but correspond to one of:
SNP, BP, CHR, A1, A2, P, Z, OR, BETA, LOG_ODDS,
SIGNED_SUMSTAT, N, N_CAS, N_CON, NSTUDY, INFO or FRQ feel free to update the MungeSumstats::sumstatsColHeaders following
the approach in the data.R file and add your mapping. Then use a pull
request on github and we will incorporate this change into the package.
A note on MungeSumstats::sumstatsColHeaders for summary statistic
files with A0/A1. The mapping in MungeSumstats::sumstatsColHeaders
converts A0 to A*, this is a special case so that the code knows to map
A0/A1 to A1/A2 (ref/alt). The special case is needed since ordinarily A1
refers to the reference not the alternative allele.
A note on MungeSumstats::sumstatsColHeaders for summary statistic
files with Effect Size (ES). By default, MSS takes BETA to be any BETA-like
value (including ES). This is coded into the mapping file -
MungeSumstats::sumstatsColHeaders. If this isn’t the case for your sumstats,
you can set the es_is_beta parameter in MungeSumstats::format_sumstats() to
FALSE to avoid this. Note this is done to try and capture most use cases of MSS.
See the Open GWAS vignette for how MungeSumstats can be used along with data from the MRC IEU Open GWAS Project and also Mungesumstats’ functionality to handle lists of summary statistics files.
## R version 4.5.0 RC (2025-04-04 r88126)
## Platform: x86_64-pc-linux-gnu
## Running under: Ubuntu 24.04.2 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so 
## LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0  LAPACK version 3.12.0
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              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       
## 
## time zone: America/New_York
## tzcode source: system (glibc)
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] MungeSumstats_1.16.0 BiocStyle_2.36.0    
## 
## loaded via a namespace (and not attached):
##  [1] tidyselect_1.2.1                           
##  [2] dplyr_1.1.4                                
##  [3] blob_1.2.4                                 
##  [4] R.utils_2.13.0                             
##  [5] Biostrings_2.76.0                          
##  [6] bitops_1.0-9                               
##  [7] fastmap_1.2.0                              
##  [8] RCurl_1.98-1.17                            
##  [9] VariantAnnotation_1.54.0                   
## [10] GenomicAlignments_1.44.0                   
## [11] XML_3.99-0.18                              
## [12] digest_0.6.37                              
## [13] lifecycle_1.0.4                            
## [14] KEGGREST_1.48.0                            
## [15] RSQLite_2.3.9                              
## [16] magrittr_2.0.3                             
## [17] compiler_4.5.0                             
## [18] rlang_1.1.6                                
## [19] sass_0.4.10                                
## [20] tools_4.5.0                                
## [21] yaml_2.3.10                                
## [22] data.table_1.17.0                          
## [23] rtracklayer_1.68.0                         
## [24] knitr_1.50                                 
## [25] S4Arrays_1.8.0                             
## [26] bit_4.6.0                                  
## [27] curl_6.2.2                                 
## [28] DelayedArray_0.34.0                        
## [29] ieugwasr_1.0.3                             
## [30] abind_1.4-8                                
## [31] BiocParallel_1.42.0                        
## [32] BiocGenerics_0.54.0                        
## [33] R.oo_1.27.0                                
## [34] grid_4.5.0                                 
## [35] stats4_4.5.0                               
## [36] SummarizedExperiment_1.38.0                
## [37] cli_3.6.4                                  
## [38] rmarkdown_2.29                             
## [39] crayon_1.5.3                               
## [40] generics_0.1.3                             
## [41] BSgenome.Hsapiens.1000genomes.hs37d5_0.99.1
## [42] httr_1.4.7                                 
## [43] rjson_0.2.23                               
## [44] DBI_1.2.3                                  
## [45] cachem_1.1.0                               
## [46] stringr_1.5.1                              
## [47] parallel_4.5.0                             
## [48] AnnotationDbi_1.70.0                       
## [49] BiocManager_1.30.25                        
## [50] XVector_0.48.0                             
## [51] restfulr_0.0.15                            
## [52] matrixStats_1.5.0                          
## [53] vctrs_0.6.5                                
## [54] Matrix_1.7-3                               
## [55] jsonlite_2.0.0                             
## [56] bookdown_0.43                              
## [57] IRanges_2.42.0                             
## [58] S4Vectors_0.46.0                           
## [59] bit64_4.6.0-1                              
## [60] GenomicFiles_1.44.0                        
## [61] GenomicFeatures_1.60.0                     
## [62] jquerylib_0.1.4                            
## [63] glue_1.8.0                                 
## [64] codetools_0.2-20                           
## [65] stringi_1.8.7                              
## [66] GenomeInfoDb_1.44.0                        
## [67] BiocIO_1.18.0                              
## [68] GenomicRanges_1.60.0                       
## [69] UCSC.utils_1.4.0                           
## [70] tibble_3.2.1                               
## [71] pillar_1.10.2                              
## [72] SNPlocs.Hsapiens.dbSNP155.GRCh37_0.99.24   
## [73] htmltools_0.5.8.1                          
## [74] GenomeInfoDbData_1.2.14                    
## [75] BSgenome_1.76.0                            
## [76] R6_2.6.1                                   
## [77] evaluate_1.0.3                             
## [78] lattice_0.22-7                             
## [79] Biobase_2.68.0                             
## [80] R.methodsS3_1.8.2                          
## [81] png_0.1-8                                  
## [82] Rsamtools_2.24.0                           
## [83] memoise_2.0.1                              
## [84] bslib_0.9.0                                
## [85] SparseArray_1.8.0                          
## [86] xfun_0.52                                  
## [87] MatrixGenerics_1.20.0                      
## [88] pkgconfig_2.0.31. Nathan G. Skene, T. E. B., Julien Bryois. Genetic identification of brain cell types underlying schizophrenia. Nature Genetics (2018). doi:10.1038/s41588-018-0129-5