View source: R/guessEnriched.R
guessEnriched | R Documentation |
As clearly enriched peptides will always have a 100% posterior probability of enrichment, BEER removes these peptides a priori to running the model. Clearly enriched peptides can be identified using edgeR estimated fold-changes or maximum likelihood estimates based on the specified prior parameters. Additional parameters for each method can be found in the details below.
guessEnriched(object, method = "mle", ...)
object |
a |
method |
one of "mle" or "edgeR", specifying which method to use to identify clearly enriched peptides |
... |
additional parameters dependent on the method used. See details for more information |
edgeR. Identification of clearly enriched peptides relies
on edgeR fold-change estimates, so edgeR
must be run on
the PhIPData
object beforehand. Additional parameters
for identifying clearly enriched peptides based on edgeR estimated
fold-changes are listed below:
object
: a PhIPData
object.
threshold
: minimum estimated fc for a peptide to be
considered super-enriched. The default value is 15.
fc.name
: assay name corresponding to the assay that stores
the edgeR estimated log2 fold-changes.
MLE. As the number of reads tends to be quite large, the estimates for the proportion of reads pulled are generally accurate. Clearly enriched peptides are identified by first comparing the observed read count to the expected read count based on the beads-only prior parameters. Peptides with observed read counts larger than 5 times the expected read counts are temporarily labeled as enriched, and attenuation constants are estimated by regressing the observed read counts on the expected read counts for all non-enriched peptides. Using this attenuation constant, peptides with fold-changes above some predefined threshold after adjusting for the attenuation constant are considered enriched. Parameters for identifying clearly enriched peptides using the MLE approach are listed below.
object
: a PhIPData
object.
threshold
: minimum estimated fc for a peptide to be
considered super-enriched.
beads.prior
: data.frame of prior parameters for beads-only
samples.
a logical matrix of the with the same dimensions as object
indicating which peptides are considered super-enriched.
sim_data <- readRDS(system.file("extdata", "sim_data.rds", package = "beer")) edgeR_out <- runEdgeR(sim_data) guessEnriched(edgeR_out, method = "edgeR", threshold = 15, fc.name = "logfc") guessEnriched(edgeR_out, method = "mle", beads.prior = getAB(edgeR_out, method = "edgeR"), threshold = 15 )
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.