NestLink 1.0.0
The following content is described in more detail in Egloff et al. (2018), (under review NMETH-A35040).
library(NestLink)
library(ExperimentHub)
eh <- ExperimentHub()
## snapshotDate(): 2019-04-29
query(eh, "NestLink")
## ExperimentHub with 8 records
## # snapshotDate(): 2019-04-29
## # $dataprovider: Functional Genomics Center Zurich (FGCZ)
## # $species: NA
## # $rdataclass: data.frame, DNAStringSet
## # additional mcols(): taxonomyid, genome, description,
## # coordinate_1_based, maintainer, rdatadateadded, preparerclass,
## # tags, rdatapath, sourceurl, sourcetype
## # retrieve records with, e.g., 'object[["EH2063"]]'
##
## title
## EH2063 | Sample NGS NB FC linkage data
## EH2064 | Flycodes tryptic digested
## EH2065 | Nanobodies tryptic digested
## EH2066 | FASTA as ground-truth for unit testing
## EH2067 | Known nanobodies
## EH2068 | Quantitaive results for SMEG and COLI
## EH2069 | F255744 Mascot Search result
## EH2070 | WU160118 Mascot Search results
# dataFolder <- file.path(path.package(package = 'NestLink'), 'extdata')
# expFile <- list.files(dataFolder, pattern='*.fastq.gz', full.names = TRUE)
expFile <- query(eh, c("NestLink", "NL42_100K.fastq.gz"))[[1]]
## see ?NestLink and browseVignettes('NestLink') for documentation
## downloading 0 resources
## loading from cache
## 'EH2063 : 2063'
scratchFolder <- tempdir()
setwd(scratchFolder)
For data QC some known NB were spiked in. Here, we load the NB DNA sequences and translate them to the corresponding AA sequences.
# knownNB_File <- list.files(dataFolder,
# pattern='knownNB.txt', full.names = TRUE)
knownNB_File <- query(eh, c("NestLink", "knownNB.txt"))[[1]]
## see ?NestLink and browseVignettes('NestLink') for documentation
## downloading 0 resources
## loading from cache
## 'EH2067 : 2067'
knownNB_data <- read.table(knownNB_File, sep='\t',
header = TRUE, row.names = 1, stringsAsFactors = FALSE)
knownNB <- Biostrings::translate(DNAStringSet(knownNB_data$Sequence))
names(knownNB) <- rownames(knownNB_data)
knownNB <- sapply(knownNB, toString)
The workflow uses the first 100 reads only for a rapid processing time.
param <- list()
param[['nReads']] <- 100 #Number of Reads from the start of fastq file to process
param[['maxMismatch']] <- 1 #Number of accepted mismatches for all pattern search steps
param[['NB_Linker1']] <- "GGCCggcggGGCC" #Linker Sequence left to nanobody
param[['NB_Linker2']] <- "GCAGGAGGA" #Linker Sequence right to nanobody
param[['ProteaseSite']] <- "TTAGTCCCAAGA" #Sequence next to flycode
param[['FC_Linker']] <- "GGCCaaggaggcCGG" #Linker Sequence next to flycode
param[['knownNB']] <- knownNB
param[['minRelBestHitFreq']] <- 0.8 #minimal fraction of the dominant nanobody for a specific flycode
param[['minConsensusScore']] <- 0.9 #minimal fraction per sequence position in nanabody consensus sequence calculation
param[['minNanobodyLength']] <- 348 #minimal nanobody length in [nt]
param[['minFlycodeLength']] <- 33 #minimal flycode length in [nt]
param[['FCminFreq']] <- 1 #minimal number of subreads for a specific flycode to keep it in the analysis
The following steps are included:
system.time(NB2FC <- runNGSAnalysis(file = expFile[1], param))
## user system elapsed
## 3.292 0.104 3.542
head(NB2FC, 2)
## NB
## 1 SQVQLVESGGGLVQAGGSLRLSCAASGFPVEAHRMYWYRQAPGKEREWVAAISSKGQQTWYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDYGWYYGDYDYWGQGTQVTVS
## 2 SQVQLVESGGGLVQAGGSLRLSCAASGFPVSWTKMYWYRQAPGKEREWVAAIWSTGSWTKYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDKGHQHAHYDYWGQGTQVTVS
## FlycodeCount
## 1 29
## 2 3
## AssociatedFlycodes
## 1 GSAAATAVTWR,GSADGQETDWR,GSADVPEAVWLTVR,GSAPTAPVSWQEGGR,GSAVDPVTVWLTVR,GSDAEGVAAWQSR,GSDAEYTTAWR,GSDDTDETDWR,GSDEAEEEGWQEGGR,GSDPGTDDEWQSR,GSDTEDWEEWQSR,GSDVWDTAVWLTVR,GSEGTDAVGWLTVR,GSEPASEVVWQEGGR,GSEPDVYTAWLTVR,GSEVLDGDEWR,GSFVASFAVWLTVR,GSGDVEGEAWQEGGR,GSGPDPPYGWLR,GSPAVDPPVWLTVR,GSPDEVEVVWLTVR,GSPDSPPAYWLTVR,GSPTVVTFLWR,GSQYTLTPTWLTVR,GSSDAASPSWLTVR,GSTGEDGVVWLTVR,GSTVVTSDPWLTVR,GSVDDQPDTWQEGGR,GSYTPGSTSWQSR
## 2 GSADFPVVAWLR,GSAEVDEADWQEGGR,GSEPDVAAGWQSR
## NB_Name
## 1
## 2
head(nanobodyFlycodeLinking.as.fasta(NB2FC))
## [1] ">NB0001 FC29 SQVQLVESGGGLVQAGGSLRLSCAASGFPVEAHRMYWYRQAPGKEREWVAAISSKGQQTWYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDYGWYYGDYDYWGQGTQVTVS\nGSAAATAVTWRGSADGQETDWRGSADVPEAVWLTVRGSAPTAPVSWQEGGRGSAVDPVTVWLTVRGSDAEGVAAWQSRGSDAEYTTAWRGSDDTDETDWRGSDEAEEEGWQEGGRGSDPGTDDEWQSRGSDTEDWEEWQSRGSDVWDTAVWLTVRGSEGTDAVGWLTVRGSEPASEVVWQEGGRGSEPDVYTAWLTVRGSEVLDGDEWRGSFVASFAVWLTVRGSGDVEGEAWQEGGRGSGPDPPYGWLRGSPAVDPPVWLTVRGSPDEVEVVWLTVRGSPDSPPAYWLTVRGSPTVVTFLWRGSQYTLTPTWLTVRGSSDAASPSWLTVRGSTGEDGVVWLTVRGSTVVTSDPWLTVRGSVDDQPDTWQEGGRGSYTPGSTSWQSR\n"
## [2] ">NB0002 FC3 SQVQLVESGGGLVQAGGSLRLSCAASGFPVSWTKMYWYRQAPGKEREWVAAIWSTGSWTKYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDKGHQHAHYDYWGQGTQVTVS\nGSADFPVVAWLRGSAEVDEADWQEGGRGSEPDVAAGWQSR\n"
## [3] ">NB0003 FC1 SQVQLVESGGGLVQAGGSLRLSCAASGFPVSWWKMYWYRQAPGKEREWVAAIWSEGWWTKYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDYGGENANYDYWGQGTQVTVS\nGSDGTTEDAWQEGGR\n"
## [4] ">NB0004 FC1 SQVQLVESGGGLVQAGGSLRLSCAASGFPVEWSWMYWYRQAPGKEREWVAAIYSQGRGTEYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDYGWWYGDYDYWGQGTQVTVS\nGSEEAADPAWR\n"
## [5] ">NB0005 FC1 SQVQLVESGGGLVQAGGSLRLSCAASGFPVEAHRMYWYRQAPGKEREWVAAISSKGQQTWYADSVKGRFTISRDNAKNTVYLQMNSLEPEDTAVYYCNVKDYGWYYGDYDYWGQGTQVTVS\nGSEEAEATWWR\n"
## [6] ">NB0006 FC2 SQVQLVESGGGLVQAGGSLRLSCAASGFPVEENFMYWYRQAPGKEREWVAAIYSHGYETEYADSVKGRFTISRDNAKNTVYLQMNSLKPEDTAVYYCNVKDQGYWWWEYDYWGQGTQVTVS\nGSGLPATPAWLRGSTDAEEGVWLTVR\n"
To analyze the expressed flycodes mass spectrometry is used.
the FASTA file containing the nanobody - flycode linkage can
be written to a file using functions such as cat
.
The exec directory provides alternative shell scripts using command line GNU tools and AWK.
Here is the output of the sessionInfo()
command.
## R version 3.6.0 (2019-04-26)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.2 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.9-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.9-bioc/R/lib/libRlapack.so
##
## locale:
## [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
## [3] LC_TIME=en_US.UTF-8 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] stats4 parallel stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] scales_1.0.0 ggplot2_3.1.1
## [3] NestLink_1.0.0 ShortRead_1.42.0
## [5] GenomicAlignments_1.20.0 SummarizedExperiment_1.14.0
## [7] DelayedArray_0.10.0 matrixStats_0.54.0
## [9] Biobase_2.44.0 Rsamtools_2.0.0
## [11] GenomicRanges_1.36.0 GenomeInfoDb_1.20.0
## [13] BiocParallel_1.18.0 protViz_0.4.0
## [15] gplots_3.0.1.1 Biostrings_2.52.0
## [17] XVector_0.24.0 IRanges_2.18.0
## [19] S4Vectors_0.22.0 ExperimentHub_1.10.0
## [21] AnnotationHub_2.16.0 BiocFileCache_1.8.0
## [23] dbplyr_1.4.0 BiocGenerics_0.30.0
## [25] BiocStyle_2.12.0
##
## loaded via a namespace (and not attached):
## [1] httr_1.4.0 bit64_0.9-7
## [3] gtools_3.8.1 shiny_1.3.2
## [5] assertthat_0.2.1 interactiveDisplayBase_1.22.0
## [7] BiocManager_1.30.4 latticeExtra_0.6-28
## [9] blob_1.1.1 GenomeInfoDbData_1.2.1
## [11] yaml_2.2.0 pillar_1.3.1
## [13] RSQLite_2.1.1 lattice_0.20-38
## [15] glue_1.3.1 digest_0.6.18
## [17] RColorBrewer_1.1-2 promises_1.0.1
## [19] colorspace_1.4-1 plyr_1.8.4
## [21] htmltools_0.3.6 httpuv_1.5.1
## [23] Matrix_1.2-17 pkgconfig_2.0.2
## [25] bookdown_0.9 zlibbioc_1.30.0
## [27] purrr_0.3.2 xtable_1.8-4
## [29] gdata_2.18.0 later_0.8.0
## [31] tibble_2.1.1 withr_2.1.2
## [33] lazyeval_0.2.2 magrittr_1.5
## [35] crayon_1.3.4 mime_0.6
## [37] memoise_1.1.0 evaluate_0.13
## [39] hwriter_1.3.2 tools_3.6.0
## [41] stringr_1.4.0 munsell_0.5.0
## [43] AnnotationDbi_1.46.0 compiler_3.6.0
## [45] caTools_1.17.1.2 rlang_0.3.4
## [47] grid_3.6.0 RCurl_1.95-4.12
## [49] rappdirs_0.3.1 labeling_0.3
## [51] bitops_1.0-6 rmarkdown_1.12
## [53] gtable_0.3.0 codetools_0.2-16
## [55] DBI_1.0.0 curl_3.3
## [57] R6_2.4.0 knitr_1.22
## [59] dplyr_0.8.0.1 bit_1.1-14
## [61] KernSmooth_2.23-15 stringi_1.4.3
## [63] Rcpp_1.0.1 tidyselect_0.2.5
## [65] xfun_0.6
Egloff, Pascal, Iwan Zimmermann, Fabian M. Arnold, Cedric A.J. Hutter, Damien Damien Morger, Lennart Opitz, Lucy Poveda, et al. 2018. “Engineered Peptide Barcodes for In-Depth Analyses of Binding Protein Ensembles.” bioRxiv. Cold Spring Harbor Laboratory. https://doi.org/10.1101/287813.