View source: R/regressIntensity.R
regressIntensity | R Documentation |
Performs linear regression on protein intensities based on selected protein (qPLEX-RIME bait)
regressIntensity(MSnSetObj, ProteinId, controlInd = NULL, plot = TRUE)
MSnSetObj |
MSnSet; an object of class MSnSet |
ProteinId |
character; Uniprot protein ID |
controlInd |
numeric; index of IgG within MSnSet |
plot |
character; Whether or not to plot the QC histograms |
This function performs regression based analysis upon protein intensities based on a selected protein. In qPLEX RIME this method could be used to regress out the effect of target protein on other interactors. This function corrects this dependency of many proteins on the target protein levels by linear regression. It sets the target protein levels as the independent variable (x) and each of the other proteins as the dependent variable (y). The resulting residuals of the linear regressions y=ax+b are the protein levels corrected for target protein dependency.
An object of class MSnSet
(see MSnSet-class
).
This consists of corrected protein levels. In addition, the function can
also plot histograms of correlation of target protein with all other
proteins before and after this correction.
data(human_anno)
data(exp3_OHT_ESR1)
MSnSet_data <- convertToMSnset(exp3_OHT_ESR1$intensities_qPLEX1,
metadata=exp3_OHT_ESR1$metadata_qPLEX1,
indExpData=c(7:16),
Sequences=2,
Accessions=6)
MSnset_P <- summarizeIntensities(MSnSet_data, sum, human_anno)
IgG_ind <- which(pData(MSnset_P)$SampleGroup == "IgG")
MSnset_reg <- regressIntensity(MSnset_P,
controlInd=IgG_ind,
ProteinId="P03372")
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.