Search for telomeres
telomeres.Rd
Searches for telomeric repeats in the whole list of tandem repeats, or narrowing the search on the ends of sequence features.
Usage
telomeres(gr, motif = "TTAGGG", narrow = c("no", "ends"))
Arguments
- gr
A
GenomicRanges::GRanges
object produced by theloadTantan
function- motif
The telomeric sequence motif, at the end of the plus strand, like in TeloBase and the Telomerase Database.
- narrow
Narrow the search on the ends.
Examples
telomeres(exampleTantan)
#> GRanges object with 8 ranges and 5 metadata columns:
#> seqnames ranges strand | score normTandem nchar
#> <Rle> <IRanges> <Rle> | <numeric> <factor> <integer>
#> [1] 1086a8205b 1-613 + | 102.66700 AACCCT 6
#> [2] 1086a8205b 137225-137239 + | 2.50000 AACCCT 6
#> [3] 1086a8205b 1005769-1005789 + | 3.50000 AACCCT 6
#> [4] 1086a8205b 1284710-1284733 - | 4.00000 AACCCT 6
#> [5] 1086a8205b 1286476-1286494 - | 3.16667 AACCCT 6
#> [6] 1086a8205b 1288480-1288495 + | 2.66667 AACCCT 6
#> [7] 1086a8205b 1505855-1505871 - | 2.83333 AACCCT 6
#> [8] 1086a8205b 1663642-1663664 + | 3.83333 AACCCT 6
#> rstart rend
#> <numeric> <numeric>
#> [1] 2232853 2233465
#> [2] 2096227 2096241
#> [3] 1227677 1227697
#> [4] 948733 948756
#> [5] 946972 946990
#> [6] 944971 944986
#> [7] 727595 727611
#> [8] 569802 569824
#> -------
#> seqinfo: 1 sequence from example genome
telomeres(exampleTantan, narrow="ends")
#> GRanges object with 1 range and 5 metadata columns:
#> seqnames ranges strand | score normTandem nchar rstart
#> <Rle> <IRanges> <Rle> | <numeric> <factor> <integer> <numeric>
#> [1] 1086a8205b 1-613 + | 102.667 AACCCT 6 2232853
#> rend
#> <numeric>
#> [1] 2233465
#> -------
#> seqinfo: 1 sequence from example genome