given table of peptide protein assigments generate matrix

prepareMatrix(
  data,
  proteinID = "proteinID",
  peptideID = "strippedSequence",
  weighting = NULL,
  sep = "|"
)

Arguments

data

generated by annotatePeptides

proteinID

protein ID column

peptideID

peptide / precursor ID column

weighting

weight type to use. Options are "one" , "AA" - amino acids, "coverage" - coverage , "inverse" - inverse peptide frequencies

sep

separator for precursor (rownames)

Value

sparse matrix

Examples

#library(prozor)
data(protpepmetashort)
library(Matrix)
colnames(protpepmetashort)
#> [1] "peptideSeq"      "lengthProtein"   "proteinID"       "peptideModSeq"  
#> [5] "precursorCharge" "m"               "lengthPeptide"  
head(protpepmetashort)
#>       peptideSeq lengthProtein               proteinID  peptideModSeq
#> 1 AAAITSDLLESLGR           838   sp|Q62009|POSTN_MOUSE AAAITSDLLESLGR
#> 2 AAAITSDLLESLGR           811 sp|Q62009-2|POSTN_MOUSE AAAITSDLLESLGR
#> 3 AAAITSDLLESLGR           810 sp|Q62009-3|POSTN_MOUSE AAAITSDLLESLGR
#> 4 AAAITSDLLESLGR           784 sp|Q62009-4|POSTN_MOUSE AAAITSDLLESLGR
#> 5 AAAITSDLLESLGR           783 sp|Q62009-5|POSTN_MOUSE AAAITSDLLESLGR
#> 6       AAAQTLER           289    sp|P23492|PNPH_MOUSE       AAAQTLER
#>   precursorCharge        m lengthPeptide
#> 1               2 708.8887            14
#> 2               2 708.8887            14
#> 3               2 708.8887            14
#> 4               2 708.8887            14
#> 5               2 708.8887            14
#> 6               2 430.2353             8
dim(protpepmetashort)
#> [1] 423   7
count = prepareMatrix( protpepmetashort, peptideID = "peptideSeq" )
dim(count)
#> [1] 149 250
inverse = prepareMatrix( protpepmetashort, peptideID = "peptideSeq" , weight = "inverse")
#aa = prepareMatrix(protpepmetashort,  peptideID = "peptideSeq" , weight = "AA")
#xx = prepareMatrix(protpepmetashort,  peptideID = "peptideSeq" , weight = "coverage")
image( as.matrix(count) )


corProt = cor( as.matrix(count) )
par(mfrow =c(1,2))
image(corProt)

#penalise peptides matching many proteins
corProtn = cor( as.matrix(inverse) )
image(corProtn)