BSmeth2Probe: A function to map WGBS methylation data to Illumina Probe IDs

View source: R/BSmeth2Probe.R

BSmeth2ProbeR Documentation

A function to map WGBS methylation data to Illumina Probe IDs

Description

A function to map WGBS methylation data to Illumina Probe IDs

Usage

BSmeth2Probe(
  probe_id_locations,
  WGBS_data,
  cutoff = 10,
  multipleMapping = FALSE
)

Arguments

probe_id_locations

Either a dataframe or GRanges object containing probe IDs and their locations. If dataframe: must contain columns named "ID", "seqnames", "Start", "End", and "Strand". If GRanges: should have locations ("seqnames", "ranges", "strand"), as well as metadata column "ID". Start and end locations should be 1-based coordinates. Note that any row with NA values will not be used.

WGBS_data

Either a GRanges object or methylKit object (methylRaw, methylBase, methylRawDB, or methylBaseDB) of CpG locations and their methylation values. Contains locations ("seqnames", "ranges", "strand") and metadata column(s) of methylation values of sample(s) (i.e. one column per sample). These methylation values must be between 0 and 1.

cutoff

The maximum number of basepairs distance to consider for probes which have not been directly covered in the WGBS data. Default value is 10.

multipleMapping

When searching for matches for probes not directly covered in WGBS data, should WGBS CpGs which have already been mapped to another probe still be considered? If TRUE, then yes. If FALSE, then no. Default value is FALSE.

Value

A dataframe with first column "IDs" for CpG IDs, then 1 or more columns for methylation values of sample(s) (same number of samples as in WGBS_data) ID for each probe which was mapped, and then methylation value(s) of the WGBS CpG to which it was matched (where either it overlapped or the gap was < cutoff). If it matched to more than one CpG, the mean methylation value is taken.

Examples

data("IlluminaMethEpicB5ProbeIDs")
load(system.file("extdata", "WGBS_GRanges.rda", package = "deconvR"))
meth_probres <- BSmeth2Probe(
  probe_id_locations = IlluminaMethEpicB5ProbeIDs,
  WGBS_data = WGBS_GRanges
)
methp_cut <- BSmeth2Probe(
  probe_id_locations = IlluminaMethEpicB5ProbeIDs,
  WGBS_data = WGBS_GRanges[5:1000],
  cutoff = 2
)

BIMSBbioinfo/deconvR documentation built on April 26, 2023, 7:05 a.m.