Modify a Ranges object

# S3 method for class 'Ranges'
mutate(.data, ...)

Arguments

.data

a Ranges object

...

Pairs of name-value expressions. The name-value pairs can either create new metadata columns or modify existing ones.

Value

a Ranges object

Examples

df <- data.frame(start = 1:10,
                 width = 5,
                 seqnames = "seq1",
                 strand = sample(c("+", "-", "*"), 10, replace = TRUE),
                 gc = runif(10))
rng <- as_granges(df)

# mutate adds new columns
rng %>%
    mutate(avg_gc = mean(gc), row_id = 1:n())
#> GRanges object with 10 ranges and 3 metadata columns:
#>        seqnames    ranges strand |        gc    avg_gc    row_id
#>           <Rle> <IRanges>  <Rle> | <numeric> <numeric> <integer>
#>    [1]     seq1       1-5      + |  0.211409  0.490834         1
#>    [2]     seq1       2-6      - |  0.463701  0.490834         2
#>    [3]     seq1       3-7      - |  0.647101  0.490834         3
#>    [4]     seq1       4-8      + |  0.960573  0.490834         4
#>    [5]     seq1       5-9      * |  0.676398  0.490834         5
#>    [6]     seq1      6-10      * |  0.445148  0.490834         6
#>    [7]     seq1      7-11      + |  0.357774  0.490834         7
#>    [8]     seq1      8-12      * |  0.455731  0.490834         8
#>    [9]     seq1      9-13      + |  0.445414  0.490834         9
#>   [10]     seq1     10-14      * |  0.245093  0.490834        10
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
# can also compute on newly created columns
rng %>%
    mutate(score = gc * width, score2 = score + 1)
#> GRanges object with 10 ranges and 3 metadata columns:
#>        seqnames    ranges strand |        gc     score    score2
#>           <Rle> <IRanges>  <Rle> | <numeric> <numeric> <numeric>
#>    [1]     seq1       1-5      + |  0.211409   1.05704   2.05704
#>    [2]     seq1       2-6      - |  0.463701   2.31851   3.31851
#>    [3]     seq1       3-7      - |  0.647101   3.23551   4.23551
#>    [4]     seq1       4-8      + |  0.960573   4.80287   5.80287
#>    [5]     seq1       5-9      * |  0.676398   3.38199   4.38199
#>    [6]     seq1      6-10      * |  0.445148   2.22574   3.22574
#>    [7]     seq1      7-11      + |  0.357774   1.78887   2.78887
#>    [8]     seq1      8-12      * |  0.455731   2.27866   3.27866
#>    [9]     seq1      9-13      + |  0.445414   2.22707   3.22707
#>   [10]     seq1     10-14      * |  0.245093   1.22546   2.22546
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
# group by partitions the data and computes within each group
rng %>%
    group_by(strand) %>%
    mutate(avg_gc = mean(gc), row_id = 1:n())
#> GRanges object with 10 ranges and 3 metadata columns:
#> Groups: strand [3]
#>        seqnames    ranges strand |        gc    avg_gc    row_id
#>           <Rle> <IRanges>  <Rle> | <numeric> <numeric> <integer>
#>    [1]     seq1       1-5      + |  0.211409  0.493792         1
#>    [2]     seq1       2-6      - |  0.463701  0.555401         1
#>    [3]     seq1       3-7      - |  0.647101  0.555401         2
#>    [4]     seq1       4-8      + |  0.960573  0.493792         2
#>    [5]     seq1       5-9      * |  0.676398  0.455593         1
#>    [6]     seq1      6-10      * |  0.445148  0.455593         2
#>    [7]     seq1      7-11      + |  0.357774  0.493792         3
#>    [8]     seq1      8-12      * |  0.455731  0.455593         3
#>    [9]     seq1      9-13      + |  0.445414  0.493792         4
#>   [10]     seq1     10-14      * |  0.245093  0.455593         4
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

# mutate can be used in conjuction with anchoring to resize ranges
rng %>%
    mutate(width = 10)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1      1-10      + |  0.211409
#>    [2]     seq1      2-11      - |  0.463701
#>    [3]     seq1      3-12      - |  0.647101
#>    [4]     seq1      4-13      + |  0.960573
#>    [5]     seq1      5-14      * |  0.676398
#>    [6]     seq1      6-15      * |  0.445148
#>    [7]     seq1      7-16      + |  0.357774
#>    [8]     seq1      8-17      * |  0.455731
#>    [9]     seq1      9-18      + |  0.445414
#>   [10]     seq1     10-19      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
# by default width modfication fixes by start
rng %>%
    anchor_start() %>%
    mutate(width = 10)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1      1-10      + |  0.211409
#>    [2]     seq1      2-11      - |  0.463701
#>    [3]     seq1      3-12      - |  0.647101
#>    [4]     seq1      4-13      + |  0.960573
#>    [5]     seq1      5-14      * |  0.676398
#>    [6]     seq1      6-15      * |  0.445148
#>    [7]     seq1      7-16      + |  0.357774
#>    [8]     seq1      8-17      * |  0.455731
#>    [9]     seq1      9-18      + |  0.445414
#>   [10]     seq1     10-19      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
# fix by end or midpoint
rng %>%
    anchor_end() %>%
    mutate(width = width + 1)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1       0-5      + |  0.211409
#>    [2]     seq1       1-6      - |  0.463701
#>    [3]     seq1       2-7      - |  0.647101
#>    [4]     seq1       3-8      + |  0.960573
#>    [5]     seq1       4-9      * |  0.676398
#>    [6]     seq1      5-10      * |  0.445148
#>    [7]     seq1      6-11      + |  0.357774
#>    [8]     seq1      7-12      * |  0.455731
#>    [9]     seq1      8-13      + |  0.445414
#>   [10]     seq1      9-14      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
rng %>%
    anchor_center() %>%
    mutate(width = width + 1)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1       0-5      + |  0.211409
#>    [2]     seq1       1-6      - |  0.463701
#>    [3]     seq1       2-7      - |  0.647101
#>    [4]     seq1       3-8      + |  0.960573
#>    [5]     seq1       4-9      * |  0.676398
#>    [6]     seq1      5-10      * |  0.445148
#>    [7]     seq1      6-11      + |  0.357774
#>    [8]     seq1      7-12      * |  0.455731
#>    [9]     seq1      8-13      + |  0.445414
#>   [10]     seq1      9-14      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
# anchoring by strand
rng %>%
    anchor_3p() %>%
    mutate(width = width * 2)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1      -4-5      + |  0.211409
#>    [2]     seq1      2-11      - |  0.463701
#>    [3]     seq1      3-12      - |  0.647101
#>    [4]     seq1      -1-8      + |  0.960573
#>    [5]     seq1       0-9      * |  0.676398
#>    [6]     seq1      1-10      * |  0.445148
#>    [7]     seq1      2-11      + |  0.357774
#>    [8]     seq1      3-12      * |  0.455731
#>    [9]     seq1      4-13      + |  0.445414
#>   [10]     seq1      5-14      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths
rng %>%
    anchor_5p() %>%
    mutate(width = width * 2)
#> GRanges object with 10 ranges and 1 metadata column:
#>        seqnames    ranges strand |        gc
#>           <Rle> <IRanges>  <Rle> | <numeric>
#>    [1]     seq1      1-10      + |  0.211409
#>    [2]     seq1      -3-6      - |  0.463701
#>    [3]     seq1      -2-7      - |  0.647101
#>    [4]     seq1      4-13      + |  0.960573
#>    [5]     seq1      5-14      * |  0.676398
#>    [6]     seq1      6-15      * |  0.445148
#>    [7]     seq1      7-16      + |  0.357774
#>    [8]     seq1      8-17      * |  0.455731
#>    [9]     seq1      9-18      + |  0.445414
#>   [10]     seq1     10-19      * |  0.245093
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths