SilacProteomicsExperiment-constructor {pulsedSilac} | R Documentation |
Constructor function for the ProteomicsExperiment class object.
It requires at minimum a SilacProteinExperiment
and a
SilacPeptideExperiment
. If the colData, metadata and metaoptions have
been already defined in those it is not necessary to give them again.
SilacProteomicsExperiment( SilacProteinExperiment, SilacPeptideExperiment, colData, linkerDf, metadata, idColProt = NA, idColPept = NA, linkedSubset = TRUE, subsetMode = "protein", conditionCol = NA, timeCol = NA, proteinCol = NA )
SilacProteinExperiment |
A |
SilacPeptideExperiment |
A |
colData |
A |
linkerDf |
A |
metadata |
A |
idColProt |
A |
idColPept |
A |
linkedSubset |
A |
subsetMode |
A |
conditionCol |
A |
timeCol |
A |
proteinCol |
A |
An object of class SilacProteomicsExperiment
.
See SilacProteomicsExperiment-class for details.
See SilacProteomicsExperiment-accessors for details.
## assays assays_protein <- list(expression = matrix(1:9, ncol = 3)) ## colData colData <- data.frame(sample = c('A1', 'A2', 'A3'), condition = c('A', 'A', 'A'), time = c(1, 2, 3)) ## rowData rowData_protein <- data.frame(prot_id = LETTERS[1:3]) ## construct the ProteinExperiment protExp <- SilacProteinExperiment(assays = assays_protein, rowData = rowData_protein, colData = colData, conditionCol = 'condition', timeCol = 'time') ## assays assays_peptide <- list(expression = matrix(1:15, ncol = 3)) ## colData colData <- data.frame(sample = c('A1', 'A2', 'A3'), condition = c('A', 'A', 'A'), time = c(1, 2, 3)) ## rowData rowData_peptide <- data.frame(pept_id = letters[1:5], prot_id = c('A', 'A', 'B', 'C', 'C')) ## construct the ProteinExperiment peptExp <- SilacPeptideExperiment(assays = assays_peptide, rowData = rowData_peptide, colData = colData, conditionCol = 'condition', timeCol = 'time') ## list with the relationships protein_to_peptide <- list(A = c('a', 'b'), B = c('c'), C = c('d', 'e')) ## function to build the data.frame linkerDf <- buildLinkerDf(protIDs = LETTERS[1:3], pepIDs = letters[1:5], protToPep = protein_to_peptide) ProteomicsExp <- SilacProteomicsExperiment(SilacProteinExperiment = protExp, SilacPeptideExperiment = peptExp, linkerDf = linkerDf)