1 Introduction

biodbUniprot is a biodb extension package that implements a connector to Uniprot database.

The UniProt Knowledge Base (Consortium 2016) can be searched using its search web service.

We present here the way to contact this web service with this package.

2 Installation

Install using Bioconductor:

if (!requireNamespace("BiocManager", quietly=TRUE))
    install.packages("BiocManager")
BiocManager::install('biodbUniprot')

3 Initialization

The first step in using biodbUniprot, is to create an instance of the biodb class BiodbMain from the main biodb package. This is done by calling the constructor of the class:

mybiodb <- biodb::newInst()

During this step the configuration is set up, the cache system is initialized and extension packages are loaded.

We will see at the end of this vignette that the biodb instance needs to be terminated with a call to the terminate() method.

4 Creating a connector to Uniprot database

In biodb the connection to a database is handled by a connector instance that you can get from the factory. biodbUniprot implements a connector to a remote database. Here is the code to instantiate a connector:

conn <- mybiodb$getFactory()$createConn('uniprot')
## Loading required package: biodbUniprot

5 Getting entries

To download entries, run the getEntry(), which returns a list of BiodbEntry objects:

entries <- conn$getEntry(c('P01011', 'P09237'))

To print the information contained in the entry objects as a data frame, run the entriesToDataframe() method attached to the BiodbMain instance:

mybiodb$entriesToDataframe(entries)
##   accession               gene.symbol kegg.genes.id molecular.mass
## 1    P01011 SERPINA3;AACT;GIG24;GIG25        hsa:12          47651
## 2    P09237          MMP7;MPSL1;PUMP1      hsa:4316          29677
##                                                                                                                            name
## 1 AACT_HUMAN;Alpha-1-antichymotrypsin;Cell growth-inhibiting gene 24/25 protein;Serpin A3;Alpha-1-antichymotrypsin His-Pro-less
## 2                             MMP7_HUMAN;Matrilysin;Matrin;Matrix metalloproteinase-7;Pump-1 protease;Uterine metalloproteinase
##   ncbi.gene.id
## 1           12
## 2         4316
##                                                                                                                                                                                                                                                                                                                                                                                                                                    aa.seq
## 1 MERMLPLLALGLLAAGFCPAVLCHPNSPLDEENLTQENQDRGTHVDLGLASANVDFAFSLYKQLVLKAPDKNVIFSPLSISTALAFLSLGAHNTTLTEILKGLKFNLTETSEAEIHQSFQHLLRTLNQSSDELQLSMGNAMFVKEQLSLLDRFTEDAKRLYGSEAFATDFQDSAAAKKLINDYVKNGTRGKITDLIKDLDSQTMMVLVNYIFFKAKWEMPFDPQDTHQSRFYLSKKKWVMVPMMSLHHLTIPYFRDEELSCTVVELKYTGNASALFILPDQDKMEEVEAMLLPETLKRWRDSLEFREIGELYLPKFSISRDYNLNDILLQLGIEEAFTSKADLSGITGARNLAVSQVVHKAVLDVFEEGTEASAATAVKITLLSALVETRTIVRFNRPFLMIIVPTDTQNIFFMSKVTNPKQA
## 2                                                                                                                                                             MRLTVLCAVCLLPGSLALPLPQEAGGMSELQWEQAQDYLKRFYLYDSETKNANSLEAKLKEMQKFFGLPITGMLNSRVIEIMQKPRCGVPDVAEYSLFPNSPKWTSKVVTYRIVSYTRDLPHITVDRLVSKALNMWGKEIPLHFRKVVWGTADIMIGFARGAHGDSYPFDGPGNTLAHAFAPGTGLGGDAHFDEDERWTDGSSLGINFLYAATHELGHSLGMGHSSDPNAVMYPTYGNGDPQNFKLSQDDIKGIQKLYGKRSNSRKK
##   aa.seq.length uniprot.id        ec expasy.enzyme.id
## 1           423     P01011      <NA>             <NA>
## 2           267     P09237 3.4.24.23        3.4.24.23

6 Using the search web service

The method wsSearch() (wsQuery() is now deprecated) implements the request to the search web service, and the parsing of its output.

To get the raw results returned by the UniProt server, run the following code:

conn$wsSearch('reviewed:true AND organism_id:9606', fields=c('accession', 'id'),
    size=2, retfmt='plain')
## [1] "Entry\tEntry Name\nA0A0C5B5G6\tMOTSC_HUMAN\nA0A1B0GTW7\tCIROP_HUMAN\n"

The first parameter is the query, as required by the web service. To learn how to write a query for UniProt, see a description of the query web service at http://www.uniprot.org/help/api_queries.

The fields parameter is the fields you want back for each entry returned by the database.

The size parameter is the maximum number of entries the server must return.

The retfmt parameter controls the type of output desired. Here "plain" states that we want the raw output from the server.

To get the output parsed by biodb and get a data frame, run:

conn$wsSearch('reviewed:true AND organism_id:9606', fields=c('accession', 'id'),
    size=2, retfmt='parsed')
##        Entry  Entry Name
## 1 A0A0C5B5G6 MOTSC_HUMAN
## 2 A0A1B0GTW7 CIROP_HUMAN

To get only the list of UniProt identifiers, run:

conn$wsSearch('reviewed:true AND organism_id:9606', fields=c('accession', 'id'),
    size=2, retfmt='ids')
## [1] "A0A0C5B5G6" "A0A1B0GTW7"

And if you are curious to see the URL request that is sent to the server, run:

conn$wsSearch('reviewed:true AND organism_id:9606', fields=c('accession', 'id'),
    size=2, retfmt='request')
## Biodb request object on https://rest.uniprot.org/uniprotkb/search?query=reviewed%3Atrue%20AND%20organism_id%3A9606&fields=accession%2Cid&format=tsv&size=2

7 Conversion of gene symbols to UniProt IDs

The method geneSymbolToUniprotIds() uses wsSearch() to search for UniProt entries that reference particular gene symbols.

For instance, if you want to get the UniProt entries that have the gene symbol G-CSF, just run:

ids <- conn$geneSymbolToUniprotIds('G-CSF')
mybiodb$entryIdsToDataframe(ids[['G-CSF']], 'uniprot', fields=c('accession',
    'gene.symbol'))
##    accession gene.symbol
## 1     Q8MKE0       G-CSF
## 2     Q9GJU0       G-CSF
## 3 A0A679AQ73       g-csf
## 4     C0STS2     G-CSF 2
## 5     C0STS3     G-CSF 1
## 6     Q4H432  GCSF;G-CSF

If you want to match also GCSF (no minus sign character), then run:

ids <- conn$geneSymbolToUniprotIds('G-CSF', ignore.nonalphanum=TRUE)
mybiodb$entryIdsToDataframe(ids[['G-CSF']], 'uniprot', fields=c('accession',
    'gene.symbol'))
##     accession        gene.symbol
## 1      P09919 CSF3;C17orf33;GCSF
## 2      P35833          CSF3;GCSF
## 3      B5L332    csf3;Gcsf;csf3b
## 4      Q8MKE0              G-CSF
## 5      B8ZHI7    gcsf;csf3;csf3a
## 6      Q9GJU0              G-CSF
## 7      Q4H432         GCSF;G-CSF
## 8  A0A679AQ73              g-csf
## 9      C0STS2            G-CSF 2
## 10     C0STS3            G-CSF 1
## 11 A0A2Z6I9R9               GCSF
## 12 A0A8M1NPY2    csf3a;csf3;gcsf
## 13 A0A8M6Z8U5    csf3a;csf3;gcsf

If you want to match G-CSFa2 too, run:

ids <- conn$geneSymbolToUniprotIds('G-CSF', partial.match=TRUE)
mybiodb$entryIdsToDataframe(ids[['G-CSF']], 'uniprot', fields=c('accession',
    'gene.symbol'))
##    accession gene.symbol
## 1     Q8MKE0       G-CSF
## 2     Q9GJU0       G-CSF
## 3 A0A679AQ73       g-csf
## 4     C0STS2     G-CSF 2
## 5     C0STS3     G-CSF 1
## 6     Q4H432  GCSF;G-CSF

The way this method works is by running wsSearch() to get a first set of entry identifiers, and then download each entry and apply a filtering on them. The downloading of the entries may quite long, wsSearch() returning potentially thousands of entries, each entry being downloaded with a single separate request and the frequency limit being 3 request per second. Entries already in cache or memory will not be downloaded again, so running the same request a second time will be faster, as it is usually the case with biodb.

8 Closing biodb instance

When done with your biodb instance you have to terminate it, in order to ensure release of resources (file handles, database connection, etc):

mybiodb$terminate()
## INFO  [15:57:49.858] Closing BiodbMain instance...
## INFO  [15:57:49.860] Connector "uniprot" deleted.

9 Session information

sessionInfo()
## R version 4.2.1 (2022-06-23)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 20.04.5 LTS
## 
## Matrix products: default
## BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
## 
## locale:
##  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
##  [3] LC_TIME=en_GB              LC_COLLATE=C              
##  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
##  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
##  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
## [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
## 
## attached base packages:
## [1] stats     graphics  grDevices utils     datasets  methods   base     
## 
## other attached packages:
## [1] biodbUniprot_1.2.1 BiocStyle_2.24.0  
## 
## loaded via a namespace (and not attached):
##  [1] progress_1.2.2      tidyselect_1.1.2    xfun_0.33          
##  [4] bslib_0.4.0         purrr_0.3.4         vctrs_0.4.1        
##  [7] generics_0.1.3      htmltools_0.5.3     BiocFileCache_2.4.0
## [10] yaml_2.3.5          utf8_1.2.2          blob_1.2.3         
## [13] XML_3.99-0.10       rlang_1.0.6         jquerylib_0.1.4    
## [16] pillar_1.8.1        withr_2.5.0         glue_1.6.2         
## [19] DBI_1.1.3           rappdirs_0.3.3      bit64_4.0.5        
## [22] dbplyr_2.2.1        lifecycle_1.0.2     plyr_1.8.7         
## [25] stringr_1.4.1       memoise_2.0.1       evaluate_0.16      
## [28] knitr_1.40          fastmap_1.1.0       curl_4.3.2         
## [31] fansi_1.0.3         biodb_1.4.2         Rcpp_1.0.9         
## [34] openssl_2.0.3       filelock_1.0.2      BiocManager_1.30.18
## [37] cachem_1.0.6        jsonlite_1.8.0      bit_4.0.4          
## [40] chk_0.8.1           askpass_1.1         hms_1.1.2          
## [43] digest_0.6.29       stringi_1.7.8       bookdown_0.29      
## [46] dplyr_1.0.10        bitops_1.0-7        cli_3.4.1          
## [49] tools_4.2.1         magrittr_2.0.3      sass_0.4.2         
## [52] RCurl_1.98-1.8      RSQLite_2.2.17      tibble_3.1.8       
## [55] crayon_1.5.1        pkgconfig_2.0.3     ellipsis_0.3.2     
## [58] prettyunits_1.1.1   assertthat_0.2.1    rmarkdown_2.16     
## [61] httr_1.4.4          lgr_0.4.4           R6_2.5.1           
## [64] compiler_4.2.1

References

Consortium, The UniProt. 2016. “UniProt: the universal protein knowledgebase.” Nucleic Acids Research 45 (D1): D158–D169. https://doi.org/10.1093/nar/gkw1099.