Skip to contents

This man page documents intra range transformations of a linkSet object.

Usage

# S4 method for class 'linkSet'
trim(x, use.names = TRUE, ...)

# S4 method for class 'linkSet'
resize(x, width, fix = "start", use.names = TRUE, ...)

# S4 method for class 'linkSet'
resizeRegions(
  x,
  width = 1000,
  fix = "start",
  use.names = TRUE,
  region = "both",
  ...
)

# S4 method for class 'linkSet'
narrow(x, start = NA, end = NA, width = NA, use.names = TRUE)

# S4 method for class 'linkSet'
narrowRegions(
  x,
  start = NA,
  end = NA,
  width = NA,
  use.names = TRUE,
  region = "both"
)

# S4 method for class 'linkSet'
shift(x, shift = 0L, use.names = TRUE)

# S4 method for class 'linkSet'
shiftRegions(x, shift = 0L, use.names = TRUE, region = "both")

# S4 method for class 'linkSet'
flank(
  x,
  width,
  start = TRUE,
  both = FALSE,
  use.names = TRUE,
  ignore.strand = FALSE
)

# S4 method for class 'linkSet'
flankRegions(
  x,
  width,
  start = TRUE,
  both = FALSE,
  use.names = TRUE,
  ignore.strand = FALSE,
  region = "both"
)

# S4 method for class 'linkSet'
promoters(x, upstream = 2000, downstream = 200, use.names = TRUE)

# S4 method for class 'linkSet'
promoterRegions(
  x,
  upstream = 2000,
  downstream = 200,
  use.names = TRUE,
  region = "both"
)

# S4 method for class 'linkSet'
width(x)

# S4 method for class 'linkSet'
reduce(x, drop.empty.ranges = FALSE, ...)

Arguments

x

A linkSet object

use.names

A logical indicating whether to use names

...

Additional arguments passed to the GenomicRanges trim method

width

The desired width of the output ranges

fix

The anchor point for resizing operations ("start", "end", or "center")

region

Which regions to modify ("both", "bait", or "oe")

start, end

The desired start and end coordinates for narrowing

shift

The number of positions to shift

both

Whether to get flanking regions on both sides

ignore.strand

TRUE or FALSE. Whether the strand of the input ranges should be ignored or not. See details below.

upstream, downstream

Number of bases upstream/downstream for promoter regions

drop.empty.ranges

Whether to drop empty ranges when reducing

Value

A linkSet object

Author

Gilbert Han

Examples

data(linkExample)
resize_bait <- resizeRegions(linkExample, width = 75, fix = "start", region = "bait")
resize_bait
#> linkSet object with 5 interactions and 1 metadata column:
#>              bait     seqnames_oe ranges_oe | anchor1.symbol
#>       <character>           <Rle> <IRanges> |    <character>
#>   [1]       Gene1 ---        chr1     50-59 |          Gene1
#>   [2]       Gene1 ---        chr2   150-159 |          Gene1
#>   [3]       Gene2 ---        chr2   250-259 |          Gene2
#>   [4]       Gene3 ---        chr4   350-359 |          Gene3
#>   [5]       Gene3 ---        chr4   450-459 |          Gene3
#>   -------
#>   regions: 10 ranges and 0 metadata columns
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

narrow_bait <- narrowRegions(linkExample, start = 1, width = 5, region = "bait")
narrow_bait
#> linkSet object with 5 interactions and 1 metadata column:
#>              bait     seqnames_oe ranges_oe | anchor1.symbol
#>       <character>           <Rle> <IRanges> |    <character>
#>   [1]       Gene1 ---        chr1     50-59 |          Gene1
#>   [2]       Gene1 ---        chr2   150-159 |          Gene1
#>   [3]       Gene2 ---        chr2   250-259 |          Gene2
#>   [4]       Gene3 ---        chr4   350-359 |          Gene3
#>   [5]       Gene3 ---        chr4   450-459 |          Gene3
#>   -------
#>   regions: 10 ranges and 0 metadata columns
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

shift_oe <- shiftRegions(linkExample, shift = 10, region = "oe")
shift_oe
#> linkSet object with 5 interactions and 1 metadata column:
#>              bait     seqnames_oe ranges_oe | anchor1.symbol
#>       <character>           <Rle> <IRanges> |    <character>
#>   [1]       Gene1 ---        chr1     60-69 |          Gene1
#>   [2]       Gene1 ---        chr2   160-169 |          Gene1
#>   [3]       Gene2 ---        chr2   260-269 |          Gene2
#>   [4]       Gene3 ---        chr4   360-369 |          Gene3
#>   [5]       Gene3 ---        chr4   460-469 |          Gene3
#>   -------
#>   regions: 10 ranges and 0 metadata columns
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

flank_bait <- flankRegions(linkExample, width = 100, start = TRUE, both = FALSE, use.names = TRUE, ignore.strand = FALSE, region = "bait")
flank_bait
#> linkSet object with 5 interactions and 1 metadata column:
#>              bait     seqnames_oe ranges_oe | anchor1.symbol
#>       <character>           <Rle> <IRanges> |    <character>
#>   [1]       Gene1 ---        chr1     50-59 |          Gene1
#>   [2]       Gene1 ---        chr2   150-159 |          Gene1
#>   [3]       Gene2 ---        chr2   250-259 |          Gene2
#>   [4]       Gene3 ---        chr4   350-359 |          Gene3
#>   [5]       Gene3 ---        chr4   450-459 |          Gene3
#>   -------
#>   regions: 10 ranges and 0 metadata columns
#>   seqinfo: 4 sequences from an unspecified genome; no seqlengths

width(linkExample)
#> $bait
#> [1] 10 10 10 10 10
#> 
#> $oe
#> [1] 10 10 10 10 10
#>