Retrieve genomic elements from regulonDB

get_dna_objects(
  regulondb,
  genome = "eschColi_K12",
  grange = GRanges("chr", IRanges(1, 5000)),
  elements = "gene"
)

Arguments

regulondb

A regulondb() object.

genome

A valid UCSC genome name.

grange

A GenomicRanges::GRanges-class() object indicating position left and right.

elements

A character vector specifying which annotation elements to plot. It can be any from: "-10 promoter box", "-35 promoter box", "gene", "promoter", "Regulatory Interaction", "sRNA interaction", or "terminator".

Value

GenomicRanges::GRanges-class() object with the elements found.

Author

Joselyn Chavez

Examples

## Connect to the RegulonDB database if necessary
if (!exists("regulondb_conn")) {
      regulondb_conn <- connect_database()
  }
#> snapshotDate(): 2021-10-20

## Build the regulondb object
e_coli_regulondb <-
    regulondb(
        database_conn = regulondb_conn,
        organism = "chr",
        database_version = "1",
        genome_version = "1"
    )

## Get all genes from E. coli
get_dna_objects(e_coli_regulondb)
#> GRanges object with 3 ranges and 4 metadata columns:
#>       seqnames    ranges strand |           id        type        name
#>          <Rle> <IRanges>  <Rle> |  <character> <character> <character>
#>   [1]      chr  337-2799      + | ECK120000987        gene        thrA
#>   [2]      chr 2801-3733      + | ECK120000988        gene        thrB
#>   [3]      chr   190-255      + | ECK120001251        gene        thrL
#>                  description
#>                  <character>
#>   [1] fused aspartate kina..
#>   [2]      homoserine kinase
#>   [3] <i>thr</i> operon le..
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

## Get genes providing Genomic Ranges
grange <- GenomicRanges::GRanges(
    "chr",
    IRanges::IRanges(5000, 10000)
)
get_dna_objects(e_coli_regulondb, grange)
#> GRanges object with 3 ranges and 4 metadata columns:
#>       seqnames    ranges strand |           id        type        name
#>          <Rle> <IRanges>  <Rle> |  <character> <character> <character>
#>   [1]      chr  337-2799      + | ECK120000987        gene        thrA
#>   [2]      chr 2801-3733      + | ECK120000988        gene        thrB
#>   [3]      chr   190-255      + | ECK120001251        gene        thrL
#>                  description
#>                  <character>
#>   [1] fused aspartate kina..
#>   [2]      homoserine kinase
#>   [3] <i>thr</i> operon le..
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths

## Get aditional elements within genomic positions
get_dna_objects(e_coli_regulondb,
    grange,
    elements = c("gene", "promoter")
)
#> GRanges object with 19 ranges and 4 metadata columns:
#>        seqnames    ranges strand |           id        type
#>           <Rle> <IRanges>  <Rle> |  <character> <character>
#>    [1]      chr  337-2799      + | ECK120000987        gene
#>    [2]      chr 2801-3733      + | ECK120000988        gene
#>    [3]      chr   190-255      + | ECK120001251        gene
#>    [4]      chr       148      + | ECK120010236    promoter
#>    [5]      chr        38      + | ECK125230824    promoter
#>    ...      ...       ...    ... .          ...         ...
#>   [15]      chr      3066      + | ECK125230834    promoter
#>   [16]      chr      3368      + | ECK125230835    promoter
#>   [17]      chr      3396      + | ECK125230836    promoter
#>   [18]      chr      3401      + | ECK125230837    promoter
#>   [19]      chr      4241      + | ECK125230838    promoter
#>                          name            description
#>                   <character>            <character>
#>    [1]                   thrA fused aspartate kina..
#>    [2]                   thrB      homoserine kinase
#>    [3]                   thrL <i>thr</i> operon le..
#>    [4] thrLp promoter with .. thrLp promoter with ..
#>    [5] TSS_1 promoter with .. TSS_1 promoter with ..
#>    ...                    ...                    ...
#>   [15] TSS_12 promoter with.. TSS_12 promoter with..
#>   [16] TSS_13 promoter with.. TSS_13 promoter with..
#>   [17] TSS_14 promoter with.. TSS_14 promoter with..
#>   [18] TSS_15 promoter with.. TSS_15 promoter with..
#>   [19] TSS_16 promoter with.. TSS_16 promoter with..
#>   -------
#>   seqinfo: 1 sequence from an unspecified genome; no seqlengths