NestLink 1.14.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(): 2022-10-24query(eh, "NestLink")## ExperimentHub with 8 records
## # snapshotDate(): 2022-10-24
## # $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## loading from cachescratchFolder <- 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## loading from cacheknownNB_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 analysisThe following steps are included:
system.time(NB2FC <- runNGSAnalysis(file = expFile[1], param))##    user  system elapsed 
##   2.503   0.256   2.760head(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        
## 2head(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 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.16-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.16-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] stats4    stats     graphics  grDevices utils     datasets  methods  
## [8] base     
## 
## other attached packages:
##  [1] scales_1.2.1                ggplot2_3.3.6              
##  [3] NestLink_1.14.0             ShortRead_1.56.0           
##  [5] GenomicAlignments_1.34.0    SummarizedExperiment_1.28.0
##  [7] Biobase_2.58.0              MatrixGenerics_1.10.0      
##  [9] matrixStats_0.62.0          Rsamtools_2.14.0           
## [11] GenomicRanges_1.50.0        BiocParallel_1.32.0        
## [13] protViz_0.7.3               gplots_3.1.3               
## [15] Biostrings_2.66.0           GenomeInfoDb_1.34.0        
## [17] XVector_0.38.0              IRanges_2.32.0             
## [19] S4Vectors_0.36.0            ExperimentHub_2.6.0        
## [21] AnnotationHub_3.6.0         BiocFileCache_2.6.0        
## [23] dbplyr_2.2.1                BiocGenerics_0.44.0        
## [25] BiocStyle_2.26.0           
## 
## loaded via a namespace (and not attached):
##  [1] nlme_3.1-160                  bitops_1.0-7                 
##  [3] bit64_4.0.5                   RColorBrewer_1.1-3           
##  [5] filelock_1.0.2                httr_1.4.4                   
##  [7] tools_4.2.1                   bslib_0.4.0                  
##  [9] utf8_1.2.2                    R6_2.5.1                     
## [11] KernSmooth_2.23-20            mgcv_1.8-41                  
## [13] colorspace_2.0-3              DBI_1.1.3                    
## [15] withr_2.5.0                   tidyselect_1.2.0             
## [17] bit_4.0.4                     curl_4.3.3                   
## [19] compiler_4.2.1                cli_3.4.1                    
## [21] DelayedArray_0.24.0           labeling_0.4.2               
## [23] bookdown_0.29                 sass_0.4.2                   
## [25] caTools_1.18.2                rappdirs_0.3.3               
## [27] stringr_1.4.1                 digest_0.6.30                
## [29] rmarkdown_2.17                jpeg_0.1-9                   
## [31] pkgconfig_2.0.3               htmltools_0.5.3              
## [33] highr_0.9                     fastmap_1.1.0                
## [35] rlang_1.0.6                   RSQLite_2.2.18               
## [37] shiny_1.7.3                   farver_2.1.1                 
## [39] jquerylib_0.1.4               generics_0.1.3               
## [41] hwriter_1.3.2.1               jsonlite_1.8.3               
## [43] gtools_3.9.3                  dplyr_1.0.10                 
## [45] RCurl_1.98-1.9                magrittr_2.0.3               
## [47] GenomeInfoDbData_1.2.9        interp_1.1-3                 
## [49] Matrix_1.5-1                  munsell_0.5.0                
## [51] Rcpp_1.0.9                    fansi_1.0.3                  
## [53] lifecycle_1.0.3               stringi_1.7.8                
## [55] yaml_2.3.6                    zlibbioc_1.44.0              
## [57] grid_4.2.1                    blob_1.2.3                   
## [59] parallel_4.2.1                promises_1.2.0.1             
## [61] crayon_1.5.2                  deldir_1.0-6                 
## [63] lattice_0.20-45               splines_4.2.1                
## [65] KEGGREST_1.38.0               magick_2.7.3                 
## [67] knitr_1.40                    pillar_1.8.1                 
## [69] codetools_0.2-18              glue_1.6.2                   
## [71] BiocVersion_3.16.0            evaluate_0.17                
## [73] latticeExtra_0.6-30           BiocManager_1.30.19          
## [75] png_0.1-7                     vctrs_0.5.0                  
## [77] httpuv_1.6.6                  purrr_0.3.5                  
## [79] gtable_0.3.1                  assertthat_0.2.1             
## [81] cachem_1.0.6                  xfun_0.34                    
## [83] mime_0.12                     xtable_1.8-4                 
## [85] later_1.3.0                   tibble_3.1.8                 
## [87] AnnotationDbi_1.60.0          memoise_2.0.1                
## [89] ellipsis_0.3.2                interactiveDisplayBase_1.36.0Egloff, 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. https://doi.org/10.1101/287813.