identifFixedModif {wrTopDownFrag} | R Documentation |
Identify Fixed Modifications
Description
Identify peptide/protein fragment based on experimental m/z values 'expMass' for given range of aa-length.
Internally all possible fragments will be predicted and their mass compared to the experimental values (argument expMass
).
Usage
identifFixedModif(
prot,
expMass,
minFragSize = 5,
maxFragSize = 60,
indexStart = 1,
suplPepTab = NULL,
internFra = TRUE,
filtChargeCatch = TRUE,
maxMod = c(p = 3, h = 1, k = 1, o = 1, m = 1, n = 1, u = 1, r = 1, s = 1),
modTy = NULL,
specModif = NULL,
knownMods = NULL,
identMeas = "ppm",
limitIdent = 5,
filtAmbiguous = FALSE,
recalibrate = FALSE,
chargeCatchFilter = TRUE,
massTy = "mono",
prefFragPat = NULL,
silent = FALSE,
debug = FALSE,
callFrom = NULL
)
Arguments
prot |
(character) amino-acid sequene of peptide or protein |
expMass |
(numeric) erperimental masses to identify peptides from |
minFragSize |
(integer) min number of AA residues for considering peptide fragments |
maxFragSize |
(integer) max number of AA residues for considering peptide fragments |
indexStart |
(integer) for starting at correct index (if not 1) |
suplPepTab |
(matrix) additional peptides to be add to theoretical peptides |
internFra |
(logical) decide whether internal fragments should be cosiered |
filtChargeCatch |
(logical) by default removing of all fragments not containing a (polar) charge-cathing residue |
maxMod |
(integer) maximum number of residue modifications to be consiered in fragments (values >1 will increase complexity and RAM consumption) |
modTy |
(character) type of fixed and variable modifications |
specModif |
(list) supplemental custom fixed or variable modifications (eg Zn++ at given residue) |
knownMods |
(character) optional custom alternative to |
identMeas |
(character) default 'ppm' |
limitIdent |
(character) thershold for identification in 'identMeas' units |
filtAmbiguous |
(logical) allows filtering/removing ambiguous results (ie same mass peptides) |
recalibrate |
(logical or numeric) may be direct recalibration-factor (numeric,length=1), if 'TRUE' fresh determination of 'recalibFact' or 'FALSE' (no action); final recalibration-factor used exported in result as $recalibFact |
chargeCatchFilter |
(logical) optionally remove all peptides not containing charge-catch AAs (K, R, H, defined via |
massTy |
(character) 'mono' or 'average' |
prefFragPat |
(numeric) pattern for preferential fragmentation (see also Haverland 2017), if |
silent |
(logical) suppress messages |
debug |
(logical) additional messages and objects exportet to current session for debugging |
callFrom |
(character) allow easier tracking of message(s) produced |
Value
list, ie result of massMatch() on 'pepTab' and 'expMass'
See Also
Examples
protP <- c(protP="PEPTIDEKR")
obsMassX <- cbind(a=c(199.1077,296.1605,397.2082,510.2922,625.3192),
b=c(227.1026,324.1554,425.2031,538.2871,653.3141),
x=c(729.2937,600.2511,503.1984,402.1507,289.0666),
y=c(703.3145,574.2719,477.2191,376.1714,263.0874))
rownames(obsMassX) <- c("E","P","T","I","D") # all 1 & 7 ions not included
identP1 <- identifFixedModif(prot=protP, expMass=as.numeric(obsMassX), minFragSize=2,
maxFragSize=7,modTy=list(basMod=c("b","y"))) # looks ok
identP2 <- identifFixedModif(prot=protP, expMass=as.numeric(obsMassX), minFragSize=2,
maxFragSize=7, modTy=list(basMod=c("a","x"), varMod=c("h","o","r","m")))
head(identP1$preMa,n=17) # predicted masses incl fixed modif
head(identP2$preMa,n=17) # predicted masses incl fixed modif