Description Usage Arguments Details Value Author(s) Examples
Normalize Associations between Genomic Signals and Target Regions into a Matrix
1 2 3 4 5 6 | normalizeToMatrix(signal, target, extend = 5000, w = max(extend)/50,
value_column = NULL, mapping_column = NULL, background = ifelse(smooth, NA, 0), empty_value = NULL,
mean_mode = c("absolute", "weighted", "w0", "coverage"), include_target = any(width(target) > 1),
target_ratio = min(c(0.4, mean(width(target))/(sum(extend) + mean(width(target))))),
k = min(c(20, min(width(target)))), smooth = FALSE, smooth_fun = default_smooth_fun,
keep = c(0, 1), trim = NULL, flip_upstream = FALSE)
|
signal |
A |
target |
A |
extend |
Extended base pairs to the upstream and/or downstream of |
w |
Window size for splitting upstream and downstream, measured in base pairs |
value_column |
Column index in |
mapping_column |
Mapping column to restrict overlapping between |
background |
Values for windows that don't overlap with |
empty_value |
Deprecated, please use |
mean_mode |
When a window is not perfectly overlapped to |
include_target |
Whether include |
target_ratio |
The ratio of |
k |
Number of windows only when |
smooth |
Whether apply smoothing on rows in the matrix? |
smooth_fun |
The smoothing function that is applied to each row in the matrix. This self-defined function accepts a numeric vector (may contain |
keep |
Percentiles in the normalized matrix to keep. The value is a vector of two percent values. Values less than the first percentile is replaces with the first pencentile and values larger than the second percentile is replaced with the second percentile. |
trim |
Deprecated, please use |
flip_upstream |
Sometimes whether the signals are on the upstream or the downstream of the targets are not important and users only want to show the relative distance to targets. If the value is set to |
In order to visualize associations between signal
and target
, the data is transformed into a matrix
and visualized as a heatmap by EnrichedHeatmap
afterwards.
Upstream and downstream also with the target body are splitted into a list of small windows and overlap
to signal
. Since regions in signal
and small windows do not always 100 percent overlap, there are four different averaging modes:
Following illustrates different settings for mean_mode
(note there is one signal region overlapping with other signals):
1 2 3 4 5 6 7 8 9 10 11 | 40 50 20 values in signal regions
++++++ +++ +++++ signal regions
30 values in signal region
++++++ signal region
================= a window (17bp), there are 4bp not overlapping to any signal regions.
4 6 3 3 overlap
absolute: (40 + 30 + 50 + 20)/4
weighted: (40*4 + 30*6 + 50*3 + 20*3)/(4 + 6 + 3 + 3)
w0: (40*4 + 30*6 + 50*3 + 20*3)/(4 + 6 + 3 + 3 + 4)
coverage: (40*4 + 30*6 + 50*3 + 20*3)/17
|
A matrix with following additional attributes:
upstream_index
column index corresponding to upstream of target
target_index
column index corresponding to target
downstream_index
column index corresponding to downstream of target
extend
extension on upstream and downstream
smooth
whether smoothing was applied on the matrix
failed_rows
index of rows which are failed after smoothing
The matrix is wrapped into a simple normalizeToMatrix
class.
Zuguang Gu <z.gu@dkfz.de>
1 2 3 4 5 6 7 8 | signal = GRanges(seqnames = "chr1",
ranges = IRanges(start = c(1, 4, 7, 11, 14, 17, 21, 24, 27),
end = c(2, 5, 8, 12, 15, 18, 22, 25, 28)),
score = c(1, 2, 3, 1, 2, 3, 1, 2, 3))
target = GRanges(seqnames = "chr1", ranges = IRanges(start = 10, end = 20))
normalizeToMatrix(signal, target, extend = 10, w = 2)
normalizeToMatrix(signal, target, extend = 10, w = 2, include_target = TRUE)
normalizeToMatrix(signal, target, extend = 10, w = 2, value_column = "score")
|
Loading required package: grid
Loading required package: ComplexHeatmap
Loading required package: GenomicRanges
Loading required package: stats4
Loading required package: BiocGenerics
Loading required package: parallel
Attaching package: 'BiocGenerics'
The following objects are masked from 'package:parallel':
clusterApply, clusterApplyLB, clusterCall, clusterEvalQ,
clusterExport, clusterMap, parApply, parCapply, parLapply,
parLapplyLB, parRapply, parSapply, parSapplyLB
The following objects are masked from 'package:stats':
IQR, mad, sd, var, xtabs
The following objects are masked from 'package:base':
Filter, Find, Map, Position, Reduce, anyDuplicated, append,
as.data.frame, basename, cbind, colMeans, colSums, colnames,
dirname, do.call, duplicated, eval, evalq, get, grep, grepl,
intersect, is.unsorted, lapply, lengths, mapply, match, mget,
order, paste, pmax, pmax.int, pmin, pmin.int, rank, rbind,
rowMeans, rowSums, rownames, sapply, setdiff, sort, table, tapply,
union, unique, unsplit, which, which.max, which.min
Loading required package: S4Vectors
Attaching package: 'S4Vectors'
The following object is masked from 'package:base':
expand.grid
Loading required package: IRanges
Loading required package: GenomeInfoDb
Loading required package: locfit
locfit 1.5-9.1 2013-03-22
Attaching package: 'EnrichedHeatmap'
The following object is masked from 'package:ComplexHeatmap':
+.AdditiveUnit
Normalize signal to target:
Upstream 10 bp (5 windows)
Downstream 10 bp (5 windows)
Include target regions (1 window)
1 signal region
Normalize signal to target:
Upstream 10 bp (5 windows)
Downstream 10 bp (5 windows)
Include target regions (1 window)
1 signal region
Normalize signal to target:
Upstream 10 bp (5 windows)
Downstream 10 bp (5 windows)
Include target regions (1 window)
1 signal region
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.