Description Usage Arguments Value Author(s) References See Also Examples
View source: R/MapRangesToGenomicIntervals.R
MapRangesToGenomicIntervals collapses all the chromosomes to their subsets covered by where.to.map and then liftover the what.to.map to this genome subset. The function is intended to prepare an estimation correlation of two features that are spatially restricted in the genome. We map them both to their location areas (the where.to.map
is common for query and for the rerefence), and then we calculate the GenometriCorrelation
between the mapped fatures.)
1 2 3 4 5 6 7 8 9 | MapRangesToGenomicIntervals(
where.to.map,
what.to.map,
chrom.suffix = "_mapped",
chromosomes.to.proceed = NA,
chromosomes.length = c(),
unmapped.chromosome.warning = TRUE,
nonnormalised.mapping.warning = TRUE
)
|
where.to.map |
The set of genomic intervals we map to. |
what.to.map |
The set of ranges that we map. |
chrom.suffix |
The suffix to be appended to all the sestination chromosome names in the mapping; default is "_mapped". |
chromosomes.to.proceed |
The default set of chromosomes to map is the intersection of the chromosomes in where.to.map and what.to.map. If we want to restrict the set, we can do it with this parameter. |
chromosomes.length |
is an alternative to seqingo() ot the GRanges of where.to.map way to pass the lengths of chromosomes to the mapping routine |
unmapped.chromosome.warning |
For each chromosome that is represented in |
nonnormalised.mapping.warning |
If the input mapping space is not normalised (e.g. contains overlapping intervals), it is normalised before the mapping. A warning is generated if |
GRanges object that is a liftover of what.to.map
to subgenome covered by where.to.map
Alexander Favorov favorov@sensi.org, Loris Mularoni, Yulia Medvedeva, Harris A. Jaffee, Ekaterina V. Zhuravleva, Veronica Busa, Leslie M. Cope, Andrey A. Mironov, Vsevolod J. Makeev, Sarah J. Wheelan.
The GenometriCorr
documentation and vignette.
1 2 3 4 | library('GenometriCorr')
intervals<-GRanges(ranges=IRanges(c(1,10001,1,10001),width=c(1000)),seqnames=c('chr1','chr1','chr2','chr2'),seqlengths=c('chr1'=300000,'chr2'=300000))
ranges=GRanges(ranges=IRanges(c(10,110,10101,100000,500,550,1055),width=c(10)),seqnames=c(rep('chr1',4),rep('chr2',3)),seqlengths=c('chr1'=300000,'chr2'=300000))
mapped<-MapRangesToGenomicIntervals(where.to.map=intervals,what.to.map=ranges)
|
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.