View source: R/convert_and_query.R
convert_and_query | R Documentation |
If it is not tabix format already
(determined by checking for a .tbi
file of the same name in the same directory),
the full summary statistics file is converted into tabix format
for super fast querying.
A query is then made using the min/max genomic positions to extract a
locus-specific summary stats file.
convert_and_query(
target_path,
target_index = paste0(target_path, ".tbi"),
target_format = NULL,
study_dir = NULL,
target_chrom_col = "CHR",
target_start_col = "POS",
target_end_col = target_start_col,
query_granges,
samples = character(),
query_save = TRUE,
query_save_path = tempfile(fileext = ".gz"),
target_genome = "GRCh37",
query_genome = "GRCh37",
convert_methods = list(sort_coordinates = "bash", run_bgzip = "Rsamtools", index =
"Rsamtools"),
query_method = c("rsamtools", "seqminer", "conda"),
conda_env = "echoR_mini",
convert_force_new = FALSE,
query_force_new = FALSE,
nThread = 1,
verbose = TRUE
)
target_path |
Path to full GWAS/QTL summary statistics file. |
target_index |
Tabix index file for |
target_format |
Format of the |
study_dir |
Path to study folder. |
target_chrom_col |
Name of the chromosome column
in the |
target_start_col |
Name of the genomic start position column
in the |
target_end_col |
Name of the genomic end position column
in the |
query_granges |
GRanges object
to be used for querying the |
samples |
[Optional] Sample names to subset the VCF by. If this option is used, the GRanges object will be converted to a ScanVcfParam for usage by readVcf. |
query_save |
Whether to save the queried data subset. |
query_save_path |
Path to save retrieved query subset to. |
target_genome |
Genome build of the VCF file. |
query_genome |
Genome build that the |
convert_methods |
A named list containing methods to run each step with. |
query_method |
Method used for querying. See query for available options. |
conda_env |
Conda environments to search in.
If |
convert_force_new |
If the |
query_force_new |
If the query subset ( |
nThread |
Number of threads to use. |
verbose |
Print messages. |
data.table or VCF
of requested subset of target_path
.
Other tabix:
index_variantannotation()
query_dat <- echodata::BST1
target_path <- echodata::example_fullSS()
query_res <- echotabix::convert_and_query(
target_path = target_path,
target_start_col = "BP",
query_granges = query_dat,
query_force_new = TRUE)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.