Extract and summarize methylation or coverage info by regions of interest

get_region_summary(
  m,
  regions = NULL,
  type = "M",
  how = "mean",
  overlap_type = "within",
  na_rm = TRUE,
  elementMetadata.col = NULL,
  verbose = TRUE,
  n_chunks = 1,
  n_cores = 1
)

Arguments

m

methrix object

regions

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

type

matrix which needs to be summarized. Coule be `M`, `C`. Default 'M'

how

mathematical function by which regions should be summarized. Can be one of the following: mean, sum, max, min. Default 'mean'

overlap_type

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.

na_rm

Remove NA's? Default TRUE

elementMetadata.col

columns in methrix@elementMetadata which needs to be summarised. Default = NULL.

verbose

Default TRUE

n_chunks

Number of chunks to split the methrix object in case it is very large. Default = 1.

n_cores

Number of parallel instances. n_cores should be less than or equal to n_chunks. If n_chunks is not specified, then n_chunks is initialized to be equal to n_cores. Default = 1.

Value

a coverage or methylation matrix

Details

Takes methrix object and summarizes regions

Examples

data('methrix_data') get_region_summary(m = methrix_data, regions = data.table(chr = 'chr21', start = 27867971, end = 27868103), type = 'M', how = 'mean')
#> -Checking for overlaps..
#> -Summarizing by average
#> -Done! Finished in:0.100s elapsed (0.120s cpu)
#> chr start end n_overlap_CpGs rid C1 C2 N1 #> 1: chr21 27867971 27868103 4 1 0.4055556 0.2378247 0.95 #> N2 #> 1: 0.9058442