NestLink 1.4.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(): 2020-04-27
query(eh, "NestLink")
## ExperimentHub with 8 records
## # snapshotDate(): 2020-04-27
## # $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 cache
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
## loading from cache
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
## 2.842 0.092 2.945
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 4.0.0 (2020-04-24)
## Platform: x86_64-pc-linux-gnu (64-bit)
## Running under: Ubuntu 18.04.4 LTS
##
## Matrix products: default
## BLAS: /home/biocbuild/bbs-3.11-bioc/R/lib/libRblas.so
## LAPACK: /home/biocbuild/bbs-3.11-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.1.0 ggplot2_3.3.0
## [3] NestLink_1.4.0 ShortRead_1.46.0
## [5] GenomicAlignments_1.24.0 SummarizedExperiment_1.18.1
## [7] DelayedArray_0.14.0 matrixStats_0.56.0
## [9] Biobase_2.48.0 Rsamtools_2.4.0
## [11] GenomicRanges_1.40.0 GenomeInfoDb_1.24.0
## [13] BiocParallel_1.22.0 protViz_0.6.8
## [15] gplots_3.0.3 Biostrings_2.56.0
## [17] XVector_0.28.0 IRanges_2.22.1
## [19] S4Vectors_0.26.0 ExperimentHub_1.14.0
## [21] AnnotationHub_2.20.0 BiocFileCache_1.12.0
## [23] dbplyr_1.4.3 BiocGenerics_0.34.0
## [25] BiocStyle_2.16.0
##
## loaded via a namespace (and not attached):
## [1] nlme_3.1-147 bitops_1.0-6
## [3] bit64_0.9-7 RColorBrewer_1.1-2
## [5] httr_1.4.1 tools_4.0.0
## [7] R6_2.4.1 KernSmooth_2.23-17
## [9] mgcv_1.8-31 DBI_1.1.0
## [11] colorspace_1.4-1 withr_2.2.0
## [13] tidyselect_1.0.0 bit_1.1-15.2
## [15] curl_4.3 compiler_4.0.0
## [17] labeling_0.3 bookdown_0.18
## [19] caTools_1.18.0 rappdirs_0.3.1
## [21] stringr_1.4.0 digest_0.6.25
## [23] rmarkdown_2.1 jpeg_0.1-8.1
## [25] pkgconfig_2.0.3 htmltools_0.4.0
## [27] fastmap_1.0.1 rlang_0.4.6
## [29] RSQLite_2.2.0 shiny_1.4.0.2
## [31] farver_2.0.3 hwriter_1.3.2
## [33] gtools_3.8.2 dplyr_0.8.5
## [35] RCurl_1.98-1.2 magrittr_1.5
## [37] GenomeInfoDbData_1.2.3 Matrix_1.2-18
## [39] Rcpp_1.0.4.6 munsell_0.5.0
## [41] lifecycle_0.2.0 stringi_1.4.6
## [43] yaml_2.2.1 zlibbioc_1.34.0
## [45] grid_4.0.0 blob_1.2.1
## [47] gdata_2.18.0 promises_1.1.0
## [49] crayon_1.3.4 lattice_0.20-41
## [51] splines_4.0.0 magick_2.3
## [53] knitr_1.28 pillar_1.4.4
## [55] codetools_0.2-16 glue_1.4.0
## [57] BiocVersion_3.11.1 evaluate_0.14
## [59] latticeExtra_0.6-29 BiocManager_1.30.10
## [61] vctrs_0.2.4 png_0.1-7
## [63] httpuv_1.5.2 gtable_0.3.0
## [65] purrr_0.3.4 assertthat_0.2.1
## [67] xfun_0.13 mime_0.9
## [69] xtable_1.8-4 later_1.0.0
## [71] tibble_3.0.1 AnnotationDbi_1.50.0
## [73] memoise_1.1.0 ellipsis_0.3.0
## [75] interactiveDisplayBase_1.26.0
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.