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

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

Value

An object of class scMethrix

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.

Examples

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")
#> Subsetting CpG sites...
#> Subsetting by contigs
#> Subsetting by samples
#> Subset in 0.12s
#> An object of class scMethrix #> n_CpGs: 147 #> n_samples: 2 #> assays: score, counts #> reduced dims: #> is_h5: FALSE #> Reference: hg19 #> Physical size: 35 Kb
#Subset to only region "chr1:1-5" subset_scMethrix(scMethrix_data, regions = regions, by = "include")
#> Subsetting CpG sites...
#> Subsetting by regions
#> Subset in 0.07s
#> An object of class scMethrix #> n_CpGs: 67 #> n_samples: 4 #> assays: score, counts #> reduced dims: #> is_h5: FALSE #> Reference: hg19 #> Physical size: 34.2 Kb
#Subset to exclude samples bed1 and bed3, and chromosome 1 subset_scMethrix(scMethrix_data, samples = samples, contigs = contigs, by = "exclude")
#> Subsetting CpG sites...
#> Subsetting by contigs
#> Subsetting CpG sites...
#> Subsetting by contigs
#> Subset in 0.09s
#> Subsetting by samples
#> Subsetting CpG sites...
#> Subsetting by samples
#> Subset in 0.06s
#> Subset in [unknown time]
#> An object of class scMethrix #> n_CpGs: 139 #> n_samples: 2 #> assays: score, counts #> reduced dims: #> is_h5: FALSE #> Reference: hg19 #> Physical size: 34.8 Kb
#Subset to exclude region "chr1:1-5" subset_scMethrix(scMethrix_data, regions = regions, by = "exclude")
#> Subsetting CpG sites...
#> Subsetting by regions
#> Subsetting CpG sites...
#> Subsetting by regions
#> Subset in 0.08s
#> Subset in [unknown time]
#> An object of class scMethrix #> n_CpGs: 219 #> n_samples: 4 #> assays: score, counts #> reduced dims: #> is_h5: FALSE #> Reference: hg19 #> Physical size: 43.1 Kb