OpenMS
|
Generate assays from a TargetedExperiment. More...
#include <OpenMS/ANALYSIS/OPENSWATH/MRMAssay.h>
Public Types | |
typedef std::vector< OpenMS::TargetedExperiment::Protein > | ProteinVectorType |
typedef std::vector< OpenMS::TargetedExperiment::Peptide > | PeptideVectorType |
typedef std::vector< OpenMS::TargetedExperiment::Compound > | CompoundVectorType |
typedef std::vector< OpenMS::ReactionMonitoringTransition > | TransitionVectorType |
typedef std::map< String, std::vector< const ReactionMonitoringTransition * > > | PeptideTransitionMapType |
typedef std::map< String, std::vector< const ReactionMonitoringTransition * > > | CompoundTransitionMapType |
typedef std::map< String, std::set< std::string > > | ModifiedSequenceMap |
Maps an unmodified sequence to all its modified sequences. More... | |
typedef boost::unordered_map< size_t, ModifiedSequenceMap > | SequenceMapT |
Stores the ModifiedSequenceMap for all SWATH windows. More... | |
typedef std::vector< std::pair< double, std::string > > | FragmentSeqMap |
Describes a fragment sequence map of : "fragment m/z" -> "modified sequence". More... | |
typedef boost::unordered_map< size_t, boost::unordered_map< String, FragmentSeqMap > > | IonMapT |
Stores a mapping : "unmodified sequence" -> FragmentSeqMap for all SWATH windows. More... | |
typedef std::vector< std::pair< std::string, double > > | IonSeries |
Describes an ion series: "ion_type" -> "fragment m/z". More... | |
typedef std::map< String, IonSeries > | PeptideMapT |
Maps a peptide sequence to an ion series: "ion_type" -> "fragment m/z". More... | |
typedef std::map< String, TargetedExperiment::Peptide > | TargetDecoyMapT |
Maps the peptide id (same for target and decoy) to the decoy peptide object. More... | |
Public Types inherited from ProgressLogger | |
enum | LogType { CMD , GUI , NONE } |
Possible log types. More... | |
Public Member Functions | |
MRMAssay () | |
Constructor. More... | |
~MRMAssay () override | |
Destructor. More... | |
void | reannotateTransitions (OpenMS::TargetedExperiment &exp, double precursor_mz_threshold, double product_mz_threshold, const std::vector< String > &fragment_types, const std::vector< size_t > &fragment_charges, bool enable_specific_losses, bool enable_unspecific_losses, int round_decPow=-4) |
Annotates and filters transitions in a TargetedExperiment. More... | |
void | restrictTransitions (OpenMS::TargetedExperiment &exp, double lower_mz_limit, double upper_mz_limit, const std::vector< std::pair< double, double > > &swathes) |
Restrict and filter transitions in a TargetedExperiment. More... | |
void | detectingTransitions (OpenMS::TargetedExperiment &exp, int min_transitions, int max_transitions) |
Select detecting fragment ions. More... | |
void | uisTransitions (OpenMS::TargetedExperiment &exp, const std::vector< String > &fragment_types, const std::vector< size_t > &fragment_charges, bool enable_specific_losses, bool enable_unspecific_losses, bool enable_ms2_precursors, double mz_threshold, const std::vector< std::pair< double, double > > &swathes, int round_decPow=-4, size_t max_num_alternative_localizations=20, int shuffle_seed=-1, bool disable_decoy_transitions=false) |
Annotate UIS / site-specific transitions. More... | |
void | filterMinMaxTransitionsCompound (OpenMS::TargetedExperiment &exp, int min_transitions, int max_transitions) |
Filters target and decoy transitions by intensity, only keeping the top N transitions. More... | |
void | filterUnreferencedDecoysCompound (OpenMS::TargetedExperiment &exp) |
Filters decoy transitions, which do not have respective target transition based on the transitionID. More... | |
Public Member Functions inherited from ProgressLogger | |
ProgressLogger () | |
Constructor. More... | |
virtual | ~ProgressLogger () |
Destructor. More... | |
ProgressLogger (const ProgressLogger &other) | |
Copy constructor. More... | |
ProgressLogger & | operator= (const ProgressLogger &other) |
Assignment Operator. More... | |
void | setLogType (LogType type) const |
Sets the progress log that should be used. The default type is NONE! More... | |
LogType | getLogType () const |
Returns the type of progress log being used. More... | |
void | setLogger (ProgressLoggerImpl *logger) |
Sets the logger to be used for progress logging. More... | |
void | startProgress (SignedSize begin, SignedSize end, const String &label) const |
Initializes the progress display. More... | |
void | setProgress (SignedSize value) const |
Sets the current progress. More... | |
void | endProgress (UInt64 bytes_processed=0) const |
void | nextProgress () const |
increment progress by 1 (according to range begin-end) More... | |
Protected Member Functions | |
std::vector< std::string > | getMatchingPeptidoforms_ (const double fragment_ion, const FragmentSeqMap &ions, const double mz_threshold) |
Check whether fragment ion are unique ion signatures in vector within threshold and return matching peptidoforms. More... | |
int | getSwath_ (const std::vector< std::pair< double, double > > &swathes, const double precursor_mz) |
Get swath index (precursor isolation window ordinal) for a particular precursor. More... | |
bool | isInSwath_ (const std::vector< std::pair< double, double > > &swathes, const double precursor_mz, const double product_mz) |
Check whether the product m/z of a transition falls into the precursor isolation window. More... | |
std::string | getRandomSequence_ (size_t sequence_size, boost::variate_generator< boost::mt19937 &, boost::uniform_int<> > pseudoRNG) |
Generates random peptide sequence. More... | |
std::vector< std::vector< size_t > > | nchoosekcombinations_ (const std::vector< size_t > &n, size_t k) |
Computes all N choose K combinations. More... | |
std::vector< OpenMS::AASequence > | addModificationsSequences_ (const std::vector< OpenMS::AASequence > &sequences, const std::vector< std::vector< size_t > > &mods_combs, const OpenMS::String &modification) |
Generate modified peptide forms based on all possible combinations. More... | |
std::vector< OpenMS::AASequence > | generateTheoreticalPeptidoforms_ (const OpenMS::AASequence &sequence) |
Generate alternative modified peptide forms according to ModificationsDB. More... | |
std::vector< OpenMS::AASequence > | generateTheoreticalPeptidoformsDecoy_ (const OpenMS::AASequence &sequence, const OpenMS::AASequence &decoy_sequence) |
Generate alternative modified peptide forms according to ModificationsDB. More... | |
void | generateTargetInSilicoMap_ (const OpenMS::TargetedExperiment &exp, const std::vector< String > &fragment_types, const std::vector< size_t > &fragment_charges, bool enable_specific_losses, bool enable_unspecific_losses, bool enable_ms2_precursors, const std::vector< std::pair< double, double > > &swathes, int round_decPow, size_t max_num_alternative_localizations, SequenceMapT &TargetSequenceMap, IonMapT &TargetIonMap, PeptideMapT &TargetPeptideMap) |
Generate target in silico map. More... | |
void | generateDecoySequences_ (const SequenceMapT &TargetSequenceMap, std::map< String, String > &DecoySequenceMap, int shuffle_seed) |
Generate decoy sequences. More... | |
void | generateDecoyInSilicoMap_ (const OpenMS::TargetedExperiment &exp, const std::vector< String > &fragment_types, const std::vector< size_t > &fragment_charges, bool enable_specific_losses, bool enable_unspecific_losses, bool enable_ms2_precursors, const std::vector< std::pair< double, double > > &swathes, int round_decPow, TargetDecoyMapT &TargetDecoyMap, PeptideMapT &TargetPeptideMap, std::map< String, String > &DecoySequenceMap, IonMapT &DecoyIonMap, PeptideMapT &DecoyPeptideMap) |
Generate decoy in silico map. More... | |
void | generateTargetAssays_ (const OpenMS::TargetedExperiment &exp, TransitionVectorType &transitions, double mz_threshold, const std::vector< std::pair< double, double > > &swathes, int round_decPow, const PeptideMapT &TargetPeptideMap, const IonMapT &TargetIonMap) |
Generate target identification transitions. More... | |
void | generateDecoyAssays_ (const OpenMS::TargetedExperiment &exp, TransitionVectorType &transitions, double mz_threshold, const std::vector< std::pair< double, double > > &swathes, int round_decPow, const PeptideMapT &DecoyPeptideMap, TargetDecoyMapT &TargetDecoyMap, const IonMapT &DecoyIonMap, const IonMapT &TargetIonMap) |
Generate decoy assays. More... | |
Additional Inherited Members | |
Protected Attributes inherited from ProgressLogger | |
LogType | type_ |
time_t | last_invoke_ |
ProgressLoggerImpl * | current_logger_ |
Static Protected Attributes inherited from ProgressLogger | |
static int | recursion_depth_ |
Generate assays from a TargetedExperiment.
Will generate assays from a raw, unfiltered TargetedExperiment, as can be generated by TargetedFileConverter.
Transitions can be selected according to a set of rules, as described in Schubert et al., 2015 (PMID 25675208).
In addition, unique ion signature (UIS) (Sherman et al., 2009; PMID 19556279) transitions can be generated based on empirically observed or in silico generated ion series. This is described in detail in the IPF paper (Rosenberger et al. 2017; PMID 28604659).
typedef std::map<String, std::vector<const ReactionMonitoringTransition*> > CompoundTransitionMapType |
typedef std::vector<OpenMS::TargetedExperiment::Compound> CompoundVectorType |
typedef std::vector<std::pair<double, std::string> > FragmentSeqMap |
Describes a fragment sequence map of : "fragment m/z" -> "modified sequence".
typedef boost::unordered_map<size_t, boost::unordered_map<String, FragmentSeqMap > > IonMapT |
Stores a mapping : "unmodified sequence" -> FragmentSeqMap for all SWATH windows.
typedef std::vector<std::pair<std::string, double> > IonSeries |
Describes an ion series: "ion_type" -> "fragment m/z".
typedef std::map<String, std::set<std::string> > ModifiedSequenceMap |
Maps an unmodified sequence to all its modified sequences.
typedef std::map<String, IonSeries > PeptideMapT |
Maps a peptide sequence to an ion series: "ion_type" -> "fragment m/z".
typedef std::map<String, std::vector<const ReactionMonitoringTransition*> > PeptideTransitionMapType |
typedef std::vector<OpenMS::TargetedExperiment::Peptide> PeptideVectorType |
typedef std::vector<OpenMS::TargetedExperiment::Protein> ProteinVectorType |
typedef boost::unordered_map<size_t, ModifiedSequenceMap> SequenceMapT |
Stores the ModifiedSequenceMap for all SWATH windows.
typedef std::map<String, TargetedExperiment::Peptide> TargetDecoyMapT |
Maps the peptide id (same for target and decoy) to the decoy peptide object.
typedef std::vector<OpenMS::ReactionMonitoringTransition> TransitionVectorType |
MRMAssay | ( | ) |
Constructor.
|
override |
Destructor.
|
protected |
Generate modified peptide forms based on all possible combinations.
sequences | template AASequences |
mods_combs | all possible combinations (e.g. from nchoosekcombinations() ) |
modification | String of the modification |
void detectingTransitions | ( | OpenMS::TargetedExperiment & | exp, |
int | min_transitions, | ||
int | max_transitions | ||
) |
Select detecting fragment ions.
exp | the input, unfiltered transitions |
min_transitions | the minimum number of transitions required per assay |
max_transitions | the maximum number of transitions required per assay |
void filterMinMaxTransitionsCompound | ( | OpenMS::TargetedExperiment & | exp, |
int | min_transitions, | ||
int | max_transitions | ||
) |
Filters target and decoy transitions by intensity, only keeping the top N transitions.
exp | the transition list which will be filtered |
min_transitions | the minimum number of transitions required per assay (targets only) |
max_transitions | the maximum number of transitions allowed per assay |
void filterUnreferencedDecoysCompound | ( | OpenMS::TargetedExperiment & | exp | ) |
Filters decoy transitions, which do not have respective target transition based on the transitionID.
References between targets and decoys will be constructed based on the transitionID and the "_decoy_" string. For example:
target: 84_CompoundName_[M+H]+_88_22 decoy: 84_CompoundName_decoy_[M+H]+_88_22
exp | the transition list which will be filtered |
|
protected |
Generate decoy assays.
Used internally by the IPF algorithm, see MRMAssay::uisTransitions()
|
protected |
Generate decoy in silico map.
For each decoy peptide sequence, compute all alternative peptidoforms using all modification-carrying residue permutations (n choose k possibilities) that are physicochemically possible according to ModificationsDB.
Used internally by the IPF algorithm, see MRMAssay::uisTransitions()
|
protected |
Generate decoy sequences.
Generate decoy sequences for IPF algorithm which share peptidoform properties with targets.
Used internally by the IPF algorithm, see MRMAssay::uisTransitions()
|
protected |
Generate target identification transitions.
This function iterates over each (theoretical) transition of each target peptide and records the identity of all peptidoforms that map to each transition. The resulting transitions are stored in transitions.
exp | The input experiment with the target peptides |
transitions | The output containing annotated transitions with potential interferences |
mz_threshold | The threshold for annotating transitions as equal |
swathes | The swath windows used |
round_decPow | round product m/z values to decimal power (default: -4) |
TargetPeptideMap | Theoretical transitions for each peptide generated before |
TargetIonMap | Theoretical transitions for each peptide generated before |
Used internally by the IPF algorithm, see MRMAssay::uisTransitions()
|
protected |
Generate target in silico map.
For each peptide, compute all alternative peptidoforms using all modification-carrying residue permutations (n choose k possibilities) that are physicochemically possible according to ModificationsDB.
Then store the ion series for each of these theoretical fragment ions in the provided maps TargetSequenceMap, TargetIonMap, TargetPeptideMap.
Used internally by the IPF algorithm, see MRMAssay::uisTransitions()
|
protected |
Generate alternative modified peptide forms according to ModificationsDB.
An input peptide sequence containing modifications is used as template to generate all modification-carrying residue permutations (n choose k possibilities) that are physicochemically possible according to ModificationsDB.
sequence | template AASequence |
|
protected |
Generate alternative modified peptide forms according to ModificationsDB.
An input peptide sequence containing modifications is used as template to generate all modification-carrying residue permutations (n choose k possibilities) that are physicochemically possible according to ModificationsDB. Instead of the target sequence, the permutations are transferred to the decoy sequence that might contain additional modifiable residues. E.g. target sequence SAS(Phospho)K could result in [SAS(Phospho)K, S(Phospho)ASK] but the responding set of the decoy sequence SSS(Phospho)K would be [SSS(Phospho)K, S(Phospho)SSK].
sequence | template AASequence |
decoy_sequence | template decoy AASequence |
|
protected |
Check whether fragment ion are unique ion signatures in vector within threshold and return matching peptidoforms.
fragment_ion | the queried fragment ion |
ions | a vector of pairs of fragment ion m/z and peptide sequences which could interfere with fragment_ion |
mz_threshold | the threshold within which to search for interferences |
|
protected |
Generates random peptide sequence.
sequence_size | length of peptide sequence |
pseudoRNG | a Boost pseudo RNG |
|
protected |
Get swath index (precursor isolation window ordinal) for a particular precursor.
swathes | the swath window settings |
precursor_mz | the query precursor m/z |
|
protected |
Check whether the product m/z of a transition falls into the precursor isolation window.
swathes | the swath window settings |
precursor_mz | the query precursor m/z |
product_mz | the query product m/z |
|
protected |
Computes all N choose K combinations.
n | vector of N indices |
k | number of K |
void reannotateTransitions | ( | OpenMS::TargetedExperiment & | exp, |
double | precursor_mz_threshold, | ||
double | product_mz_threshold, | ||
const std::vector< String > & | fragment_types, | ||
const std::vector< size_t > & | fragment_charges, | ||
bool | enable_specific_losses, | ||
bool | enable_unspecific_losses, | ||
int | round_decPow = -4 |
||
) |
Annotates and filters transitions in a TargetedExperiment.
exp | the input, unfiltered transitions |
precursor_mz_threshold | the precursor m/z threshold in Th for annotation |
product_mz_threshold | the product m/z threshold in Th for annotation |
fragment_types | the fragment types to consider for annotation |
fragment_charges | the fragment charges to consider for annotation |
enable_specific_losses | whether specific neutral losses should be considered |
enable_unspecific_losses | whether unspecific neutral losses (H2O1, H3N1, C1H2N2, C1H2N1O1) should be considered |
round_decPow | round product m/z values to decimal power (default: -4) |
void restrictTransitions | ( | OpenMS::TargetedExperiment & | exp, |
double | lower_mz_limit, | ||
double | upper_mz_limit, | ||
const std::vector< std::pair< double, double > > & | swathes | ||
) |
Restrict and filter transitions in a TargetedExperiment.
exp | the input, unfiltered transitions |
lower_mz_limit | the lower product m/z limit in Th |
upper_mz_limit | the upper product m/z limit in Th |
swathes | the swath window settings (to exclude fragment ions falling into the precursor isolation window) |
void uisTransitions | ( | OpenMS::TargetedExperiment & | exp, |
const std::vector< String > & | fragment_types, | ||
const std::vector< size_t > & | fragment_charges, | ||
bool | enable_specific_losses, | ||
bool | enable_unspecific_losses, | ||
bool | enable_ms2_precursors, | ||
double | mz_threshold, | ||
const std::vector< std::pair< double, double > > & | swathes, | ||
int | round_decPow = -4 , |
||
size_t | max_num_alternative_localizations = 20 , |
||
int | shuffle_seed = -1 , |
||
bool | disable_decoy_transitions = false |
||
) |
Annotate UIS / site-specific transitions.
Performs the following actions:
The IPF algorithm uses the concept of "identification transitions" that are used to discriminate different peptidoforms, these are generated in this function. In brief, the algorithm takes the existing set of peptides and transitions and then appends these "identification transitions" for targets and decoys. The novel transitions are set to be non-detecting and non-quantifying and are annotated with the set of peptidoforms to which they map.
exp | the input, unfiltered transitions |
fragment_types | the fragment types to consider for annotation |
fragment_charges | the fragment charges to consider for annotation |
enable_specific_losses | whether specific neutral losses should be considered |
enable_unspecific_losses | whether unspecific neutral losses (H2O1, H3N1, C1H2N2, C1H2N1O1) should be considered |
enable_ms2_precursors | whether MS2 precursors should be considered |
mz_threshold | the product m/z threshold in Th for annotation |
swathes | the swath window settings (to exclude fragment ions falling |
round_decPow | round product m/z values to decimal power (default: -4) |
max_num_alternative_localizations | maximum number of allowed peptide sequence permutations |
shuffle_seed | set seed for shuffle (-1: select seed based on time) |
disable_decoy_transitions | whether to disable generation of decoy UIS transitions |