OpenMS
MRMAssay Class Reference

Generate assays from a TargetedExperiment. More...

#include <OpenMS/ANALYSIS/OPENSWATH/MRMAssay.h>

Inheritance diagram for MRMAssay:
[legend]
Collaboration diagram for MRMAssay:
[legend]

Public Types

typedef std::vector< OpenMS::TargetedExperiment::ProteinProteinVectorType
 
typedef std::vector< OpenMS::TargetedExperiment::PeptidePeptideVectorType
 
typedef std::vector< OpenMS::TargetedExperiment::CompoundCompoundVectorType
 
typedef std::vector< OpenMS::ReactionMonitoringTransitionTransitionVectorType
 
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, ModifiedSequenceMapSequenceMapT
 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, IonSeriesPeptideMapT
 Maps a peptide sequence to an ion series: "ion_type" -> "fragment m/z". More...
 
typedef std::map< String, TargetedExperiment::PeptideTargetDecoyMapT
 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...
 
ProgressLoggeroperator= (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::AASequenceaddModificationsSequences_ (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::AASequencegenerateTheoreticalPeptidoforms_ (const OpenMS::AASequence &sequence)
 Generate alternative modified peptide forms according to ModificationsDB. More...
 
std::vector< OpenMS::AASequencegenerateTheoreticalPeptidoformsDecoy_ (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_
 
ProgressLoggerImplcurrent_logger_
 
- Static Protected Attributes inherited from ProgressLogger
static int recursion_depth_
 

Detailed Description

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).

Member Typedef Documentation

◆ CompoundTransitionMapType

typedef std::map<String, std::vector<const ReactionMonitoringTransition*> > CompoundTransitionMapType

◆ CompoundVectorType

◆ FragmentSeqMap

typedef std::vector<std::pair<double, std::string> > FragmentSeqMap

Describes a fragment sequence map of : "fragment m/z" -> "modified sequence".

◆ IonMapT

typedef boost::unordered_map<size_t, boost::unordered_map<String, FragmentSeqMap > > IonMapT

Stores a mapping : "unmodified sequence" -> FragmentSeqMap for all SWATH windows.

◆ IonSeries

typedef std::vector<std::pair<std::string, double> > IonSeries

Describes an ion series: "ion_type" -> "fragment m/z".

◆ ModifiedSequenceMap

typedef std::map<String, std::set<std::string> > ModifiedSequenceMap

Maps an unmodified sequence to all its modified sequences.

◆ PeptideMapT

typedef std::map<String, IonSeries > PeptideMapT

Maps a peptide sequence to an ion series: "ion_type" -> "fragment m/z".

◆ PeptideTransitionMapType

typedef std::map<String, std::vector<const ReactionMonitoringTransition*> > PeptideTransitionMapType

◆ PeptideVectorType

◆ ProteinVectorType

◆ SequenceMapT

typedef boost::unordered_map<size_t, ModifiedSequenceMap> SequenceMapT

Stores the ModifiedSequenceMap for all SWATH windows.

◆ TargetDecoyMapT

Maps the peptide id (same for target and decoy) to the decoy peptide object.

◆ TransitionVectorType

Constructor & Destructor Documentation

◆ MRMAssay()

MRMAssay ( )

Constructor.

◆ ~MRMAssay()

~MRMAssay ( )
override

Destructor.

Member Function Documentation

◆ addModificationsSequences_()

std::vector<OpenMS::AASequence> addModificationsSequences_ ( const std::vector< OpenMS::AASequence > &  sequences,
const std::vector< std::vector< size_t > > &  mods_combs,
const OpenMS::String modification 
)
protected

Generate modified peptide forms based on all possible combinations.

Parameters
sequencestemplate AASequences
mods_combsall possible combinations (e.g. from nchoosekcombinations() )
modificationString of the modification
Returns
a vector of all modified peptides.

◆ detectingTransitions()

void detectingTransitions ( OpenMS::TargetedExperiment exp,
int  min_transitions,
int  max_transitions 
)

Select detecting fragment ions.

Parameters
expthe input, unfiltered transitions
min_transitionsthe minimum number of transitions required per assay
max_transitionsthe maximum number of transitions required per assay

◆ filterMinMaxTransitionsCompound()

void filterMinMaxTransitionsCompound ( OpenMS::TargetedExperiment exp,
int  min_transitions,
int  max_transitions 
)

Filters target and decoy transitions by intensity, only keeping the top N transitions.

Parameters
expthe transition list which will be filtered
min_transitionsthe minimum number of transitions required per assay (targets only)
max_transitionsthe maximum number of transitions allowed per assay

◆ filterUnreferencedDecoysCompound()

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

Parameters
expthe transition list which will be filtered

◆ generateDecoyAssays_()

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 
)
protected

Generate decoy assays.

Used internally by the IPF algorithm, see MRMAssay::uisTransitions()

◆ generateDecoyInSilicoMap_()

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 
)
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()

◆ generateDecoySequences_()

void generateDecoySequences_ ( const SequenceMapT TargetSequenceMap,
std::map< String, String > &  DecoySequenceMap,
int  shuffle_seed 
)
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()

◆ generateTargetAssays_()

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 
)
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.

Parameters
expThe input experiment with the target peptides
transitionsThe output containing annotated transitions with potential interferences
mz_thresholdThe threshold for annotating transitions as equal
swathesThe swath windows used
round_decPowround product m/z values to decimal power (default: -4)
TargetPeptideMapTheoretical transitions for each peptide generated before
TargetIonMapTheoretical transitions for each peptide generated before

Used internally by the IPF algorithm, see MRMAssay::uisTransitions()

◆ generateTargetInSilicoMap_()

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 
)
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()

◆ generateTheoreticalPeptidoforms_()

std::vector<OpenMS::AASequence> generateTheoreticalPeptidoforms_ ( const OpenMS::AASequence sequence)
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.

Parameters
sequencetemplate AASequence
Returns
a vector of all alternative modified peptides.

◆ generateTheoreticalPeptidoformsDecoy_()

std::vector<OpenMS::AASequence> generateTheoreticalPeptidoformsDecoy_ ( const OpenMS::AASequence sequence,
const OpenMS::AASequence decoy_sequence 
)
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].

Parameters
sequencetemplate AASequence
decoy_sequencetemplate decoy AASequence
Returns
a vector of all alternative modified peptides.

◆ getMatchingPeptidoforms_()

std::vector<std::string> getMatchingPeptidoforms_ ( const double  fragment_ion,
const FragmentSeqMap ions,
const double  mz_threshold 
)
protected

Check whether fragment ion are unique ion signatures in vector within threshold and return matching peptidoforms.

Parameters
fragment_ionthe queried fragment ion
ionsa vector of pairs of fragment ion m/z and peptide sequences which could interfere with fragment_ion
mz_thresholdthe threshold within which to search for interferences
Returns
a vector of strings containing all peptidoforms with which fragment_ion overlaps

◆ getRandomSequence_()

std::string getRandomSequence_ ( size_t  sequence_size,
boost::variate_generator< boost::mt19937 &, boost::uniform_int<> >  pseudoRNG 
)
protected

Generates random peptide sequence.

Parameters
sequence_sizelength of peptide sequence
pseudoRNGa Boost pseudo RNG
Returns
random peptide sequence

◆ getSwath_()

int getSwath_ ( const std::vector< std::pair< double, double > > &  swathes,
const double  precursor_mz 
)
protected

Get swath index (precursor isolation window ordinal) for a particular precursor.

Parameters
swathesthe swath window settings
precursor_mzthe query precursor m/z
Returns
index of swath where precursor_mz falls into

◆ isInSwath_()

bool isInSwath_ ( const std::vector< std::pair< double, double > > &  swathes,
const double  precursor_mz,
const double  product_mz 
)
protected

Check whether the product m/z of a transition falls into the precursor isolation window.

Parameters
swathesthe swath window settings
precursor_mzthe query precursor m/z
product_mzthe query product m/z
Returns
whether product m/z falls into precursor isolation window

◆ nchoosekcombinations_()

std::vector<std::vector<size_t> > nchoosekcombinations_ ( const std::vector< size_t > &  n,
size_t  k 
)
protected

Computes all N choose K combinations.

Parameters
nvector of N indices
knumber of K
Returns
a vector of all N index combinations

◆ reannotateTransitions()

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.

Parameters
expthe input, unfiltered transitions
precursor_mz_thresholdthe precursor m/z threshold in Th for annotation
product_mz_thresholdthe product m/z threshold in Th for annotation
fragment_typesthe fragment types to consider for annotation
fragment_chargesthe fragment charges to consider for annotation
enable_specific_losseswhether specific neutral losses should be considered
enable_unspecific_losseswhether unspecific neutral losses (H2O1, H3N1, C1H2N2, C1H2N1O1) should be considered
round_decPowround product m/z values to decimal power (default: -4)

◆ restrictTransitions()

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.

Parameters
expthe input, unfiltered transitions
lower_mz_limitthe lower product m/z limit in Th
upper_mz_limitthe upper product m/z limit in Th
swathesthe swath window settings (to exclude fragment ions falling into the precursor isolation window)

◆ uisTransitions()

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.

Parameters
expthe input, unfiltered transitions
fragment_typesthe fragment types to consider for annotation
fragment_chargesthe fragment charges to consider for annotation
enable_specific_losseswhether specific neutral losses should be considered
enable_unspecific_losseswhether unspecific neutral losses (H2O1, H3N1, C1H2N2, C1H2N1O1) should be considered
enable_ms2_precursorswhether MS2 precursors should be considered
mz_thresholdthe product m/z threshold in Th for annotation
swathesthe swath window settings (to exclude fragment ions falling
round_decPowround product m/z values to decimal power (default: -4)
max_num_alternative_localizationsmaximum number of allowed peptide sequence permutations
shuffle_seedset seed for shuffle (-1: select seed based on time)
disable_decoy_transitionswhether to disable generation of decoy UIS transitions