Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
MultiplexFiltering Class Reference

base class for filtering centroided and profile data for peak patterns More...

#include <OpenMS/TRANSFORMATIONS/FEATUREFINDER/MultiplexFiltering.h>

Inheritance diagram for MultiplexFiltering:
ProgressLogger MultiplexFilteringCentroided MultiplexFilteringProfile

Classes

struct  BlackListEntry
 structure for peak blacklisting More...
 
struct  PeakReference
 structure for peak position in neighbouring spectra More...
 

Public Member Functions

 MultiplexFiltering (const PeakMap &exp_picked, const std::vector< MultiplexIsotopicPeakPattern > patterns, int peaks_per_peptide_min, int peaks_per_peptide_max, bool missing_peaks, double intensity_cutoff, double mz_tolerance, bool mz_tolerance_unit, double peptide_similarity, double averagine_similarity, double averagine_similarity_scaling, String averagine_type="peptide")
 constructor More...
 
- Public Member Functions inherited from ProgressLogger
 ProgressLogger ()
 Constructor. More...
 
 ~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 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 () const
 Ends the progress display. More...
 

Protected Member Functions

int positionsAndBlacklistFilter_ (const MultiplexIsotopicPeakPattern &pattern, int spectrum, const std::vector< double > &peak_position, int peak, std::vector< double > &mz_shifts_actual, std::vector< int > &mz_shifts_actual_indices) const
 position and blacklist filter More...
 
bool monoIsotopicPeakIntensityFilter_ (const MultiplexIsotopicPeakPattern &pattern, int spectrum_index, const std::vector< int > &mz_shifts_actual_indices) const
 mono-isotopic peak intensity filter More...
 
bool zerothPeakFilter_ (const MultiplexIsotopicPeakPattern &pattern, const std::vector< double > &intensities_actual) const
 zeroth peak filter More...
 
bool peptideSimilarityFilter_ (const MultiplexIsotopicPeakPattern &pattern, const std::vector< double > &intensities_actual, int peaks_found_in_all_peptides_spline) const
 peptide similarity filter More...
 
bool averagineSimilarityFilter_ (const MultiplexIsotopicPeakPattern &pattern, const std::vector< double > &intensities_actual, int peaks_found_in_all_peptides_spline, double mz) const
 averagine similarity filter More...
 
void blacklistPeaks_ (const MultiplexIsotopicPeakPattern &pattern, int spectrum, const std::vector< int > &mz_shifts_actual_indices, int peaks_found_in_all_peptides_spline)
 blacklist peaks More...
 
int getPeakIndex_ (const std::vector< double > &peak_position, int start, double mz, double scaling) const
 returns the index of a peak at m/z (finds not only a valid peak, i.e. within certain m/z deviation, but the best of the valid peaks) More...
 
double getPatternSimilarity_ (const std::vector< double > &pattern1, const std::vector< double > &pattern2) const
 returns similarity of two isotope patterns (simple Pearson correlation coefficient) More...
 
double getAveragineSimilarity_ (const std::vector< double > &pattern, double m) const
 returns similarity of an isotope pattern and an averagine pattern at mass m More...
 

Protected Attributes

PeakMap exp_picked_
 centroided experimental data More...
 
std::vector< std::vector< PeakReference > > registry_
 auxiliary structs for navigation and blacklisting More...
 
std::vector< std::vector< BlackListEntry > > blacklist_
 
std::vector< MultiplexIsotopicPeakPatternpatterns_
 list of peak patterns More...
 
int peaks_per_peptide_min_
 minimum number of isotopic peaks per peptide More...
 
int peaks_per_peptide_max_
 maximum number of isotopic peaks per peptide More...
 
bool missing_peaks_
 flag for missing peaks More...
 
double intensity_cutoff_
 intensity cutoff More...
 
double mz_tolerance_
 m/z shift tolerance More...
 
bool mz_tolerance_unit_
 unit for m/z shift tolerance (ppm - true, Da - false) More...
 
double peptide_similarity_
 peptide similarity More...
 
double averagine_similarity_
 averagine similarity More...
 
double averagine_similarity_scaling_
 averagine similarity scaling More...
 
String averagine_type_
 type of averagine to use More...
 
- Protected Attributes inherited from ProgressLogger
LogType type_
 
time_t last_invoke_
 
ProgressLoggerImplcurrent_logger_
 

Additional Inherited Members

- Public Types inherited from ProgressLogger
enum  LogType { CMD, GUI, NONE }
 Possible log types. More...
 
- Static Protected Member Functions inherited from ProgressLogger
static String logTypeToFactoryName_ (LogType type)
 Return the name of the factory product used for this log type. More...
 
- Static Protected Attributes inherited from ProgressLogger
static int recursion_depth_
 

Detailed Description

base class for filtering centroided and profile data for peak patterns

The algorithm searches for patterns of multiple peptides in the data. The peptides appear as characteristic patterns of isotopic peaks in MS1 spectra. We first search the centroided data, and optionally in a second step the spline interpolated profile data. For each peak pattern the algorithm generates a filter result.

The algorithm differs slightly for centroided and profile input data. This base class comprises code common to both. The two child classes MultiplexFilteringCentroided and MultiplexFilteringProfile contain specific functions and the primary filter() method.

See also
MultiplexIsotopicPeakPattern
MultiplexFilterResult
MultiplexFilteringCentroided
MultiplexFilteringProfile

Constructor & Destructor Documentation

◆ MultiplexFiltering()

MultiplexFiltering ( const PeakMap exp_picked,
const std::vector< MultiplexIsotopicPeakPattern patterns,
int  peaks_per_peptide_min,
int  peaks_per_peptide_max,
bool  missing_peaks,
double  intensity_cutoff,
double  mz_tolerance,
bool  mz_tolerance_unit,
double  peptide_similarity,
double  averagine_similarity,
double  averagine_similarity_scaling,
String  averagine_type = "peptide" 
)

constructor

Parameters
exp_pickedexperimental data in centroid mode
patternspatterns of isotopic peaks to be searched for
peaks_per_peptide_minminimum number of isotopic peaks in peptides
peaks_per_peptide_maxmaximum number of isotopic peaks in peptides
missing_peaksflag for missing peaks
intensity_cutoffintensity cutoff
mz_toleranceerror margin in m/z for matching expected patterns to experimental data
mz_tolerance_unitunit for mz_tolerance, ppm (true), Da (false)
peptide_similaritysimilarity score for two peptides in the same multiplet
averagine_similaritysimilarity score for peptide isotope pattern and averagine model
averagine_similarity_scalingscaling factor x for the averagine similarity parameter p when detecting peptide singlets. With p' = p + x(1-p).

Member Function Documentation

◆ averagineSimilarityFilter_()

bool averagineSimilarityFilter_ ( const MultiplexIsotopicPeakPattern pattern,
const std::vector< double > &  intensities_actual,
int  peaks_found_in_all_peptides_spline,
double  mz 
) const
protected

averagine similarity filter

Checks similarity of the isotopic distribution with the expected averagine distribution. Does the isotope distribution look like a peptide?

Parameters
patternpattern of isotopic peaks to be searched for
intensities_actualspline-interpolated intensities at the actual m/z shift positions
peaks_found_in_all_peptides_splinenumber of isotopic peaks seen for each peptide (profile)
mzm/z at which the averagine distribution is calculated
Returns
true if isotope distribution looks like an average peptide

◆ blacklistPeaks_()

void blacklistPeaks_ ( const MultiplexIsotopicPeakPattern pattern,
int  spectrum,
const std::vector< int > &  mz_shifts_actual_indices,
int  peaks_found_in_all_peptides_spline 
)
protected

blacklist peaks

If a datapoint passes all filters, the corresponding peak in this and the two neighbouring spectra is blacklisted.

Parameters
patternpattern of isotopic peaks to be searched for
spectrumindex of the spectrum in exp_picked_ and boundaries_
peaks_found_in_all_peptides_splinenumber of isotopic peaks seen for each peptide (profile)

◆ getAveragineSimilarity_()

double getAveragineSimilarity_ ( const std::vector< double > &  pattern,
double  m 
) const
protected

returns similarity of an isotope pattern and an averagine pattern at mass m

Parameters
patternisotope pattern
mmass at which the averagine distribution is calculated
Returns
similarity (+1 best, -1 worst)

◆ getPatternSimilarity_()

double getPatternSimilarity_ ( const std::vector< double > &  pattern1,
const std::vector< double > &  pattern2 
) const
protected

returns similarity of two isotope patterns (simple Pearson correlation coefficient)

Parameters
pattern1isotope pattern 1
pattern2isotope pattern 2
Returns
similarity (+1 best, -1 worst)

◆ getPeakIndex_()

int getPeakIndex_ ( const std::vector< double > &  peak_position,
int  start,
double  mz,
double  scaling 
) const
protected

returns the index of a peak at m/z (finds not only a valid peak, i.e. within certain m/z deviation, but the best of the valid peaks)

Parameters
peak_positionm/z position of the peaks
startindex in peak_position for starting the search
mzm/z position of the peak
scalingrescaling of limits
Returns
index of the peak in spectrum

◆ monoIsotopicPeakIntensityFilter_()

bool monoIsotopicPeakIntensityFilter_ ( const MultiplexIsotopicPeakPattern pattern,
int  spectrum_index,
const std::vector< int > &  mz_shifts_actual_indices 
) const
protected

mono-isotopic peak intensity filter

Quick check if the intensities of the mono-isotopic peaks are above the intensity cutoff.

Parameters
patternpattern of isotopic peaks to be searched for
spectrum_indexindex of the spectrum in exp_picked_ and boundaries_
mz_shifts_actual_indicesindices of peaks corresponding to the pattern
Returns
true if all intensities above threshold

◆ peptideSimilarityFilter_()

bool peptideSimilarityFilter_ ( const MultiplexIsotopicPeakPattern pattern,
const std::vector< double > &  intensities_actual,
int  peaks_found_in_all_peptides_spline 
) const
protected

peptide similarity filter

The algorithm takes only MS1 spectra into account i.e. we have no knowledge of the peptide sequences. But we do know that peptides in a pair should have the same sequence and hence the same isotopic distributions. The filter checks the similarity of the lightest peptide with all of the other peptides of the pattern. (In high-complexity samples two peptides can have the correct mass shift by chance. Such accidental pairs show different isotopic distributions and are therefore filtered out.)

Parameters
patternpattern of isotopic peaks to be searched for
intensities_actualspline-interpolated intensities at the actual m/z shift positions
peaks_found_in_all_peptides_splinenumber of isotopic peaks seen for each peptide (profile)
Returns
true if peptide isotope patterns are similar

◆ positionsAndBlacklistFilter_()

int positionsAndBlacklistFilter_ ( const MultiplexIsotopicPeakPattern pattern,
int  spectrum,
const std::vector< double > &  peak_position,
int  peak,
std::vector< double > &  mz_shifts_actual,
std::vector< int > &  mz_shifts_actual_indices 
) const
protected

position and blacklist filter

Checks if there are peaks at positions corresponding to the pattern and that these peaks are not blacklisted.

Parameters
patternpattern of isotopic peaks to be searched for
spectrumindex of the spectrum in exp_picked_ and boundaries_
peak_positionm/z positions of the peaks in spectrum
peakindex of the peak in peak_position
mz_shifts_actualoutput for actual m/z shifts seen in the spectrum (will differ slightly from expected m/z shifts in pattern)
mz_shifts_actual_indicesoutput for indices of peaks corresponding to the pattern
Returns
number of isotopic peaks seen for each peptide

◆ zerothPeakFilter_()

bool zerothPeakFilter_ ( const MultiplexIsotopicPeakPattern pattern,
const std::vector< double > &  intensities_actual 
) const
protected

zeroth peak filter

The mono-isotopic peak is the first peak of each peptide. A peak one m/z shift to the left (e.g. 0.5Th for 2+) is called zeroth peak. High-intensity zeroth peaks indicate incorrect pattern matches. A different pattern is likely to be a better fit.

Parameters
patternpattern of isotopic peaks to be searched for
intensities_actualspline-interpolated intensities at the actual m/z shift positions
Returns
true if there are high-intensity zeroth peaks

Member Data Documentation

◆ averagine_similarity_

double averagine_similarity_
protected

averagine similarity

◆ averagine_similarity_scaling_

double averagine_similarity_scaling_
protected

averagine similarity scaling

◆ averagine_type_

String averagine_type_
protected

type of averagine to use

◆ blacklist_

std::vector<std::vector<BlackListEntry> > blacklist_
protected

◆ exp_picked_

PeakMap exp_picked_
protected

centroided experimental data

◆ intensity_cutoff_

double intensity_cutoff_
protected

intensity cutoff

◆ missing_peaks_

bool missing_peaks_
protected

flag for missing peaks

◆ mz_tolerance_

double mz_tolerance_
protected

m/z shift tolerance

◆ mz_tolerance_unit_

bool mz_tolerance_unit_
protected

unit for m/z shift tolerance (ppm - true, Da - false)

◆ patterns_

std::vector<MultiplexIsotopicPeakPattern> patterns_
protected

list of peak patterns

◆ peaks_per_peptide_max_

int peaks_per_peptide_max_
protected

maximum number of isotopic peaks per peptide

◆ peaks_per_peptide_min_

int peaks_per_peptide_min_
protected

minimum number of isotopic peaks per peptide

◆ peptide_similarity_

double peptide_similarity_
protected

peptide similarity

◆ registry_

std::vector<std::vector<PeakReference> > registry_
protected

auxiliary structs for navigation and blacklisting


OpenMS / TOPP release 2.3.0 Documentation generated on Tue Jan 9 2018 18:22:11 using doxygen 1.8.13