subset_scMethrix: Subsets an 'scMethrix' object based on 'regions', 'contigs'...

Description Usage Arguments Details Value Examples

View source: R/scmethrix_operations.R

Description

Subsets an scMethrix object based on regions, contigs and/or samples.

Usage

1
2
3
4
5
6
7
8
9
subset_scMethrix(
  scm = NULL,
  regions = NULL,
  contigs = NULL,
  samples = NULL,
  by = c("include", "exclude"),
  overlap_type = c("within", "start", "end", "any", "equal"),
  verbose = TRUE
)

Arguments

scm

scMethrix; the single cell methylation experiment

regions

genomic regions to subset by. Could be a data.table with 3 columns (chr, start, end) or a GenomicRanges object

contigs

string; array of chromosome names to subset by

samples

string; array of sample names to subset by

by

string to decide whether to "include" or "exclude" the given criteria from the subset

overlap_type

string; defines the type of the overlap of the CpG sites with the target region. Default value is within. For detailed description, see the findOverlaps function of the IRanges package.

verbose

boolean; Flag for outputting function status messages. Default = TRUE

Details

Takes scMethrix object and filters CpGs based on region, contig and/or sample. Can either subset (include) to or filter (exclude) the specified parameters.

Value

An object of class scMethrix

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
data('scMethrix_data')

contigs <- c("chr1","chr3")
regions <- GenomicRanges::GRanges(seqnames = "chr1", ranges = IRanges(1,100000000)) 
samples <- c("C1","C2")

#Subset to only samples bed1 and bed3, and chromosome 1
subset_scMethrix(scMethrix_data, samples = samples, contigs = contigs, by = "include")

#Subset to only region "chr1:1-5"
subset_scMethrix(scMethrix_data, regions = regions, by = "include")

#Subset to exclude samples bed1 and bed3, and chromosome 1
subset_scMethrix(scMethrix_data, samples = samples, contigs = contigs, by = "exclude")

#Subset to exclude region "chr1:1-5"
subset_scMethrix(scMethrix_data, regions = regions, by = "exclude")

CompEpigen/scMethrix documentation built on Nov. 6, 2021, 3:09 p.m.