-
Notifications
You must be signed in to change notification settings - Fork 4
/
Copy pathget_signalp5.R
545 lines (500 loc) · 19 KB
/
get_signalp5.R
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
#' Query SignalP 5.0 web server.
#'
#' The SignalP 5.0 server predicts the presence of signal peptides and the location of their cleavage sites in proteins from Archaea, Gram-positive Bacteria, Gram-negative Bacteria and Eukarya.
#'
#' @aliases get_signalp5 get_signalp5.default get_signalp5.character get_signalp5.data.frame get_signalp5.list get_signalp5.AAStringSet
#' @param data A data frame with protein amino acid sequences as strings in one column and corresponding id's in another. Alternatively a path to a .fasta file with protein sequences. Alternatively a list with elements of class \code{\link[seqinr]{SeqFastaAA}} resulting from \code{\link[seqinr]{read.fasta}} call. Alternatively an \code{\link[Biostrings]{AAStringSet}} object. Should be left blank if vectors are provided to sequence and id arguments.
#' @param sequence A vector of strings representing protein amino acid sequences, or the appropriate column name if a data.frame is supplied to data argument. If .fasta file path, or list with elements of class "SeqFastaAA" provided to data, this should be left blank.
#' @param id A vector of strings representing protein identifiers, or the appropriate column name if a data.frame is supplied to data argument. If .fasta file path, or list with elements of class "SeqFastaAA" provided to data, this should be left blank.
#' @param org_type One of "euk", "gram-", "gram+" or "archea". Default is "euk". Are the protein sequences from Eukarya, Gram-negative Bacteria, Gram-positive Bacteria or Archaea.
#' @param splitter An integer indicating the number of sequences to be in each .fasta file that is to be sent to the server. Default is 2500. Change only in case of a server side error. Accepted values are in range of 1 to 5000.
#' @param attempts Integer, number of attempts if server unresponsive, at default set to 2.
#' @param progress Boolean, whether to show messages of the job id for each batch. Default is FALSE.
#' @param ... currently no additional arguments are accepted apart the ones documented bellow.
#'
#' @return if org_type is "euk" a data frame with columns:
#' \describe{
#' \item{id}{Character, as from input}
#' \item{Prediction}{Integer, The type of signal peptide predicted: Sec/SPI, Tat/SPI, Sec/SPII or Other if no signal peptide predicted}
#' \item{SP.Sec.SPI}{Numeric, marginal probability that the protein contains a Sec N-terminal signal peptide (Sec/SPI).}
#' \item{Other}{Numeric, the probability that the sequence does not have any kind of signal peptid.}
#' \item{CS_pos}{Character, cleavage site position.}
#' \item{Pr}{Numeric, probability of the predicted cleavage site position.}
#' \item{cleave.site}{Numeric,llocal amino acid sequence arround the predicted cleavage site.}
#' \item{is.signalp}{Logical, did SignalP5 predict the presence of a signal peptide.}
#' \item{sp.length}{Integer, length of the predicted signal peptide.}
#' }
#'
#' if org_type is one of "gram-", "gram+" or "archea" the returned data frame will have two additional columns between `SP.Sec.SPI` and `Other`:
#'
#' \describe{
#' \item{TAT.Tat.SPI}{Numeric, marginal probability that the protein contains a Tat N-terminal signal peptide (Tat/SPI).}
#' \item{LIPO.Sec.SPII}{Numeric, marginal probability that the protein contains a Lipoprotein N-terminal signal peptide (Sec/SPII).}
#' }
#'
#'
#' @note This function creates temporary files in the working directory. If something goes wrong during communication with the server and progress was set to TRUE, predictions can be obtained using `file.path("http://www.cbs.dtu.dk/services/SignalP-5.0/tmp", jobid, "output_protein_type.txt")` eg `read.delim(file.path(...), header = TRUE, skip = 1)`.
#'
#'
#' @source \url{https://services.healthtech.dtu.dk/service.php?SignalP-5.0}
#' @references Almagro Armenteros JJ, Tsirigos KD, Sønderby CK, Petersen TN, Winther O, Brunak S, von Heijne G, Nielsen H. (2019) SignalP 5.0 improves signal peptide predictions using deep neural networks. Nature Biotechnology, 37:420-423, doi:10.1038/s41587-019-0036-z
#
#' @seealso \code{\link[ragp]{get_signalp}} \code{\link[ragp]{get_targetp}}
#'
#' @examples
#'
#' library(ragp)
#' sp5_pred <- get_signalp5(data = at_nsp[1:10,],
#' sequence,
#' Transcript.id)
#' sp5_pred
#'
#' @import seqinr
#' @import httr
#' @import xml2
#' @export
get_signalp5 <- function(data, ...){
if (missing(data) || is.null(data)) get_signalp5.default(...)
else UseMethod("get_signalp5")
}
#' @rdname get_signalp5
#' @method get_signalp5 character
#' @export
get_signalp5.character <- function(data,
org_type = c("euk", "gram-", "gram+", "archea"),
splitter = 2500L,
attempts = 2,
progress = FALSE,
...){
if (missing(splitter)) {
splitter <- 2500L
}
if (length(splitter) > 1){
splitter <- 2500L
warning("splitter should be of length 1, setting to default: splitter = 2500",
call. = FALSE)
}
if (!is.numeric(splitter)){
splitter <- as.numeric(splitter)
warning("splitter is not numeric, converting using 'as.numeric'",
call. = FALSE)
}
if (is.na(splitter)){
splitter <- 2500L
warning("splitter was set to NA, setting to default: splitter = 2500",
call. = FALSE)
}
if (is.numeric(splitter)) {
splitter <- floor(splitter)
}
if (!(splitter %in% 1:5000)) {
splitter <- 2500L
warning("Illegal splitter input, splitter will be set to 2500",
call. = FALSE)
}
if (length(attempts) > 1){
attempts <- 2L
warning("attempts should be of length 1, setting to default: attempts = 2",
call. = FALSE)
}
if (!is.numeric(attempts)){
attempts<- as.numeric(attempts)
warning("attempts is not numeric, converting using 'as.numeric'",
call. = FALSE)
}
if (is.na(attempts)){
attempts <- 2L
warning("attempts was set to NA, setting to default: attempts = 2",
call. = FALSE)
}
if (is.numeric(attempts)) {
attempts <- floor(attempts)
}
if (attempts < 1) {
attempts <- 2L
warning("attempts was set to less then 1, setting to default: attempts = 2",
call. = FALSE)
}
if (missing(progress)) {
progress <- FALSE
}
if (length(progress) > 1){
progress <- FALSE
warning("progress should be of length 1, setting to default: progress = FALSE",
call. = FALSE)
}
if (!is.logical(progress)){
progress <- as.logical(progress)
warning("progress is not logical, converting using 'as.logical'",
call. = FALSE)
}
if (is.na(progress)){
progress <- FALSE
warning("progress was set to NA, setting to default: progress = FALSE",
call. = FALSE)
}
if (missing(org_type)) {
org_type <- "euk"
}
if (!org_type %in% c("euk", "gram-", "gram+", "archea")) {
stop("org_type should be one of: 'euk', 'gram-', 'gram+', 'archea'",
call. = FALSE)
}
if (length(org_type) > 1){
stop("org_type should be one of: 'euk', 'gram-', 'gram+', 'archea'",
call. = FALSE)
}
if(length(data) > 1){
stop("one fasta file per function call can be supplied",
call. = FALSE)
}
if (file.exists(data)){
file_name <- data
} else {
stop("cannot find file in the specified path",
call. = FALSE)
}
organism <- c("Eukarya",
"Gram-negative bacteria",
"Gram-positive bacteria",
"Archaea")
names(organism) <- c("euk", "gram-", "gram+", "archea")
organism <- organism[names(organism) == org_type]
url <- "https://services.healthtech.dtu.dk/cgi-bin/webface2.fcgi"
cfg_file <- "/var/www/html/services/SignalP-5.0/webface.cf"
file_list <- ragp::split_fasta(path_in = file_name,
path_out = "tmp_signalp5_",
num_seq = splitter)
if(grepl("temp_", file_name)){
unlink(file_name)
}
output <- vector("list", length(file_list))
for(k in seq_along(file_list)){
x <- file_list[k]
file_up <- httr::upload_file(x)
res <- httr::POST(url = url,
encode = "multipart",
body = list(configfile = cfg_file,
uploadfile = file_up,
organism = organism,
format = "short"))
if(!grepl("jobid=", res$url)){
stop("something went wrong on server side")
}
res <- sub("https://services.healthtech.dtu.dk/cgi-bin/webface2.cgi?jobid=",
"",
res$url,
fixed = TRUE)
res <- sub("&wait=20",
"",
res,
fixed = TRUE)
jobid <- res
if(progress){
message(paste("batch", k, "jobid is:", jobid))
}
time1 <- Sys.time()
repeat {
res2 <- httr::GET(url = file.path("https://services.healthtech.dtu.dk/services/SignalP-5.0/tmp",
jobid,
"output_protein_type.txt"))
code <- res2$status_code
Sys.sleep(5)
if (code == 200L) {
res2 <- read.delim(
file.path("https://services.healthtech.dtu.dk/services/SignalP-5.0/tmp",
jobid,
"output_protein_type.txt"),
header = TRUE,
skip = 1)
break
}
time2 <- Sys.time()
max.time <- as.difftime(pmax(50, splitter),
units = "secs")
if ((time2 - time1) > max.time) {
code <- 404L
if(progress) message("file",
x,
"took longer then expected")
break
}
}
if (code == 404L) {
tms <- 0
while(tms < attempts && code == 404L){
if(progress) message(
"reattempting file",
x)
file_up <- httr::upload_file(x)
res <- httr::POST(url = url,
encode = "multipart",
body = list(configfile = cfg_file,
uploadfile = file_up,
organism = organism,
format = "short"))
if(!grepl("jobid=", res$url)){
stop("something went wrong on server side")
}
res <- sub("https://services.healthtech.dtu.dk/cgi-bin/webface2.cgi?jobid=",
"",
res$url,
fixed = TRUE)
res <- sub("&wait=20",
"",
res,
fixed = TRUE)
jobid <- res
if(progress){
message(paste("batch", k, "reatempt jobid is:", jobid))
}
time1 <- Sys.time()
repeat {
res2 <- httr::GET(url = file.path("https://services.healthtech.dtu.dk/services/SignalP-5.0/tmp",
jobid,
"output_protein_type.txt"))
code <- res2$status_code
Sys.sleep(5)
if (code == 200L) {
res2 <- read.delim(
file.path("https://services.healthtech.dtu.dk/services/SignalP-5.0/tmp",
jobid,
"output_protein_type.txt"),
header = TRUE,
skip = 1,
stringsAsFactors = FALSE)
break
}
time2 <- Sys.time()
max.time <- as.difftime(pmax(100, splitter*1.5),
units = "secs")
if ((time2 - time1) > max.time) {
code <- 404L
if(progress) message(
"file",
x,
"took longer then expected")
break
}
}
tms <- tms + 1
}
}
if (code == 404L){
output <- do.call(rbind,
output)
output$is.signalp <- output$is.sp == "Y"
warning(
"maximum attempts reached at",
x,
"returning finished queries",
call. = FALSE)
return(output)
}
unlink(x)
res2 <- data.frame(res2[,colnames(res2) != "CS.Position"],
CS_pos = gsub("CS pos: (\\d+)-(\\d+)\\..*$", "\\1-\\2", res2[,"CS.Position"]),
Pr = gsub(".*Pr: ", "", res2[,"CS.Position"]),
substr = gsub("^.*?([ARNDCEQGHILKMFPSTWYV]{3}-[ARNDCEQGHILKMFPSTWYV]{2}).*$" ,
"\\1", res2[,"CS.Position"]),
stringsAsFactors = FALSE)
if(org_type == "euk"){
colnames(res2) <- c("id",
"Prediction",
"SP.Sec.SPI",
"Other",
"CS_pos",
"Pr",
"cleave.site")
res2$SP.Sec.SPI <- as.numeric(res2$SP.Sec.SPI)
}
if(org_type != "euk"){
colnames(res2) <- c("id",
"Prediction",
"SP.Sec.SPI",
"TAT.Tat.SPI",
"LIPO.Sec.SPII",
"Other",
"CS_pos",
"Pr",
"cleave.site")
res2$SP.Sec.SPI <- as.numeric(res2$SP.Sec.SPI)
res2$TAT.Tat.SPI <- as.numeric(res2$TAT.Tat.SPI)
res2$LIPO.Sec.SPII <- as.numeric(res2$LIPO.Sec.SPII)
}
res2$Other <- as.numeric(res2$Other)
res2$Pr <- as.numeric(res2$Pr)
output[[k]] <- res2
}
output <- do.call(rbind,
output)
output$is.signalp <- !is.na(output$Pr)
output$sp.length <- as.integer(gsub("(\\d+)-\\d+", "\\1", output[,"CS_pos"]))
return(output)
}
#' @rdname get_signalp5
#' @method get_signalp5 data.frame
#' @export
get_signalp5.data.frame <- function(data,
sequence,
id,
...){
if(missing(sequence)){
stop("the column name with the sequences must be specified",
call. = FALSE)
}
if(missing(id)){
stop("the column name with the sequence id's must be specified",
call. = FALSE)
}
id <- as.character(substitute(id))
sequence <- as.character(substitute(sequence))
if (length(id) != 1L){
stop("only one column name for 'id' must be specifed",
call. = FALSE)
}
if (length(sequence) != 1L){
stop("only one column name for 'sequence' must be specifed",
call. = FALSE)
}
id <- if(id %in% colnames(data)){
data[[id]]
} else {
stop("specified 'id' not found in data",
call. = FALSE)
}
id <- as.character(id)
sequence <- if(sequence %in% colnames(data)){
data[[sequence]]
} else {
stop("specified 'sequence' not found in data",
call. = FALSE)
}
sequence <- toupper(as.character(sequence))
sequence <- sub("\\*$",
"",
sequence)
aa_regex <- "[^ARNDCQEGHILKMFPSTWYV]"
if (any(grepl(aa_regex, sequence))){
warning(paste("sequences: ",
paste(id[grepl(aa_regex,
sequence)],
collapse = ", "),
" contain symbols not corresponding to amino acids",
sep = ""),
call. = FALSE)
}
file_name <- paste("temp_",
gsub("^X",
"",
make.names(Sys.time())),
".fasta",
sep = "")
seqinr::write.fasta(sequence = strsplit(sequence, ""),
name = id,
file = file_name)
res <- get_signalp5.character(data = file_name, ...)
return(res)
}
#' @rdname get_signalp5
#' @method get_signalp5 list
#' @export
get_signalp5.list <- function(data,
...){
if(class(data[[1]]) == "SeqFastaAA"){
dat <- lapply(data,
paste0,
collapse ="")
id <- names(dat)
sequence <- toupper(as.character(unlist(dat)))
sequence <- sub("\\*$",
"",
sequence)
aa_regex <- "[^ARNDCQEGHILKMFPSTWYV]"
if (any(grepl(aa_regex, sequence))){
warning(paste("sequences: ",
paste(id[grepl(aa_regex,
sequence)],
collapse = ", "),
" contain symbols not corresponding to amino acids",
sep = ""),
call. = FALSE)
}
file_name <- paste("temp_",
gsub("^X",
"",
make.names(Sys.time())),
".fasta",
sep = "")
seqinr::write.fasta(sequence = strsplit(sequence, ""),
name = id,
file = file_name)
} else {
stop("only lists containing objects of class SeqFastaAA are supported")
}
res <- get_signalp5.character(data = file_name, ...)
return(res)
}
#' @rdname get_signalp5
#' @method get_signalp5 default
#' @export
get_signalp5.default <- function(data = NULL,
sequence,
id,
...){
if (missing(sequence)){
stop("protein sequence must be provided to obtain predictions",
call. = FALSE)
}
if (missing(id)){
stop("protein id must be provided to obtain predictions",
call. = FALSE)
}
id <- as.character(id)
sequence <- toupper(as.character(sequence))
if (length(sequence) != length(id)){
stop("id and sequence vectors are not of same length",
call. = FALSE)
}
sequence <- sub("\\*$",
"",
sequence)
aa_regex <- "[^ARNDCQEGHILKMFPSTWYV]"
if (any(grepl(aa_regex, sequence))){
warning(paste("sequences: ",
paste(id[grepl(aa_regex,
sequence)],
collapse = ", "),
" contain symbols not corresponding to amino acids",
sep = ""),
call. = FALSE)
}
file_name <- paste("temp_",
gsub("^X",
"",
make.names(Sys.time())),
".fasta",
sep = "")
seqinr::write.fasta(sequence = strsplit(sequence, ""),
name = id,
file = file_name)
res <- get_signalp5.character(data = file_name, ...)
return(res)
}
#' @rdname get_signalp5
#' @method get_signalp5 AAStringSet
#' @export
get_signalp5.AAStringSet <- function(data,
...){
sequence <- as.character(data)
id <- names(sequence)
sequence <- unname(sequence)
sequence <- toupper(sequence)
sequence <- sub("\\*$",
"",
sequence)
res <- get_signalp5.default(sequence = sequence,
id = id,
...)
return(res)
}