Browse code

Make writeMSData a method

- Import writeMSData generic from ProtGenerics.
- Convert writeMSData from a function to a method.

jotsetung authored on 20/09/2017 10:19:40
Showing 6 changed files

... ...
@@ -18,7 +18,7 @@ Description: mzR provides a unified API to the common file formats and
18 18
 License: Artistic-2.0
19 19
 LazyLoad: yes
20 20
 Depends: Rcpp (>= 0.10.1), methods, utils
21
-Imports: Biobase, BiocGenerics (>= 0.13.6), ProtGenerics
21
+Imports: Biobase, BiocGenerics (>= 0.13.6), ProtGenerics (>= 1.9.1)
22 22
 Suggests: msdata (>= 0.15.1), RUnit, mzID, BiocStyle (>= 2.5.19), knitr, XML
23 23
 VignetteBuilder: knitr
24 24
 LinkingTo: Rcpp, zlibbioc
... ...
@@ -13,7 +13,6 @@ export(openMSfile,
13 13
        openIDfile,
14 14
        pwiz.version,
15 15
        nChrom,
16
-       writeMSData,
17 16
        copyWriteMSData)
18 17
 
19 18
 exportMethods(close,
... ...
@@ -47,7 +46,8 @@ exportMethods(close,
47 46
               para,
48 47
               tic,
49 48
               chromatogram,
50
-              chromatograms)
49
+              chromatograms,
50
+              writeMSData)
51 51
 
52 52
 exportClasses("mzR",
53 53
               "mzRramp",
... ...
@@ -52,45 +52,43 @@ openIDfile <- function(filename, verbose = FALSE) {
52 52
                fileName=filename))
53 53
 }
54 54
 
55
-writeMSData <- function(filename, header, data, backend = "pwiz",
56
-                        outformat = "mzml",
57
-                        rtime_seconds = TRUE,
58
-                        software_processing) {
59
-    backend <- match.arg(backend)
60
-    ## supp_formats <- c("mzml", "mgf", "mzxml")
61
-    supp_formats <- c("mzml", "mzxml")
62
-    outformat <- match.arg(tolower(outformat), supp_formats)
63
-    if (missing(filename))
64
-        stop("'filename' is a required parameter")
65
-    if (missing(header) | missing(data))
66
-        stop("'header' and 'data' are required")
67
-    ## Other checks:
68
-    header <- .validateHeader(header)
69
-    if (is(header, "character"))
70
-        stop("Error checking parameter 'header': ", header)
71
-    is_ok <- .validSpectrumList(data)
72
-    if (is(is_ok, "character"))
73
-        stop("Error checking parameter 'data'. First error was: ", is_ok)
74
-    ## Check software_processing:
75
-    software_processing <- .check_software_processing(software_processing)
76
-    ## Add mzR processing:
77
-    mzR <- c("mzR", paste(packageVersion("mzR"), collapse = "."), "MS:-1")
78
-    if (outformat == "mzml")
79
-        mzR <- c(mzR, "MS:1000544")
80
-    if (outformat == "mzxml")
81
-        mzR <- c(mzR, "MS:1000545")
82
-    software_processing <- c(software_processing, list(mzR))
83
-    if (backend == "pwiz") {
84
-        if (outformat == "mzxml" & any(header$injectionTime > 0))
85
-            warning("mzXML export does not support writing ion injection time")
86
-        pwizModule <- new(Pwiz)
87
-        pwizModule$writeSpectrumList(filename, outformat,
88
-                                     header, data, rtime_seconds,
89
-                                     software_processing)
90
-    }
91
-}
55
+setMethod("writeMSData", signature(object = "list", file = "character"),
56
+          function(object, file, header, backend = "pwiz",
57
+                   outformat = "mzml", rtime_seconds = TRUE,
58
+                   software_processing) {
59
+              backend <- match.arg(backend)
60
+              ## supp_formats <- c("mzml", "mgf", "mzxml")
61
+              supp_formats <- c("mzml", "mzxml")
62
+              outformat <- match.arg(tolower(outformat), supp_formats)
63
+              if (missing(header))
64
+                  stop("'header' is required")
65
+              ## Other checks:
66
+              header <- .validateHeader(header)
67
+              if (is(header, "character"))
68
+                  stop("Error checking parameter 'header': ", header)
69
+              is_ok <- .validSpectrumList(object)
70
+              if (is(is_ok, "character"))
71
+                  stop("Error checking parameter 'object'. First error was: ", is_ok)
72
+              ## Check software_processing:
73
+              software_processing <- .check_software_processing(software_processing)
74
+              ## Add mzR processing:
75
+              mzR <- c("mzR", paste(packageVersion("mzR"), collapse = "."), "MS:-1")
76
+              if (outformat == "mzml")
77
+                  mzR <- c(mzR, "MS:1000544")
78
+              if (outformat == "mzxml")
79
+                  mzR <- c(mzR, "MS:1000545")
80
+              software_processing <- c(software_processing, list(mzR))
81
+              if (backend == "pwiz") {
82
+                  if (outformat == "mzxml" & any(header$injectionTime > 0))
83
+                      warning("mzXML export does not support writing ion injection time")
84
+                  pwizModule <- new(Pwiz)
85
+                  pwizModule$writeSpectrumList(file, outformat,
86
+                                               header, object, rtime_seconds,
87
+                                               software_processing)
88
+              }
89
+          })
92 90
 
93
-copyWriteMSData <- function(filename, original_file, header, data,
91
+copyWriteMSData <- function(object, file, original_file, header,
94 92
                             backend = "pwiz",
95 93
                             outformat = "mzml",
96 94
                             rtime_seconds = TRUE,
... ...
@@ -99,21 +97,21 @@ copyWriteMSData <- function(filename, original_file, header, data,
99 97
     ## supp_formats <- c("mzml", "mgf", "mzxml")
100 98
     supp_formats <- c("mzml", "mzxml")
101 99
     outformat <- match.arg(tolower(outformat), supp_formats)
102
-    if (missing(filename))
103
-        stop("'filename' is a required parameter")
100
+    if (missing(file))
101
+        stop("'file' is a required parameter")
104 102
     if (missing(original_file))
105 103
         stop("'original_file' is a required parameter")
106
-    if (missing(header) | missing(data))
107
-        stop("'header' and 'data' are required")
104
+    if (missing(header) | missing(object))
105
+        stop("'header' and 'object' are required")
108 106
     if (!file.exists(original_file))
109 107
         stop("Original file ", original_file, " not found")
110 108
     ## Other checks:
111 109
     header <- .validateHeader(header)
112 110
     if (is(header, "character"))
113 111
         stop("Error checking parameter 'header': ", header)
114
-    is_ok <- .validSpectrumList(data)
112
+    is_ok <- .validSpectrumList(object)
115 113
     if (is(is_ok, "character"))
116
-        stop("Error checking parameter 'data'. First error was: ", is_ok)
114
+        stop("Error checking parameter 'object'. First error was: ", is_ok)
117 115
     ## Check software_processing:
118 116
     software_processing <- .check_software_processing(software_processing)
119 117
     ## Add mzR processing:
... ...
@@ -129,8 +127,8 @@ copyWriteMSData <- function(filename, original_file, header, data,
129 127
             header$injectionTime = 0
130 128
         }
131 129
         pwizModule <- new(Pwiz)
132
-        pwizModule$copyWriteMSfile(filename, outformat, original_file,
133
-                                   header, data, rtime_seconds,
130
+        pwizModule$copyWriteMSfile(file, outformat, original_file,
131
+                                   header, object, rtime_seconds,
134 132
                                    software_processing)
135 133
     }
136 134
 }
... ...
@@ -8,7 +8,7 @@ test_copyWriteMSData <- function() {
8 8
 
9 9
     ## INPUT: mzXML
10 10
     orig_file <- system.file("threonine", "threonine_i2_e35_pH_tree.mzXML",
11
-                        package = "msdata")
11
+                             package = "msdata")
12 12
     mzxml <- openMSfile(orig_file, backend = "pwiz")
13 13
     pks <- peaks(mzxml)
14 14
     hdr <- header(mzxml)
... ...
@@ -17,8 +17,8 @@ test_copyWriteMSData <- function() {
17 17
 
18 18
     ## OUTPUT: mzML
19 19
     fnew <- paste0(test_folder, "test_copyWrite.mzML")
20
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
21
-                          header = hdr, data = pks, backend = "pwiz")
20
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
21
+                         header = hdr, object = pks, backend = "pwiz")
22 22
     ## Check content is same
23 23
     mzml_new <- openMSfile(fnew, backend = "pwiz")
24 24
     pks_new <- peaks(mzml_new)
... ...
@@ -31,9 +31,9 @@ test_copyWriteMSData <- function() {
31 31
     
32 32
     ## OUTPUT: mzXML
33 33
     fnew <- paste0(test_folder, "test_copyWrite.mzXML")
34
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
35
-                          header = hdr, data = pks, backend = "pwiz",
36
-                          outformat = "mzxml")
34
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
35
+                         header = hdr, object = pks, backend = "pwiz",
36
+                         outformat = "mzxml")
37 37
     ## Check content is same
38 38
     mzml_new <- openMSfile(fnew, backend = "pwiz")
39 39
     pks_new <- peaks(mzml_new)
... ...
@@ -64,15 +64,15 @@ test_copyWriteMSData <- function() {
64 64
     pks_sub <- pks[c(1, 3, 5)]
65 65
     fnew <- paste0(test_folder, "test_copyWrite.mzML")
66 66
     ## index is not OK after subsetting
67
-    checkException(mzR:::copyWriteMSData(filename = fnew,
67
+    checkException(mzR:::copyWriteMSData(file = fnew,
68 68
                                          original_file = orig_file,
69
-                                         header = hdr_sub, data = pks_sub,
69
+                                         header = hdr_sub, object = pks_sub,
70 70
                                          backend = "pwiz"))
71 71
     hdr_sub$seqNum <- seq_len(nrow(hdr_sub))
72 72
     ## mzML
73
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
74
-                          header = hdr_sub, data = pks_sub, backend = "pwiz",
75
-                          outformat = "mzml")
73
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
74
+                         header = hdr_sub, object = pks_sub, backend = "pwiz",
75
+                         outformat = "mzml")
76 76
     ## Check content is same
77 77
     mzml_new <- openMSfile(fnew, backend = "pwiz")
78 78
     pks_new <- peaks(mzml_new)
... ...
@@ -83,9 +83,9 @@ test_copyWriteMSData <- function() {
83 83
     checkEquals(pks_new, pks_sub)
84 84
     checkEquals(hdr_new, hdr_sub)
85 85
     ## mzXML
86
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
87
-                          header = hdr_sub, data = pks_sub, backend = "pwiz",
88
-                          outformat = "mzxml")
86
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
87
+                         header = hdr_sub, object = pks_sub, backend = "pwiz",
88
+                         outformat = "mzxml")
89 89
     ## Check content is same
90 90
     mzml_new <- openMSfile(fnew, backend = "pwiz")
91 91
     pks_new <- peaks(mzml_new)
... ...
@@ -108,27 +108,27 @@ test_copyWriteMSData <- function() {
108 108
     ## wrong header.
109 109
     ## wrong spectra.
110 110
     ## wrong data processing.
111
-    checkException(mzR:::copyWriteMSData(filename = fnew,
112
-                                         original_file = orig_file,
113
-                                         header = pks, data = hdr,
114
-                                         backend = "pwiz"))
115
-    checkException(mzR:::copyWriteMSData(filename = fnew,
116
-                                         original_file = orig_file,
117
-                                         header = hdr, data = hdr,
118
-                                         backend = "pwiz"))
119
-    checkException(mzR:::copyWriteMSData(filename = fnew,
120
-                                         original_file = orig_file,
121
-                                         header = hdr, data = pks,
122
-                                         backend = "Ramp"))
123
-    checkException(mzR:::copyWriteMSData(filename = fnew,
124
-                                         original_file = "somefile",
125
-                                         header = hdr, data = pks,
126
-                                         backend = "pwiz"))
127
-    checkException(mzR:::copyWriteMSData(filename = fnew,
128
-                                         original_file = orig_file,
129
-                                         header = hdr, data = pks,
130
-                                         backend = "pwiz",
131
-                                         software_processing = c("other")))
111
+    checkException(mzR::copyWriteMSData(file = fnew,
112
+                                        original_file = orig_file,
113
+                                        header = pks, object = hdr,
114
+                                        backend = "pwiz"))
115
+    checkException(mzR::copyWriteMSData(file = fnew,
116
+                                        original_file = orig_file,
117
+                                        header = hdr, object = hdr,
118
+                                        backend = "pwiz"))
119
+    checkException(mzR::copyWriteMSData(file = fnew,
120
+                                        original_file = orig_file,
121
+                                        header = hdr, object = pks,
122
+                                        backend = "Ramp"))
123
+    checkException(mzR::copyWriteMSData(file = fnew,
124
+                                        original_file = "somefile",
125
+                                        header = hdr, object = pks,
126
+                                        backend = "pwiz"))
127
+    checkException(mzR::copyWriteMSData(file = fnew,
128
+                                        original_file = orig_file,
129
+                                        header = hdr, object = pks,
130
+                                        backend = "pwiz",
131
+                                        software_processing = c("other")))
132 132
     
133 133
     ## INPUT: mzML
134 134
     orig_file <- system.file("proteomics",
... ...
@@ -142,8 +142,8 @@ test_copyWriteMSData <- function() {
142 142
 
143 143
     ## OUTPUT: mzML
144 144
     fnew <- paste0(test_folder, "test_copyWrite.mzML")
145
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
146
-                          header = hdr, data = pks, backend = "pwiz")
145
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
146
+                         header = hdr, object = pks, backend = "pwiz")
147 147
     ## Check content is same
148 148
     mzml_new <- openMSfile(fnew, backend = "pwiz")
149 149
     pks_new <- peaks(mzml_new)
... ...
@@ -156,9 +156,9 @@ test_copyWriteMSData <- function() {
156 156
     
157 157
     ## OUTPUT: mzXML
158 158
     fnew <- paste0(test_folder, "test_copyWrite.mzXML")
159
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
160
-                          header = hdr, data = pks, backend = "pwiz",
161
-                          outformat = "mzxml")
159
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
160
+                         header = hdr, object = pks, backend = "pwiz",
161
+                         outformat = "mzxml")
162 162
     ## Check content is same
163 163
     mzml_new <- openMSfile(fnew, backend = "pwiz")
164 164
     pks_new <- peaks(mzml_new)
... ...
@@ -194,9 +194,9 @@ test_copyWriteMSData <- function() {
194 194
     pks_sub <- pks[idx]
195 195
     hdr_sub$seqNum <- 1:nrow(hdr_sub)
196 196
     fnew <- paste0(test_folder, "test_copyWrite.mzML")
197
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
198
-                          header = hdr_sub, data = pks_sub, backend = "pwiz",
199
-                          outformat = "mzml")
197
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
198
+                         header = hdr_sub, object = pks_sub, backend = "pwiz",
199
+                         outformat = "mzml")
200 200
     ## Check content is same
201 201
     mzml_new <- openMSfile(fnew, backend = "pwiz")
202 202
     pks_new <- peaks(mzml_new)
... ...
@@ -211,9 +211,9 @@ test_copyWriteMSData <- function() {
211 211
 
212 212
     ## Subset with mzXML
213 213
     fnew <- paste0(test_folder, "test_copyWrite.mzXML")
214
-    mzR:::copyWriteMSData(filename = fnew, original_file = orig_file,
215
-                          header = hdr_sub, data = pks_sub, backend = "pwiz",
216
-                          outformat = "mzxml")
214
+    mzR::copyWriteMSData(file = fnew, original_file = orig_file,
215
+                         header = hdr_sub, object = pks_sub, backend = "pwiz",
216
+                         outformat = "mzxml")
217 217
     ## Check content is same
218 218
     mzml_new <- openMSfile(fnew, backend = "pwiz")
219 219
     pks_new <- peaks(mzml_new)
... ...
@@ -244,9 +244,9 @@ test_copyWriteMSData <- function() {
244 244
     
245 245
     ## mzML
246 246
     out_file <- paste0(test_folder, "test_copyWrite.mzML")
247
-    mzR:::copyWriteMSData(filename = out_file, original_file = test_file,
248
-                          header = hdr, data = pks,
249
-                          software_processing = c("MSnbase", "2.3.8"))
247
+    mzR::copyWriteMSData(file = out_file, original_file = test_file,
248
+                         header = hdr, object = pks,
249
+                         software_processing = c("MSnbase", "2.3.8"))
250 250
     in_file <- openMSfile(out_file, backend = "pwiz")
251 251
     hdr_2 <- header(in_file)
252 252
     pks_2 <- peaks(in_file)
... ...
@@ -258,9 +258,9 @@ test_copyWriteMSData <- function() {
258 258
     
259 259
     ## mzXML output:
260 260
     out_file <- paste0(test_folder, "test_copyWrite.mzXML")
261
-    mzR:::copyWriteMSData(filename = out_file, original_file = test_file,
262
-                          header = hdr, data = pks, outformat = "mzXML",
263
-                          software_processing = c("MSnbase", "2.3.8"))
261
+    mzR::copyWriteMSData(file = out_file, original_file = test_file,
262
+                         header = hdr, object = pks, outformat = "mzXML",
263
+                         software_processing = c("MSnbase", "2.3.8"))
264 264
     in_file <- openMSfile(out_file, backend = "pwiz")
265 265
     hdr_2 <- header(in_file)
266 266
     pks_2 <- peaks(in_file)
... ...
@@ -287,7 +287,7 @@ test_writeMSData <- function() {
287 287
 
288 288
     ## mzML
289 289
     out_file <- paste0(test_folder, "test_write.mzML")
290
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks)
290
+    writeMSData(file = out_file, header = hdr, object = pks)
291 291
     in_file <- openMSfile(out_file, backend = "pwiz")
292 292
     hdr_2 <- header(in_file)
293 293
     pks_2 <- peaks(in_file)
... ...
@@ -299,7 +299,7 @@ test_writeMSData <- function() {
299 299
     hdr_sub <- hdr[c(1, 3, 5), ]
300 300
     hdr_sub$seqNum <- 1:nrow(hdr_sub)
301 301
     pks_sub <- pks[c(1, 3, 5)]
302
-    mzR:::writeMSData(out_file, header = hdr_sub, data = pks_sub)
302
+    writeMSData(pks_sub, out_file, header = hdr_sub)
303 303
     in_file <- openMSfile(out_file)
304 304
     hdr_sub_2 <- header(in_file)
305 305
     pks_sub_2 <- peaks(in_file)
... ...
@@ -319,8 +319,8 @@ test_writeMSData <- function() {
319 319
     
320 320
     ## mzXML output:
321 321
     out_file <- paste0(test_folder, "test_write.mzXML")
322
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks,
323
-                      outformat = "mzXML")
322
+    writeMSData(file = out_file, header = hdr, object = pks,
323
+                outformat = "mzXML")
324 324
     in_file <- openMSfile(out_file, backend = "pwiz")
325 325
     hdr_2 <- header(in_file)
326 326
     pks_2 <- peaks(in_file)
... ...
@@ -331,8 +331,8 @@ test_writeMSData <- function() {
331 331
     hdr_sub <- hdr[c(1, 3, 5), ]
332 332
     hdr_sub$seqNum <- 1:nrow(hdr_sub)
333 333
     pks_sub <- pks[c(1, 3, 5)]
334
-    mzR:::writeMSData(out_file, header = hdr_sub, data = pks_sub,
335
-                      outformat = "mzXML")
334
+    writeMSData(file = out_file, header = hdr_sub, object = pks_sub,
335
+                outformat = "mzXML")
336 336
     in_file <- openMSfile(out_file)
337 337
     hdr_sub_2 <- header(in_file)
338 338
     pks_sub_2 <- peaks(in_file)
... ...
@@ -373,7 +373,7 @@ test_writeMSData <- function() {
373 373
     
374 374
     ## mzML
375 375
     out_file <- paste0(test_folder, "test_write.mzML")
376
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks)
376
+    writeMSData(pks, file = out_file, header = hdr)
377 377
     in_file <- openMSfile(out_file, backend = "pwiz")
378 378
     hdr_2 <- header(in_file)
379 379
     pks_2 <- peaks(in_file)
... ...
@@ -390,8 +390,8 @@ test_writeMSData <- function() {
390 390
 
391 391
     ## mzXML output:
392 392
     out_file <- paste0(test_folder, "test_write.mzXML")
393
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks,
394
-                      outformat = "mzXML")
393
+    writeMSData(file = out_file, header = hdr, object = pks,
394
+                outformat = "mzXML")
395 395
     in_file <- openMSfile(out_file, backend = "pwiz")
396 396
     hdr_2 <- header(in_file)
397 397
     pks_2 <- peaks(in_file)
... ...
@@ -424,9 +424,9 @@ test_writeMSData <- function() {
424 424
     hdr_sub$seqNum <- 1:nrow(hdr_sub)
425 425
     pks_sub <- pks[idx]
426 426
     fnew <- paste0(test_folder, "test_copyWrite.mzML")
427
-    mzR:::writeMSData(filename = fnew,
428
-                      header = hdr_sub, data = pks_sub, backend = "pwiz",
429
-                      outformat = "mzml")
427
+    writeMSData(file = fnew,
428
+                header = hdr_sub, object = pks_sub, backend = "pwiz",
429
+                outformat = "mzml")
430 430
     ## Check content is same
431 431
     mzml_new <- openMSfile(fnew, backend = "pwiz")
432 432
     pks_new <- peaks(mzml_new)
... ...
@@ -450,8 +450,8 @@ test_writeMSData <- function() {
450 450
     hdr_sub$seqNum <- 1:nrow(hdr_sub)
451 451
     pks_sub <- pks[idx]
452 452
     fnew <- paste0(test_folder, "test_copyWrite.mzXML")
453
-    mzR:::writeMSData(filename = fnew, header = hdr_sub, data = pks_sub,
454
-                      backend = "pwiz", outformat = "mzxml")
453
+    writeMSData(file = fnew, header = hdr_sub, object = pks_sub,
454
+                backend = "pwiz", outformat = "mzxml")
455 455
     ## Check content is same
456 456
     mzml_new <- openMSfile(fnew, backend = "pwiz")
457 457
     pks_new <- peaks(mzml_new)
... ...
@@ -482,7 +482,7 @@ test_writeMSData <- function() {
482 482
     
483 483
     ## mzML
484 484
     out_file <- paste0(test_folder, "test_write.mzML")
485
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks)
485
+    writeMSData(file = out_file, header = hdr, object = pks)
486 486
     in_file <- openMSfile(out_file, backend = "pwiz")
487 487
     hdr_2 <- header(in_file)
488 488
     pks_2 <- peaks(in_file)
... ...
@@ -492,8 +492,8 @@ test_writeMSData <- function() {
492 492
 
493 493
     ## mzXML output:
494 494
     out_file <- paste0(test_folder, "test_write.mzXML")
495
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks,
496
-                      outformat = "mzXML")
495
+    writeMSData(file = out_file, header = hdr, object = pks,
496
+                outformat = "mzXML")
497 497
     in_file <- openMSfile(out_file, backend = "pwiz")
498 498
     hdr_2 <- header(in_file)
499 499
     pks_2 <- peaks(in_file)
... ...
@@ -511,7 +511,7 @@ test_writeMSData <- function() {
511 511
 
512 512
     ## mzML
513 513
     out_file <- paste0(test_folder, "test_write.mzML")
514
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks)
514
+    writeMSData(file = out_file, header = hdr, object = pks)
515 515
     in_file <- openMSfile(out_file, backend = "pwiz")
516 516
     hdr_2 <- header(in_file)
517 517
     pks_2 <- peaks(in_file)
... ...
@@ -521,8 +521,8 @@ test_writeMSData <- function() {
521 521
 
522 522
     ## mzXML output:
523 523
     out_file <- paste0(test_folder, "test_write.mzXML")
524
-    mzR:::writeMSData(filename = out_file, header = hdr, data = pks,
525
-                      outformat = "mzXML")
524
+    writeMSData(file = out_file, header = hdr, object = pks,
525
+                outformat = "mzXML")
526 526
     in_file <- openMSfile(out_file, backend = "pwiz")
527 527
     hdr_2 <- header(in_file)
528 528
     pks_2 <- peaks(in_file)
... ...
@@ -7,12 +7,20 @@
7 7
 }
8 8
 
9 9
 \usage{
10
-  copyWriteMSData(filename, original_file, header, data, backend =
10
+  copyWriteMSData(object, file, original_file, header, backend =
11 11
     "pwiz", outformat = "mzml", rtime_seconds = TRUE, software_processing)
12 12
 }
13 13
 
14 14
 \arguments{
15
-  \item{filename}{
15
+  
16
+  \item{object}{
17
+    \code{list} containing for each spectrum one \code{matrix} with
18
+    columns \code{mz} (first column) and \code{intensity} (second
19
+    column). See also \code{\link{peaks}} for the method that reads such
20
+    data from an MS file.
21
+  }
22
+
23
+  \item{file}{
16 24
     \code{character(1)} defining the name of the file.
17 25
   }
18 26
 
... ...
@@ -26,14 +34,7 @@
26 34
     the format as the \code{data.frame} returned by the
27 35
     \code{\link{header}} method.
28 36
   }
29
-  
30
-  \item{data}{
31
-    \code{list} containing for each spectrum one \code{matrix} with
32
-    columns \code{mz} (first column) and \code{intensity} (second
33
-    column). See also \code{\link{peaks}} for the method that reads such
34
-    data from an MS file.
35
-  }
36
-  
37
+    
37 38
   \item{backend}{
38 39
     \code{character(1)} defining the backend that should be used for
39 40
     writing. Currently only \code{"pwiz"} backend is supported.
... ...
@@ -113,6 +114,6 @@ pks <- lapply(pks, function(z) {
113 114
 ## Copy metadata and additional information from the originating file
114 115
 ## and save it, along with the modified data, to a new mzML file.
115 116
 out_file <- tempfile()
116
-copyWriteMSData(filename = out_file, original_file = fl,
117
-    header = hdr, data = pks)
117
+copyWriteMSData(pks, file = out_file, original_file = fl,
118
+    header = hdr)
118 119
 }
... ...
@@ -1,15 +1,24 @@
1 1
 \name{writeMSData}
2 2
 \alias{writeMSData}
3
+\alias{writeMSData,list,character-method}
3 4
 
4 5
 \title{
5 6
   Write MS spectrum data to an MS file
6 7
 }
7 8
 \usage{
8
-writeMSData(filename, header, data, backend = "pwiz", outformat =
9
-   "mzml", rtime_seconds = TRUE, software_processing)
9
+
10
+\S4method{writeMSData}{list,character}(object, file, header,
11
+    backend = "pwiz", outformat = "mzml", rtime_seconds = TRUE,
12
+    software_processing)
10 13
 }
11 14
 \arguments{
12
-  \item{filename}{
15
+  \item{object}{
16
+    \code{list} containing for each spectrum one \code{matrix} with
17
+    columns \code{mz} (first column) and \code{intensity} (second
18
+    column). See also \code{\link{peaks}} for the method that reads such
19
+    data from an MS file.
20
+  }
21
+  \item{file}{
13 22
     \code{character(1)} defining the name of the file.
14 23
   }
15 24
   \item{header}{
... ...
@@ -17,12 +26,6 @@ writeMSData(filename, header, data, backend = "pwiz", outformat =
17 26
     the format as the \code{data.frame} returned by the
18 27
     \code{\link{header}} method.
19 28
   }
20
-  \item{data}{
21
-    \code{list} containing for each spectrum one \code{matrix} with
22
-    columns \code{mz} (first column) and \code{intensity} (second
23
-    column). See also \code{\link{peaks}} for the method that reads such
24
-    data from an MS file.
25
-  }
26 29
   \item{backend}{
27 30
     \code{character(1)} defining the backend that should be used for
28 31
     writing. Currently only \code{"pwiz"} backend is supported.
... ...
@@ -84,5 +87,5 @@ pks <- lapply(pks, function(z) {
84 87
 
85 88
 ## Write the data to a mzML file.
86 89
 out_file <- tempfile()
87
-writeMSData(filename = out_file, header = hdr, data = pks)
90
+writeMSData(object = pks, file = out_file, header = hdr)
88 91
 }