| import {methimpute} | R Documentation |
This page provides an overview of all methimpute data import functions.
importBSMAP(file, chrom.lengths = NULL, skip = 1, contexts = c(CG = "NNCGN", CHG = "NNCHG", CHH = "NNCHH")) importMethylpy(file, chrom.lengths = NULL, skip = 1, contexts = c(CG = "CGN", CHG = "CHG", CHH = "CHH")) importBSSeeker(file, chrom.lengths = NULL, skip = 0) importBismark(file, chrom.lengths = NULL, skip = 0)
file |
The file to import. |
chrom.lengths |
A data.frame with chromosome names in the first, and chromosome lengths in the second column. Only chromosomes named in here will be returned. Alternatively a tab-separated file with such a data.frame (with headers). |
skip |
The number of lines to skip. Usually 1 if the file contains a header and 0 otherwise. |
contexts |
A character vector of the contexts that are to be assigned. Since some programs report 5-letter contexts, this parameter can be used to obtain a reduced number of contexts. Will yield contexts CG, CHG, CHH by default. Set |
A methimputeData object.
importBSMAP: Import a BSMAP methylation extractor file.
importMethylpy: Import a Methylpy methylation extractor file.
importBSSeeker: Import a BSSeeker methylation extractor file.
importBismark: Import a Bismark methylation extractor file.
## Get an example file in BSSeeker format
file <- system.file("extdata","arabidopsis_bsseeker.txt.gz", package="methimpute")
data(arabidopsis_chromosomes)
bsseeker.data <- importBSSeeker(file, chrom.lengths=arabidopsis_chromosomes)
## Get an example file in Bismark format
file <- system.file("extdata","arabidopsis_bismark.txt", package="methimpute")
data(arabidopsis_chromosomes)
arabidopsis_chromosomes$chromosome <- sub('chr', '', arabidopsis_chromosomes$chromosome)
bismark.data <- importBismark(file, chrom.lengths=arabidopsis_chromosomes)
## Get an example file in BSMAP format
file <- system.file("extdata","arabidopsis_BSMAP.txt", package="methimpute")
data(arabidopsis_chromosomes)
bsmap.data <- importBSMAP(file, chrom.lengths=arabidopsis_chromosomes)
## Get an example file in Methylpy format
file <- system.file("extdata","arabidopsis_methylpy.txt", package="methimpute")
data(arabidopsis_chromosomes)
arabidopsis_chromosomes$chromosome <- sub('chr', '', arabidopsis_chromosomes$chromosome)
methylpy.data <- importMethylpy(file, chrom.lengths=arabidopsis_chromosomes)