Description Usage Arguments Details Value Examples
View source: R/computeMatrix.R
Computes a metagene matrix from scratch.
1 2 3 4 5 6 7 8 9 10 11 12 | metagene_matrix(
bw_plus,
bw_minus,
anno,
anchor = "TSS",
upstream = 1000,
downstream = 1000,
collapse_to_win_size = NULL,
collapse_to_n_wins = NULL,
collapse_avg_fun = rowMeans,
negate_neg_strand_values = FALSE
)
|
bw_plus |
bigwig filename for plus strand signal |
bw_minus |
bigwig filename for minus strand signal |
anno |
Bed filename for annotations to use OR a GRanges object (see Details). |
anchor |
'TSS', 'TES', 'center' or NULL (default='TSS') |
upstream |
bp upstream of anchor (default=1000) |
downstream |
bp downstream of anchor (default=1000) |
collapse_to_win_size |
size for binning in bp (default=1, ie no binning) |
collapse_to_n_wins |
number of bins to collapse to (default=NULL, ie no binning) |
collapse_avg_fun |
function for averaging values inside windows (Default = rowMeans, see Details). |
negate_neg_strand_values |
negate values from minus strand bigwig (default=FALSE). |
Loads the bed file. If anchor is given, extracts the position of the anchor point and creates ranges from x bp upstream to x bp downstream. If anchor is NULL, collects signal in specified region from TSS-upstream to TES+downstream. Queries from both bigwigs and the correct strand. Summarizes signal of window_size bins and returns a tidy data.frame. negate_neg_strand_values can be set to TRUE to deal with minus strand bigwigs which contain all negative values as used by some labs. When using a GRanges object for anno, this assumes the bigwigs will be queries from start to end! All ranges must be exactly same size!
Tidy data frame with columns: gene (name, ie col4 from *bed* file), rel_pos (position in bp relative to *anchor*), value (averaged value from bigwig file).
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 | bedfile <- system.file("extdata", "Chen_PROMPT_TSSs_liftedTohg38.bed", package = "RMetaTools")
bw_plus <- system.file("extdata", "GSM1573841_mNET_8WG16_siLuc_plus_hg38.bw", package = "RMetaTools")
bw_minus <- system.file("extdata", "GSM1573841_mNET_8WG16_siLuc_minus_hg38.bw", package = "RMetaTools")
anno=bedfile;anchor='center';upstream=1000;downstream=1000;window_size=10
collapse_fun = rowMeans;negate_neg_strand_values=FALSE
metamat <- metagene_matrix(bw_plus, bw_minus, bedfile, anchor, upstream, downstream, window_size)
metamat <- metagene_matrix(bw_plus, bw_minus, bedfile, anchor, upstream, downstream, collapse_to_n_wins=5, collapse_avg_fun=mean)
bedfile <- system.file("extdata", "snRNA_MS31_hg38.bed", package = "RMetaTools")
metamat <- metagene_matrix(bw_plus, bw_minus, bedfile, anchor, upstream, downstream, window_size)
metamat <- metagene_matrix(bw_plus, bw_minus, bedfile, anchor=NULL, 10, 10, collapse_to_n_wins=5, collapse_avg_fun=mean)
# for pre-loaded regions (useful for ie modifying regions by hand)
# this most useful
regions <- meta_regions(bedfile, 'TSS', 1000, 5000)
start(regions) <- start(regions) + 1000
end(regions) <- end(regions) - 1000
regions <- regions[width(regions)>0,]
metamat <- metagene_matrix(bw_plus, bw_minus, regions, anchor=NULL, upstream=0, downstream=0, collapse_to_n_wins=5, collapse_avg_fun=mean)
metamat <- metagene_matrix(bw_plus, bw_minus, regions, anchor='TSS', upstream=100, downstream=100, collapse_to_win_size=50, collapse_avg_fun=rowMeans)
metamat <- metagene_matrix(bw_plus, bw_minus, bedfile, anchor='TSS', upstream=100, downstream=100, collapse_to_win_size=50, collapse_avg_fun=rowMeans)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.