Title: | Statistical analysis of peptide microarrays |
---|---|
Description: | Statistical analysis of peptide microarrays |
Authors: | Raphael Gottardo, Gregory C Imholte, Renan Sauteraud, Mike Jiang |
Maintainer: | Gregory C Imholte <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.41.0 |
Built: | 2024-10-30 09:24:48 UTC |
Source: | https://github.com/bioc/pepStat |
Correct intensities by substracting PRE visit sample intensities.
baseline_correct(pSet, verbose = FALSE)
baseline_correct(pSet, verbose = FALSE)
pSet |
A |
verbose |
A |
The function will try to pair as many sample as possible. The remaining subjects with a POST and no PRE will use the average expression of all baseline samples. Subjects with baseline only will not be represented in the resulting matrix.
A matrix
of the baseline corrected intensities, with as many
columns as there are samples POST visit
Renan Sauteraud
Correct intensities by substracting PRE visit sample intensities.
baselineCorrect.pSet(pSet, verbose = FALSE)
baselineCorrect.pSet(pSet, verbose = FALSE)
pSet |
A |
verbose |
A |
If samples are not PAIRED (One PRE and POST for each ptid), then the average expression of all PRE visit samples is substracted from each sample.
A matrix
of the baseline corrected intensities, with as many
columns as there are samples POST visit
Raphael Gottardo, Gregory Imholte
Constructor to create peptide collection to be used in
summarizePeptides
.
create_db(position)
create_db(position)
position |
A |
position
can have additional columns. These columns will be kept in
the resulting peptide collection. This is especially useful to include
clades and grouping parameters for the makeCalls
function.
If the input contains all the z-scores (z1 to z5), then they will not be re-calculated. If some (but not all) z-scores are missing, a warning message will be sent and the z-scores are re-calculated.
Renan Sauteraud
#construct data.frame object AA <- c("A", "C", "D", "E", "F", "G", "H", "I", "K", "L", "M", "N", "P", "Q","R", "S", "T", "V", "W", "Y") starts <- seq(1, 30, 3) ends <- starts + 14 peptides <- sapply(1:10, function(x) { paste0(AA[floor(runif(15, 1, 20))], collapse = "") }) data <- data.frame(start = starts, end = ends, peptide = peptides) #from data.frame new_pep <- create_db(data) #from GRanges new_pep <- create_db(new_pep)
#construct data.frame object AA <- c("A", "C", "D", "E", "F", "G", "H", "I", "K", "L", "M", "N", "P", "Q","R", "S", "T", "V", "W", "Y") starts <- seq(1, 30, 3) ends <- starts + 14 peptides <- sapply(1:10, function(x) { paste0(AA[floor(runif(15, 1, 20))], collapse = "") }) data <- data.frame(start = starts, end = ends, peptide = peptides) #from data.frame new_pep <- create_db(data) #from GRanges new_pep <- create_db(new_pep)
After normalization and data smoothing, this last step makes the call for each peptide of the peptideSet after baseline correcting the peptide intenstities.
makeCalls(peptideSet, cutoff = 1.2, method = "absolute", freq = TRUE, group = NULL, verbose = FALSE)
makeCalls(peptideSet, cutoff = 1.2, method = "absolute", freq = TRUE, group = NULL, verbose = FALSE)
peptideSet |
A |
cutoff |
A |
method |
A |
freq |
A |
group |
A |
verbose |
A |
This function requires specific variables ptid and visit in pData(peptideSet).
The variable ptid
should indicate subjects, and the variable visit
should be a factor with levels pre and post.
If slides are paired for subjects, intensities corresponding to post-visit are substracted from pre. If slides are not paired, slides with pre have intensities averaged by peptides, and averaged peptide intensities are subtracted from slides that have entry post. Calls are made on these baseline corrected intensities.
When method = FDR, a left-tail method is used to generate a threshold controlling
the False Discovery Rate at level cutoff
. When method = absolute, Intensities
exceeding the threshold are labelled as positive.
When freq = TRUE a group variable may be specified. The argument group indicates the name of a variable in pData(peptideSet) by which positive calls should be grouped. The call frequency for each peptide is calculated within groups.
If freq = TRUE, a numeric
matrix
with peptides as rows and
groups as columns where the values are the frequency of response in the group. If
freq = FALSE, a logical
matrix
indicating binding events for each
peptide in each subject.
Greg Imholte
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
This function reads GenePix results (.gpr) files and creates a peptideSet object gathering experiment information.
makePeptideSet(files = NULL, path = NULL, mapping.file = NULL, use.flags = FALSE, rm.control.list = NULL, empty.control.list = NULL, bgCorrect.method = "normexp", log = TRUE, check.row.order = FALSE, verbose = FALSE)
makePeptideSet(files = NULL, path = NULL, mapping.file = NULL, use.flags = FALSE, rm.control.list = NULL, empty.control.list = NULL, bgCorrect.method = "normexp", log = TRUE, check.row.order = FALSE, verbose = FALSE)
files |
A |
path |
A |
mapping.file |
A |
use.flags |
A |
rm.control.list |
A |
empty.control.list |
A |
bgCorrect.method |
A |
log |
A |
check.row.order |
A |
verbose |
A |
GenePix results files (.gpr) are read when found in either the files or path arguments. By default, the foreground and background median intensities of the "red" channels, "F635 Median" and "B635 Median", are read. The background correction specified in bgCorrect.method is passed to the backgroundCorrect method in the limma package.
The mapping.file can be either a filename or a data.frame. In any case, it should
contain at least three columns labeled "filename", "ptid" and "visit". The
filenames given here should match those read from the path or files argument,
or be a subset of it. For each ptid (patient ID), the visit column should have at
least one "pre" and one "post" sample. Any additional column will be kept into
the resulting peptideSet
and can be used later on for groupping.
If check.row.order = TRUE, the final set of probes is taken to be those with IDs found in all arrays that were read.
Row, Column and Block spatial array position for each probe are stored in the
featureRanges
slot of the returned object.
A peptideSet
object that contain the intensities, peptide
sequences and annotations available in the mapping file.
Raphael Gottardo, Gregory Imholte
peptideSet
, read.maimages
,
backgroundCorrect
# Read gpr files library(pepDat) mapFile <- system.file("extdata/mapping.csv", package = "pepDat") dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) # Specify controls to be removed and empty controls # to be used for background correction. pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log = TRUE, rm.control.list = c("JPT-control", "Ig", "Cy3"), empty.control.list= c("empty", "blank control"))
# Read gpr files library(pepDat) mapFile <- system.file("extdata/mapping.csv", package = "pepDat") dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) # Specify controls to be removed and empty controls # to be used for background correction. pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log = TRUE, rm.control.list = c("JPT-control", "Ig", "Cy3"), empty.control.list= c("empty", "blank control"))
This function is used to normalize the peptide microarray data using sequence information.
normalizeArray(peptideSet, method = "ZpepQuad", robust = TRUE, centered = TRUE)
normalizeArray(peptideSet, method = "ZpepQuad", robust = TRUE, centered = TRUE)
peptideSet |
A |
method |
A |
robust |
A |
centered |
A |
The available methods are "Zpep" and "ZpepQuad". These methods fit a linear model using either linear or linear and quadratic terms (respectively), regressing intensity on the peptides' five Z-scale scores. A peptide Z-scale score is obtained by summing over the Z-scale values in Sandburg et al (1998) of the amino acids the peptide comprises.
Peptide Z-scale scores may be provided in the featureRange slot of peptideSet.
This slot is a GRanges
object x, and the function will seek five columns
labelled z1 through z5 in values(x). If these are not found, the function attempts
to calculate Z-scales from sequence information found in peptide(peptideSet)
If robust = TRUE the linear model is fit with t_4 distributed errors. The method returns the residuals of each peptide intensity in the fitted linear model. If centered = TRUE the fitted intercept term is added back to the residuals of the fit.
A peptideSet
object with updated normalized intensity values.
Raphael Gottardo, Gregory Imholte
Sandberg, M., Eriksson, L., Jonsson, J., Sjostrom, M., and Wold, S. (1998). New chemical descriptors relevant for the design of biologically active peptides. A multivariate characterization of 87 amino acids. Journal of Medicinal Chemistry 41, 2481-2491.
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
This class gathers all information from gpr files, annotation data and sequence data
See ?`peptideSet-methods`
for a list of accessors and method associated
with the class.
A GRanges
object. The ranges and sequences of
the peptides and their associated annotation.
An AnnotatedDataFrame
. Annotation for the samples.
Slots inherited from ExpressionSet
.
Greg Imholte
ExpressionSet
, peptideSet-methods
Methods for handling peptideSet objects
nrow(x)
:The number of peptides in x.
ncol(x)
:The number of samples in x.
start(x)
:Get the starts of the peptides.
end(x)
:Get the ends of the peptides.
width(x)
:Get the widths of the peptides.
position(x)
:Get the coordinates of the central amino-acid of
each peptide, given by: round((start(x) + end(x))/2)
.
ranges(x)
:Returns a GRanges
object that contains
the annotations for the peptides.
ranges(x)<- value
Set annotations for the peptides.
values(x)
:Returns a SplitDataFrameList
. Accessor for the
values of the featureRange slot.
clade(x)
:If available, returns the clade information for each
peptide as a matrix
.
peptide(x)
:Get the sequence of the peptides.
peptide(x) <- value
Set the sequence of the peptides.
featureID(x)
:Get the ID of the peptides.
pepZscore(x)
:If available, returns a matrix
of the zScores
for each peptide.
pepZscore(x) <- value
Set the zScores for each peptide
show(object)
:Display a peptideSet object.
summary(object)
:Summarize a peptideSet object.
x[i, j]
:Subset x by peptides (i), or samples (j).
subset(x, subset, drop=FALSE)
:Subset x given an expression 'subset'.
Plot a color image of the intensities on a slide. These plots are helpful to diagnose spatial abnormalities.
plotArrayImage(peptideSet, array.index = NULL, low = "white", high = "steelblue", ask = dev.interactive(orNone = TRUE) & 1 < length(array.index)) plotArrayResiduals(peptideSet, array.index = NULL, smooth = FALSE, low = "blue", high = "red", ask = dev.interactive(orNone = TRUE) & 1 < length(array.index))
plotArrayImage(peptideSet, array.index = NULL, low = "white", high = "steelblue", ask = dev.interactive(orNone = TRUE) & 1 < length(array.index)) plotArrayResiduals(peptideSet, array.index = NULL, smooth = FALSE, low = "blue", high = "red", ask = dev.interactive(orNone = TRUE) & 1 < length(array.index))
peptideSet |
A |
array.index |
A vector subsetting |
smooth |
A |
low |
A |
high |
A |
ask |
A |
The most coherent results are achieved when the peptideSet
object is
read with makePeptideSet
with empty.control.list = NULL and rm.control.list
= NULL
Gregory Imholte
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
Tabulate the results of a peptide microarray analysis.
restab(peptideSet, calls)
restab(peptideSet, calls)
peptideSet |
A |
calls |
A |
The peptideSet should be the one used in the function call to makeCalls
that generated the calls used. They should have identical peptides.
A data.frame
with the peptides and some information from the
peptideSet
as well as the frequency of binding for each group of the
calls.
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
Launches the pepStat
Shiny application, providing an interactive
interface for constructing peptide sets, normalizing intensities, generating
calls. Quality control is also facilitated through interactive plotting
features.
shinyPepStat()
shinyPepStat()
if (interactive()) { shinyPepStat() }
if (interactive()) { shinyPepStat() }
This function applies a sliding mean window to intensities to reduce noise generated by experimental variation, as well as take advantage of the overlapping nature of array peptides to share signal.
slidingMean(peptideSet, width = 9, verbose = FALSE, split.by.clade = TRUE)
slidingMean(peptideSet, width = 9, verbose = FALSE, split.by.clade = TRUE)
peptideSet |
A |
width |
A |
verbose |
A |
split.by.clade |
A |
Peptide membership in the sliding mean window is determined by its position and the width argument. Two peptides are in the same window if the difference in their positions is less than or equal to width/2. A peptide's position is taken to be position(peptideSet).
A peptide's intensity is replaced by the mean of all peptide intensities within the peptide's sliding mean window.
When split.by.clade = TRUE, peptides are smoothed within clades defined by the clade column of the GRanges object occupying the featureRange slot of peptideSet. If set to FALSE, a peptide at a given position will borrow information from the neighboring peptides as well as the ones from other clades around this position.
A peptideSet
object with smoothed intensities.
Gregory Imholte
summarizePeptides
, normalizeArray
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
This function merges the replicates and adds information from a peptide collection to a peptideSet. This collection can include coordinates, alignment information, Z-scales, and other peptide information.
summarizePeptides(peptideSet, summary = "median", position = NULL)
summarizePeptides(peptideSet, summary = "median", position = NULL)
peptideSet |
A |
summary |
A |
position |
A |
The object in the position argument will be passed to create_db
, it
can either be a GRanges
object with a peptide as a metadata column, or
a data.frame
that can be used to create such GRanges
.
Some peptide collections can be found in the pepDat
package.
An object of class peptideSet
with added columns and updated ranges.
Raphael Gottardo, Greory Imholte
makePeptideSet
, create_db
,
create_db
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }
## This example curated from the vignette -- please see vignette("pepStat") ## for more information if (require("pepDat")) { ## Get example GPR files + associated mapping file dirToParse <- system.file("extdata/gpr_samples", package = "pepDat") mapFile <- system.file("extdata/mapping.csv", package = "pepDat") ## Make a peptide set pSet <- makePeptideSet(files = NULL, path = dirToParse, mapping.file = mapFile, log=TRUE) ## Plot array images -- useful for quality control plotArrayImage(pSet, array.index = 1) plotArrayResiduals(pSet, array.index = 1, smooth = TRUE) ## Summarize peptides, using pep_hxb2 as the position database data(pep_hxb2) psSet <- summarizePeptides(pSet, summary = "mean", position = pep_hxb2) ## Normalize the peptide set pnSet <- normalizeArray(psSet) ## Smooth psmSet <- slidingMean(pnSet, width = 9) ## Make calls calls <- makeCalls(psmSet, freq = TRUE, group = "treatment", cutoff = .1, method = "FDR", verbose = TRUE) ## Produce a summary of the results summary <- restab(psmSet, calls) }