scMethrix object based on regions, contigs and/or samples.R/scmethrix_operations.R
subset_scMethrix.RdSubsets 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 )
| scm |
|
|---|---|
| regions | genomic regions to subset by. Could be a data.table with 3 columns (chr, start, end) or a |
| 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 |
| verbose | boolean; Flag for outputting function status messages. Default = TRUE |
An object of class scMethrix
Takes scMethrix object and filters CpGs based on region, contig and/or sample. Can
either subset (include) to or filter (exclude) the specified parameters.
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")#>#>#>#>#> 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")#>#>#>#> 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")#>#>#>#>#>#>#>#>#>#>#> 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")#>#>#>#>#>#>#> 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