OpenMS
MRMFeatureFinderScoring Class Reference

The MRMFeatureFinder finds and scores peaks of transitions that co-elute. More...

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

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

Public Types

typedef OpenSwath::LightTransition TransitionType
 Type definitions. More...
 
typedef OpenSwath::LightTargetedExperiment TargetedExpType
 
typedef OpenSwath::LightCompound PeptideType
 
typedef OpenSwath::LightProtein ProteinType
 
typedef OpenSwath::LightModification ModificationType
 
typedef MRMTransitionGroup< MSChromatogram, TransitionTypeMRMTransitionGroupType
 
typedef std::map< String, MRMTransitionGroupTypeTransitionGroupMapType
 
- Public Types inherited from ProgressLogger
enum  LogType { CMD , GUI , NONE }
 Possible log types. More...
 

Public Member Functions

 MRMFeatureFinderScoring ()
 Constructor. More...
 
 ~MRMFeatureFinderScoring () override
 Destructor. More...
 
void pickExperiment (const PeakMap &chromatograms, FeatureMap &output, const TargetedExperiment &transition_exp, const TransformationDescription &trafo, const PeakMap &swath_map)
 Picker and prepare functions. More...
 
void pickExperiment (const OpenSwath::SpectrumAccessPtr &input, FeatureMap &output, const OpenSwath::LightTargetedExperiment &transition_exp, const TransformationDescription &trafo, const std::vector< OpenSwath::SwathMap > &swath_maps, TransitionGroupMapType &transition_group_map)
 Pick and score features in a single experiment from chromatograms. More...
 
void prepareProteinPeptideMaps_ (const OpenSwath::LightTargetedExperiment &transition_exp)
 Prepares the internal mappings of peptides and proteins. More...
 
void scorePeakgroups (MRMTransitionGroupType &transition_group, const TransformationDescription &trafo, const std::vector< OpenSwath::SwathMap > &swath_maps, FeatureMap &output, bool ms1only=false) const
 Score all peak groups of a transition group. More...
 
void setStrictFlag (bool f)
 Set the flag for strict mapping. More...
 
void setMS1Map (OpenSwath::SpectrumAccessPtr ms1_map)
 Add an MS1 map containing spectra. More...
 
void mapExperimentToTransitionList (const OpenSwath::SpectrumAccessPtr &input, const OpenSwath::LightTargetedExperiment &transition_exp, TransitionGroupMapType &transition_group_map, TransformationDescription trafo, double rt_extraction_window)
 Map the chromatograms to the transitions. More...
 
- Public Member Functions inherited from DefaultParamHandler
 DefaultParamHandler (const String &name)
 Constructor with name that is displayed in error messages. More...
 
 DefaultParamHandler (const DefaultParamHandler &rhs)
 Copy constructor. More...
 
virtual ~DefaultParamHandler ()
 Destructor. More...
 
DefaultParamHandleroperator= (const DefaultParamHandler &rhs)
 Assignment operator. More...
 
virtual bool operator== (const DefaultParamHandler &rhs) const
 Equality operator. More...
 
void setParameters (const Param &param)
 Sets the parameters. More...
 
const ParamgetParameters () const
 Non-mutable access to the parameters. More...
 
const ParamgetDefaults () const
 Non-mutable access to the default parameters. More...
 
const StringgetName () const
 Non-mutable access to the name. More...
 
void setName (const String &name)
 Mutable access to the name. More...
 
const std::vector< String > & getSubsections () const
 Non-mutable access to the registered subsections. 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...
 

Private Member Functions

void splitTransitionGroupsDetection_ (const MRMTransitionGroupType &transition_group, MRMTransitionGroupType &transition_group_detection) const
 Splits combined transition groups into detection transition groups. More...
 
void splitTransitionGroupsIdentification_ (const MRMTransitionGroupType &transition_group, MRMTransitionGroupType &transition_group_identification, MRMTransitionGroupType &transition_group_identification_decoy) const
 Splits combined transition groups into identification transition groups. More...
 
OpenSwath_Ind_Scores scoreIdentification_ (MRMTransitionGroupType &transition_group_identification, OpenSwathScoring &scorer, const size_t feature_idx, const std::vector< std::string > &native_ids_detection, const double det_intensity_ratio_score, const double det_mi_ratio_score, const std::vector< OpenSwath::SwathMap > &swath_maps) const
 Provides scoring for target and decoy identification against detecting transitions. More...
 
void prepareFeatureOutput_ (OpenMS::MRMFeature &mrmfeature, bool ms1only, int charge) const
 
void updateMembers_ () override
 Synchronize members with param class. More...
 

Private Attributes

double rt_extraction_window_
 
double quantification_cutoff_
 
int stop_report_after_feature_
 
bool write_convex_hull_
 
bool strict_
 
bool use_ms1_ion_mobility_
 
String scoring_model_
 
double rt_normalization_factor_
 
int add_up_spectra_
 
String spectrum_addition_method_
 
double spacing_for_spectra_resampling_
 
double uis_threshold_sn_
 
double uis_threshold_peak_area_
 
double sn_win_len_
 
unsigned int sn_bin_count_
 
bool write_log_messages_
 
double im_extra_drift_
 
std::map< OpenMS::String, const PeptideType * > PeptideRefMap_
 
OpenSwath_Scores_Usage su_
 
OpenMS::DIAScoring diascoring_
 
OpenMS::SONARScoring sonarscoring_
 
OpenMS::EmgScoring emgscoring_
 
OpenSwath::SpectrumAccessPtr ms1_map_
 

Additional Inherited Members

- Static Public Member Functions inherited from DefaultParamHandler
static void writeParametersToMetaValues (const Param &write_this, MetaInfoInterface &write_here, const String &key_prefix="")
 Writes all parameters to meta values. More...
 
- Protected Member Functions inherited from DefaultParamHandler
void defaultsToParam_ ()
 Updates the parameters after the defaults have been set in the constructor. More...
 
- Protected Attributes inherited from DefaultParamHandler
Param param_
 Container for current parameters. More...
 
Param defaults_
 Container for default parameters. This member should be filled in the constructor of derived classes! More...
 
std::vector< Stringsubsections_
 Container for registered subsections. This member should be filled in the constructor of derived classes! More...
 
String error_name_
 Name that is displayed in error messages during the parameter checking. More...
 
bool check_defaults_
 If this member is set to false no checking if parameters in done;. More...
 
bool warn_empty_defaults_
 If this member is set to false no warning is emitted when defaults are empty;. More...
 
- 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

The MRMFeatureFinder finds and scores peaks of transitions that co-elute.

It does so using an internal peakpicker for each chromatogram and then creating consensus / meta-peaks (MRMFeatures) that contain the information of all corresponding chromatograms at the peak-position. It then goes on to score those MRMFeatures using different criteria described in the MRMScoring class.

Internally, all peak group detection is performed in MRMTransitionGroupPicker which segments the data and determines consensus peaks across traces (MRMFeatures). All scoring is delegated to the OpenSwathScoring class which implements i) chromatographic scores, ii) library based scores and iii) full spectrum (DIA) scores. These scores are retrieved from the OpenSwathScoring class and added to the MRMFeatures found in this algorithm. Note that the OpenSwathScoring is a facade that can be used to communicate with the underlying actual scoring engines and assembles the scores inside a scoring object called OpenSwath_Scores where they are easy to retrieve.

Parameters of this class are:

NameTypeDefaultRestrictionsDescription
stop_report_after_feature int-1  Stop reporting after feature (ordered by quality; -1 means do not stop).
rt_extraction_window float-1.0  Only extract RT around this value (-1 means extract over the whole range, a value of 500 means to extract around +/- 500 s of the expected elution). For this to work, the TraML input file needs to contain normalized RT values.
rt_normalization_factor float1.0  The normalized RT is expected to be between 0 and 1. If your normalized RT has a different range, pass this here (e.g. it goes from 0 to 100, set this value to 100)
quantification_cutoff float0.0 min: 0.0Cutoff in m/z below which peaks should not be used for quantification any more
write_convex_hull stringfalse true, falseWhether to write out all points of all features into the featureXML
spectrum_addition_method stringsimple simple, resampleFor spectrum addition, either use simple concatenation or use peak resampling
add_up_spectra int1 min: 1Add up spectra around the peak apex (needs to be a non-even integer)
spacing_for_spectra_resampling float5.0e-03 min: 0.0If spectra are to be added, use this spacing to add them up
uis_threshold_sn int-1  S/N threshold to consider identification transition (set to -1 to consider all)
uis_threshold_peak_area int0  Peak area threshold to consider identification transition (set to -1 to consider all)
scoring_model stringdefault default, single_transitionScoring model to use
im_extra_drift float0.0 min: 0.0Extra drift time to extract for IM scoring (as a fraction, e.g. 0.25 means 25% extra on each side)
strict stringtrue true, falseWhether to error (true) or skip (false) if a transition in a transition group does not have a corresponding chromatogram.
use_ms1_ion_mobility stringtrue  Performs ion mobility extraction in MS1. Set to false if MS1 spectra do not contain ion mobility
TransitionGroupPicker:stop_after_feature int-1  Stop finding after feature (ordered by intensity; -1 means do not stop).
TransitionGroupPicker:stop_after_intensity_ratio float1.0e-04  Stop after reaching intensity ratio
TransitionGroupPicker:min_peak_width float1.0e-03  Minimal peak width (s), discard all peaks below this value (-1 means no action).
TransitionGroupPicker:peak_integration stringoriginal original, smoothedCalculate the peak area and height either the smoothed or the raw chromatogram data.
TransitionGroupPicker:background_subtraction stringnone none, original, exactRemove background from peak signal using estimated noise levels. The 'original' method is only provided for historical purposes, please use the 'exact' method and set parameters using the PeakIntegrator: settings. The same original or smoothed chromatogram specified by peak_integration will be used for background estimation.
TransitionGroupPicker:recalculate_peaks stringfalse true, falseTries to get better peak picking by looking at peak consistency of all picked peaks. Tries to use the consensus (median) peak border if the variation within the picked peaks is too large.
TransitionGroupPicker:use_precursors stringfalse true, falseUse precursor chromatogram for peak picking (note that this may lead to precursor signal driving the peak picking)
TransitionGroupPicker:use_consensus stringtrue true, falseUse consensus peak boundaries when computing transition group picking (if false, compute independent peak boundaries for each transition)
TransitionGroupPicker:recalculate_peaks_max_z float1.0  Determines the maximal Z-Score (difference measured in standard deviations) that is considered too large for peak boundaries. If the Z-Score is above this value, the median is used for peak boundaries (default value 1.0).
TransitionGroupPicker:minimal_quality float-1.0e04  Only if compute_peak_quality is set, this parameter will not consider peaks below this quality threshold
TransitionGroupPicker:resample_boundary float15.0  For computing peak quality, how many extra seconds should be sample left and right of the actual peak
TransitionGroupPicker:compute_peak_quality stringfalse true, falseTries to compute a quality value for each peakgroup and detect outlier transitions. The resulting score is centered around zero and values above 0 are generally good and below -1 or -2 are usually bad.
TransitionGroupPicker:compute_peak_shape_metrics stringfalse true, falseCalculates various peak shape metrics (e.g., tailing) that can be used for downstream QC/QA.
TransitionGroupPicker:compute_total_mi stringfalse true, falseCompute mutual information metrics for individual transitions that can be used for OpenSWATH/IPF scoring.
TransitionGroupPicker:boundary_selection_method stringlargest largest, widestMethod to use when selecting the best boundaries for peaks.
TransitionGroupPicker:PeakPickerChromatogram:sgolay_frame_length int15  The number of subsequent data points used for smoothing.
This number has to be uneven. If it is not, 1 will be added.
TransitionGroupPicker:PeakPickerChromatogram:sgolay_polynomial_order int3  Order of the polynomial that is fitted.
TransitionGroupPicker:PeakPickerChromatogram:gauss_width float50.0  Gaussian width in seconds, estimated peak size.
TransitionGroupPicker:PeakPickerChromatogram:use_gauss stringtrue false, trueUse Gaussian filter for smoothing (alternative is Savitzky-Golay filter)
TransitionGroupPicker:PeakPickerChromatogram:peak_width float-1.0  Force a certain minimal peak_width on the data (e.g. extend the peak at least by this amount on both sides) in seconds. -1 turns this feature off.
TransitionGroupPicker:PeakPickerChromatogram:signal_to_noise float1.0 min: 0.0Signal-to-noise threshold at which a peak will not be extended any more. Note that setting this too high (e.g. 1.0) can lead to peaks whose flanks are not fully captured.
TransitionGroupPicker:PeakPickerChromatogram:sn_win_len float1000.0  Signal to noise window length.
TransitionGroupPicker:PeakPickerChromatogram:sn_bin_count int30  Signal to noise bin count.
TransitionGroupPicker:PeakPickerChromatogram:write_sn_log_messages stringfalse true, falseWrite out log messages of the signal-to-noise estimator in case of sparse windows or median in rightmost histogram bin
TransitionGroupPicker:PeakPickerChromatogram:remove_overlapping_peaks stringfalse false, trueTry to remove overlapping peaks during peak picking
TransitionGroupPicker:PeakPickerChromatogram:method stringcorrected legacy, corrected, crawdadWhich method to choose for chromatographic peak-picking (OpenSWATH legacy on raw data, corrected picking on smoothed chromatogram or Crawdad on smoothed chromatogram).
TransitionGroupPicker:PeakIntegrator:integration_type stringintensity_sum intensity_sum, simpson, trapezoidThe integration technique to use in integratePeak() and estimateBackground() which uses either the summed intensity, integration by Simpson's rule or trapezoidal integration.
TransitionGroupPicker:PeakIntegrator:baseline_type stringbase_to_base base_to_base, vertical_division, vertical_division_min, vertical_division_maxThe baseline type to use in estimateBackground() based on the peak boundaries. A rectangular baseline shape is computed based either on the minimal intensity of the peak boundaries, the maximum intensity or the average intensity (base_to_base).
TransitionGroupPicker:PeakIntegrator:fit_EMG stringfalse false, trueFit the chromatogram/spectrum to the EMG peak model.
DIAScoring:dia_extraction_window float0.05 min: 0.0DIA extraction window in Th or ppm.
DIAScoring:dia_extraction_unit stringTh Th, ppmDIA extraction window unit
DIAScoring:dia_centroided stringfalse true, falseUse centroided DIA data.
DIAScoring:dia_byseries_intensity_min float300.0 min: 0.0DIA b/y series minimum intensity to consider.
DIAScoring:dia_byseries_ppm_diff float10.0 min: 0.0DIA b/y series minimal difference in ppm to consider.
DIAScoring:dia_nr_isotopes int4 min: 0DIA number of isotopes to consider.
DIAScoring:dia_nr_charges int4 min: 0DIA number of charges to consider.
DIAScoring:peak_before_mono_max_ppm_diff float20.0 min: 0.0DIA maximal difference in ppm to count a peak at lower m/z when searching for evidence that a peak might not be monoisotopic.
EMGScoring:interpolation_step float0.2  Sampling rate for the interpolation of the model function.
EMGScoring:tolerance_stdev_bounding_box float3.0  Bounding box has range [minimim of data, maximum of data] enlarged by tolerance_stdev_bounding_box times the standard deviation of the data.
EMGScoring:max_iteration int500  Maximum number of iterations using by Levenberg-Marquardt algorithm.
EMGScoring:init_mom stringfalse true, falseInitialize parameters using method of moments estimators.
EMGScoring:statistics:mean float1.0  Centroid position of the model.
EMGScoring:statistics:variance float1.0  Variance of the model.
Scores:use_shape_score stringtrue true, falseUse the shape score (this score measures the similarity in shape of the transitions using a cross-correlation)
Scores:use_coelution_score stringtrue true, falseUse the coelution score (this score measures the similarity in coelution of the transitions using a cross-correlation)
Scores:use_rt_score stringtrue true, falseUse the retention time score (this score measure the difference in retention time)
Scores:use_library_score stringtrue true, falseUse the library score
Scores:use_elution_model_score stringtrue true, falseUse the elution model (EMG) score (this score fits a gaussian model to the peak and checks the fit)
Scores:use_intensity_score stringtrue true, falseUse the intensity score
Scores:use_nr_peaks_score stringtrue true, falseUse the number of peaks score
Scores:use_total_xic_score stringtrue true, falseUse the total XIC score
Scores:use_total_mi_score stringfalse true, falseUse the total MI score
Scores:use_sn_score stringtrue true, falseUse the SN (signal to noise) score
Scores:use_mi_score stringfalse true, falseUse the MI (mutual information) score
Scores:use_dia_scores stringtrue true, falseUse the DIA (SWATH) scores. If turned off, will not use fragment ion spectra for scoring.
Scores:use_ms1_correlation stringfalse true, falseUse the correlation scores with the MS1 elution profiles
Scores:use_sonar_scores stringfalse true, falseUse the scores for SONAR scans (scanning swath)
Scores:use_ion_mobility_scores stringfalse true, falseUse the scores for Ion Mobility scans
Scores:use_ms1_fullscan stringfalse true, falseUse the full MS1 scan at the peak apex for scoring (ppm accuracy of precursor and isotopic pattern)
Scores:use_ms1_mi stringfalse true, falseUse the MS1 MI score
Scores:use_uis_scores stringfalse true, falseUse UIS scores for peptidoform identification
Scores:use_ionseries_scores stringtrue true, falseUse MS2-level b/y ion-series scores for peptidoform identification
Scores:use_ms2_isotope_scores stringtrue true, falseUse MS2-level isotope scores (pearson & manhattan) across product transitions (based on ID if annotated or averagine)

Note:
  • If a section name is documented, the documentation is displayed as tooltip.
  • Advanced parameter names are italic.

Member Typedef Documentation

◆ ModificationType

◆ MRMTransitionGroupType

◆ PeptideType

◆ ProteinType

◆ TargetedExpType

◆ TransitionGroupMapType

◆ TransitionType

Type definitions.

Constructor & Destructor Documentation

◆ MRMFeatureFinderScoring()

Constructor.

◆ ~MRMFeatureFinderScoring()

~MRMFeatureFinderScoring ( )
override

Destructor.

Member Function Documentation

◆ mapExperimentToTransitionList()

void mapExperimentToTransitionList ( const OpenSwath::SpectrumAccessPtr input,
const OpenSwath::LightTargetedExperiment transition_exp,
TransitionGroupMapType transition_group_map,
TransformationDescription  trafo,
double  rt_extraction_window 
)

Map the chromatograms to the transitions.

Map an input chromatogram experiment (mzML) and transition list (TraML) onto each other when they share identifiers, e.g. if the transition id is the same as the chromatogram native id.

Parameters
inputThe input chromatograms
transition_expThe transition list describing the experiment
transition_group_mapMapping of transition groups
trafoOptional transformation of the experimental retention time to the normalized retention time space used in the transition list.
rt_extraction_windowThe used retention time extraction window

◆ pickExperiment() [1/2]

void pickExperiment ( const OpenSwath::SpectrumAccessPtr input,
FeatureMap output,
const OpenSwath::LightTargetedExperiment transition_exp,
const TransformationDescription trafo,
const std::vector< OpenSwath::SwathMap > &  swath_maps,
TransitionGroupMapType transition_group_map 
)

Pick and score features in a single experiment from chromatograms.

Parameters
inputThe input chromatograms
outputThe output features with corresponding scores
transition_expThe transition list describing the experiment
trafoOptional transformation of the experimental retention time to the normalized retention time space used in the transition list.
swath_mapsOptional SWATH-MS (DIA) map corresponding from which the chromatograms were extracted. Use empty map if no data is available.
transition_group_mapOutput mapping of transition groups

◆ pickExperiment() [2/2]

void pickExperiment ( const PeakMap chromatograms,
FeatureMap output,
const TargetedExperiment transition_exp,
const TransformationDescription trafo,
const PeakMap swath_map 
)

Picker and prepare functions.

Pick and score features in a single experiment from chromatograms

Function for wrapping in Python, only uses OpenMS datastructures and does not return the map.

Parameters
chromatogramsThe input chromatograms
outputThe output features with corresponding scores
transition_expThe transition list describing the experiment
trafoOptional transformation of the experimental retention time to the normalized retention time space used in the transition list
swath_mapOptional SWATH-MS (DIA) map corresponding from which the chromatograms were extracted

◆ prepareFeatureOutput_()

void prepareFeatureOutput_ ( OpenMS::MRMFeature mrmfeature,
bool  ms1only,
int  charge 
) const
private

◆ prepareProteinPeptideMaps_()

void prepareProteinPeptideMaps_ ( const OpenSwath::LightTargetedExperiment transition_exp)

Prepares the internal mappings of peptides and proteins.

Calling this method is required before calling scorePeakgroups.

Parameters
transition_expThe transition list describing the experiment

◆ scoreIdentification_()

OpenSwath_Ind_Scores scoreIdentification_ ( MRMTransitionGroupType transition_group_identification,
OpenSwathScoring scorer,
const size_t  feature_idx,
const std::vector< std::string > &  native_ids_detection,
const double  det_intensity_ratio_score,
const double  det_mi_ratio_score,
const std::vector< OpenSwath::SwathMap > &  swath_maps 
) const
private

Provides scoring for target and decoy identification against detecting transitions.

The function is used twice, for target and decoy identification transitions. The results are reported analogously to the ones for detecting transitions but must be stored separately.

Parameters
transition_group_identificationContaining all detecting and identifying transitions
scorerAn instance of OpenSwathScoring
feature_idxThe index of the current feature
native_ids_detectionThe native IDs of the detecting transitions
det_intensity_ratio_scoreThe intensity score of the detection transitions for normalization
det_mi_ratio_scoreThe MI score of the detection transitions for normalization
swath_mapsOptional SWATH-MS (DIA) map corresponding from which the chromatograms were extracted. Use empty map if no data is available.
Returns
a struct of type OpenSwath_Ind_Scores containing either target or decoy values

◆ scorePeakgroups()

void scorePeakgroups ( MRMTransitionGroupType transition_group,
const TransformationDescription trafo,
const std::vector< OpenSwath::SwathMap > &  swath_maps,
FeatureMap output,
bool  ms1only = false 
) const

Score all peak groups of a transition group.

Iterate through all features found along the chromatograms of the transition group and score each one individually.

Parameters
transition_groupThe MRMTransitionGroup to be scored (input)
trafoOptional transformation of the experimental retention time to the normalized retention time space used in the transition list.
swath_mapsOptional SWATH-MS (DIA) map corresponding from which the chromatograms were extracted. Use empty map if no data is available.
outputThe output features with corresponding scores (the found features will be added to this FeatureMap).
ms1onlyWhether to only do MS1 scoring and skip all MS2 scoring

◆ setMS1Map()

void setMS1Map ( OpenSwath::SpectrumAccessPtr  ms1_map)
inline

Add an MS1 map containing spectra.

For DIA (SWATH-MS), an optional MS1 map can be supplied which can be used to extract precursor ion signal and provides additional scores. If no MS1 map is provided, the respective scores are not calculated.

Parameters
ms1_mapThe raw mass spectrometric MS1 data

◆ setStrictFlag()

void setStrictFlag ( bool  f)
inline

Set the flag for strict mapping.

◆ splitTransitionGroupsDetection_()

void splitTransitionGroupsDetection_ ( const MRMTransitionGroupType transition_group,
MRMTransitionGroupType transition_group_detection 
) const
private

Splits combined transition groups into detection transition groups.

For standard assays, transition_group_detection is identical to transition_group and the others are empty.

Parameters
transition_groupContaining all detecting, identifying transitions
transition_group_detectionTo be filled with detecting transitions

◆ splitTransitionGroupsIdentification_()

void splitTransitionGroupsIdentification_ ( const MRMTransitionGroupType transition_group,
MRMTransitionGroupType transition_group_identification,
MRMTransitionGroupType transition_group_identification_decoy 
) const
private

Splits combined transition groups into identification transition groups.

For standard assays, transition_group_identification is empty. When UIS scoring is enabled, it contains the corresponding identification transitions.

Parameters
transition_groupContaining all detecting, identifying transitions
transition_group_identificationTo be filled with identifying transitions
transition_group_identification_decoyTo be filled with identifying decoy transitions

◆ updateMembers_()

void updateMembers_ ( )
overrideprivatevirtual

Synchronize members with param class.

Reimplemented from DefaultParamHandler.

Member Data Documentation

◆ add_up_spectra_

int add_up_spectra_
private

◆ diascoring_

OpenMS::DIAScoring diascoring_
private

◆ emgscoring_

OpenMS::EmgScoring emgscoring_
private

◆ im_extra_drift_

double im_extra_drift_
private

◆ ms1_map_

OpenSwath::SpectrumAccessPtr ms1_map_
private

◆ PeptideRefMap_

std::map<OpenMS::String, const PeptideType*> PeptideRefMap_
private

◆ quantification_cutoff_

double quantification_cutoff_
private

◆ rt_extraction_window_

double rt_extraction_window_
private

◆ rt_normalization_factor_

double rt_normalization_factor_
private

◆ scoring_model_

String scoring_model_
private

◆ sn_bin_count_

unsigned int sn_bin_count_
private

◆ sn_win_len_

double sn_win_len_
private

◆ sonarscoring_

OpenMS::SONARScoring sonarscoring_
private

◆ spacing_for_spectra_resampling_

double spacing_for_spectra_resampling_
private

◆ spectrum_addition_method_

String spectrum_addition_method_
private

◆ stop_report_after_feature_

int stop_report_after_feature_
private

◆ strict_

bool strict_
private

◆ su_

◆ uis_threshold_peak_area_

double uis_threshold_peak_area_
private

◆ uis_threshold_sn_

double uis_threshold_sn_
private

◆ use_ms1_ion_mobility_

bool use_ms1_ion_mobility_
private

◆ write_convex_hull_

bool write_convex_hull_
private

◆ write_log_messages_

bool write_log_messages_
private