OpenMS
SpectralDeconvolution Class Reference

FLASHDeconv algorithm: ultrafast mass deconvolution algorithm for top down mass spectrometry dataset From MSSpectrum, this class outputs DeconvolvedSpectrum. Deconvolution takes three steps: i) decharging and select candidate masses - speed up via binning ii) collecting isotopes from the candidate masses and deisotope - peak groups are defined here iii) scoring and filter out low scoring masses (i.e., peak groups) More...

#include <OpenMS/ANALYSIS/TOPDOWN/SpectralDeconvolution.h>

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

Public Types

typedef FLASHHelperClasses::PrecalculatedAveragine PrecalculatedAveragine
 
typedef FLASHHelperClasses::LogMzPeak LogMzPeak
 

Public Member Functions

 SpectralDeconvolution ()
 default constructor More...
 
 SpectralDeconvolution (const SpectralDeconvolution &)=default
 copy constructor More...
 
 SpectralDeconvolution (SpectralDeconvolution &&other)=default
 move constructor More...
 
SpectralDeconvolutionoperator= (const SpectralDeconvolution &fd)=default
 assignment operator More...
 
SpectralDeconvolutionoperator= (SpectralDeconvolution &&fd)=default
 move assignment operator More...
 
 ~SpectralDeconvolution ()=default
 destructor More...
 
void performSpectrumDeconvolution (const MSSpectrum &spec, int scan_number, const PeakGroup &precursor_peak_group)
 main deconvolution function that generates the deconvolved target and dummy spectrum based on the original spectrum. More...
 
DeconvolvedSpectrumgetDeconvolvedSpectrum ()
 return deconvolved spectrum More...
 
const PrecalculatedAveraginegetAveragine ()
 get calculated averagine. Call after calculateAveragine is called. More...
 
void setAveragine (const PrecalculatedAveragine &avg)
 set calculated averagine More...
 
void setTargetMasses (const std::vector< double > &masses, bool exclude=false)
 set targeted or excluded masses for targeted deconvolution. Masses are targeted or excluded in all ms levels. More...
 
void calculateAveragine (bool use_RNA_averagine)
 precalculate averagine (for predefined mass bins) to speed up averagine generation More...
 
void setToleranceEstimation ()
 when estimating tolerance, max_mass_dalton_tolerance_ should be large More...
 
void setTargetDecoyType (PeakGroup::TargetDecoyType target_decoy_type, const DeconvolvedSpectrum &target_dspec_for_decoy_calcualtion)
 
- 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...
 

Static Public Member Functions

static int getNominalMass (double mass)
 convert double to nominal mass More...
 
static float getCosine (const std::vector< float > &a, int a_start, int a_end, const IsotopeDistribution &b, int offset, int min_iso_len)
 
static float getIsotopeCosineAndIsoOffset (double mono_mass, const std::vector< float > &per_isotope_intensities, int &offset, const PrecalculatedAveragine &avg, int iso_int_shift, int window_width, const std::vector< double > &excluded_masses)
 Examine intensity distribution over isotope indices. Also determines the most plausible isotope index or, monoisotopic mono_mass. More...
 
- 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...
 

Static Public Attributes

static const int min_iso_size = 2
 minimum isotopologue count in a peak group More...
 

Protected Member Functions

void updateMembers_ () override
 This method is used to update extra member variables at the end of the setParameters() method. More...
 
- Protected Member Functions inherited from DefaultParamHandler
void defaultsToParam_ ()
 Updates the parameters after the defaults have been set in the constructor. More...
 

Private Member Functions

void updateLogMzPeaks_ ()
 generate log mz peaks from the input spectrum More...
 
void binLogMzPeaks_ (Size bin_number, std::vector< float > &binned_log_mz_peak_intensities)
 generate mz bins and intensity per mz bin from log mz peaks More...
 
double getMassFromMassBin_ (Size mass_bin, double bin_mul_factor) const
 get mass value for input mass bin: used for debugging More...
 
double getMzFromMzBin_ (Size mass_bin, double bin_mul_factor) const
 get mz value for input mz bin: used for debugging More...
 
void generatePeakGroupsFromSpectrum_ ()
 Generate peak groups from the input spectrum. More...
 
Matrix< int > updateMassBins_ (const std::vector< float > &mz_intensities)
 Update binned_log_masses_. It select candidate mass bins using the universal pattern, eliminate possible harmonic masses. This function does not perform deisotoping. More...
 
Matrix< int > filterMassBins_ (const std::vector< float > &mass_intensities)
 Subfunction of updateMassBins_. More...
 
void updateCandidateMassBins_ (std::vector< float > &mass_intensities, const std::vector< float > &mz_intensities)
 Subfunction of updateMassBins_. It select candidate masses and update binned_log_masses_ using the universal pattern, eliminate possible harmonic masses. More...
 
void getCandidatePeakGroups_ (const Matrix< int > &per_mass_abs_charge_ranges)
 For selected masses in binned_log_masses_, select the peaks from the original spectrum. Also isotopic peaks are clustered in this function. More...
 
void setFilters_ ()
 Make the universal pattern. More...
 
void scoreAndFilterPeakGroups_ ()
 function for peak group scoring and filtering More...
 
void removeChargeErrorPeakGroups_ (DeconvolvedSpectrum &dspec, const PeakGroup::TargetDecoyType &target_decoy_type) const
 filter out charge error masses More...
 
void removeExcludedMasses_ (DeconvolvedSpectrum &dspec, std::vector< double > excluded_masses, double tol) const
 filter out excluded masses More...
 
void setTargetPrecursorCharge_ ()
 
void prepareSignalDecoyExclusions_ ()
 prepare signal decoy exclusions for decoy runs More...
 
void prepareNoiseDecoySpectrum_ (const MSSpectrum &spec)
 prepare noise decoy spectrum by filtering signal peaks More...
 
void registerPrecursorForMSn_ (const PeakGroup &precursor_peak_group)
 register precursor information for MSn spectra (n > 1) More...
 
bool isPeakGroupInExcludedMassForDecoyRuns_ (const PeakGroup &peak_group, double tol, int offset=0) const
 

Static Private Member Functions

static double getBinValue_ (Size bin, double min_value, double bin_mul_factor)
 static function that converts bin to value More...
 
static Size getBinNumber_ (double value, double min_value, double bin_mul_factor)
 static function that converts value to bin More...
 
static void removeOverlappingPeakGroups_ (DeconvolvedSpectrum &dspec, double tol, PeakGroup::TargetDecoyType target_decoy_type=PeakGroup::TargetDecoyType::target)
 filter out overlapping masses More...
 

Private Attributes

int allowed_iso_error_ = -1
 FLASHDeconv parameters. More...
 
int min_abs_charge_
 min charge and max charge subject to analysis, set by users More...
 
int max_abs_charge_
 
bool is_positive_
 is positive mode More...
 
double min_mass_
 mass ranges of deconvolution, set by users More...
 
double max_mass_
 
int current_max_charge_
 current_max_charge_: controlled by precursor charge for MSn n>1; otherwise just max_abs_charge_ More...
 
double current_max_mass_
 max mass is controlled by precursor mass for MSn n>1; otherwise just max_mass More...
 
double current_min_mass_
 max mass is max_mass for MS1 and 50 for MS2 More...
 
double noise_iso_delta_ = .9444
 isotope distance for noise decoy. This distance has been determined so i) it is close to the mass of 13C - 12C yet is not close to any correct isotope distances of different charges (from 1 to 10). More...
 
DoubleList tolerance_
 tolerance in ppm for each MS level More...
 
DoubleList bin_mul_factors_
 bin multiplication factor (log mz * bin_mul_factors_ = bin number) - for fast convolution, binning is used More...
 
DoubleList min_isotope_cosine_
 Isotope cosine threshold for each MS level. More...
 
DoubleList min_snr_
 snr threshold for each MS level More...
 
double min_qscore_ = .2
 minimum Qscore of a deconvolved mass More...
 
const DeconvolvedSpectrumtarget_dspec_for_decoy_calculation_
 the peak group vector from normal run. This is used when dummy masses are generated. More...
 
PeakGroup::TargetDecoyType target_decoy_type_ = PeakGroup::TargetDecoyType::target
 PeakGroup::TargetDecoyType values. More...
 
FLASHHelperClasses::PrecalculatedAveragine avg_
 precalculated averagine distributions for fast averagine generation More...
 
boost::dynamic_bitset target_mass_bins_
 mass bins that are targeted for FLASHIda global targeting mode More...
 
std::vector< double > target_mono_masses_
 
std::vector< double > excluded_masses_
 mass bins that are excluded for FLASHIda global targeting mode More...
 
boost::dynamic_bitset excluded_mass_bins_for_decoy_runs_
 mass bins that are previously deconvolved and excluded for decoy mass generation More...
 
std::vector< double > excluded_peak_masses_for_decoy_runs_
 
std::vector< double > excluded_masses_for_decoy_runs_
 
std::vector< LogMzPeaklog_mz_peaks_
 Stores log mz peaks. More...
 
DeconvolvedSpectrum deconvolved_spectrum_
 selected_peak_groups_ stores the deconvolved mass peak groups More...
 
boost::dynamic_bitset binned_log_masses_
 binned_log_masses_ stores the selected bins for this spectrum + overlapped spectrum (previous a few spectra). More...
 
boost::dynamic_bitset binned_log_mz_peaks_
 binned_log_mz_peaks_ stores the binned log mz peaks More...
 
std::vector< double > universal_pattern_
 This stores the "universal pattern". More...
 
Matrix< double > harmonic_pattern_matrix_
 This stores the patterns for harmonic reduction. More...
 
std::vector< int > binned_universal_pattern_
 This stores the "universal pattern" in binned dimension. More...
 
Matrix< int > binned_harmonic_patterns
 This stores the patterns for harmonic reduction in binned dimension. More...
 
double mass_bin_min_value_
 minimum mass and mz values representing the first bin of massBin and mzBin, respectively: to save memory space More...
 
double mz_bin_min_value_
 
uint ms_level_
 current ms Level More...
 
double iso_da_distance_
 isotope dalton distance More...
 
int target_precursor_charge_ = 0
 for precursor targetting. More...
 
double target_precursor_mz_ = 0
 
double max_mass_dalton_tolerance_ = .16
 this is additional mass tolerance in Da to get more high signal-to-ratio peaks in this candidate peakgroup finding More...
 

Static Private Attributes

static const int min_support_peak_count_ = 2
 minimum number of peaks supporting a mass More...
 

Additional Inherited Members

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

Detailed Description

FLASHDeconv algorithm: ultrafast mass deconvolution algorithm for top down mass spectrometry dataset From MSSpectrum, this class outputs DeconvolvedSpectrum. Deconvolution takes three steps: i) decharging and select candidate masses - speed up via binning ii) collecting isotopes from the candidate masses and deisotope - peak groups are defined here iii) scoring and filter out low scoring masses (i.e., peak groups)

Member Typedef Documentation

◆ LogMzPeak

◆ PrecalculatedAveragine

Constructor & Destructor Documentation

◆ SpectralDeconvolution() [1/3]

default constructor

◆ SpectralDeconvolution() [2/3]

copy constructor

◆ SpectralDeconvolution() [3/3]

move constructor

◆ ~SpectralDeconvolution()

~SpectralDeconvolution ( )
default

destructor

Member Function Documentation

◆ binLogMzPeaks_()

void binLogMzPeaks_ ( Size  bin_number,
std::vector< float > &  binned_log_mz_peak_intensities 
)
private

generate mz bins and intensity per mz bin from log mz peaks

Parameters
bin_numbernumber of mz bins
binned_log_mz_peak_intensitiesintensity per mz bin

◆ calculateAveragine()

void calculateAveragine ( bool  use_RNA_averagine)

precalculate averagine (for predefined mass bins) to speed up averagine generation

Parameters
use_RNA_averagineif set, averagine for RNA (nucleotides) is calculated

◆ filterMassBins_()

Matrix<int> filterMassBins_ ( const std::vector< float > &  mass_intensities)
private

Subfunction of updateMassBins_.

Parameters
mass_intensitiesper mass bin intensity
Returns
a matrix containing charge ranges for all found masses

◆ generatePeakGroupsFromSpectrum_()

void generatePeakGroupsFromSpectrum_ ( )
private

Generate peak groups from the input spectrum.

◆ getAveragine()

const PrecalculatedAveragine& getAveragine ( )

get calculated averagine. Call after calculateAveragine is called.

◆ getBinNumber_()

static Size getBinNumber_ ( double  value,
double  min_value,
double  bin_mul_factor 
)
staticprivate

static function that converts value to bin

Parameters
valuevalue
min_valueminimum value (corresponding to bin number = 0)
bin_mul_factorbin multiplication factor: bin_number = (bin_value * bin_mul_factors_ - min_value)
Returns
bin corresponding to value

◆ getBinValue_()

static double getBinValue_ ( Size  bin,
double  min_value,
double  bin_mul_factor 
)
staticprivate

static function that converts bin to value

Parameters
binbin number
min_valueminimum value (corresponding to bin number = 0)
bin_mul_factorbin multiplication factor: bin_value = (min_value + bin_number/ bin_mul_factors_)
Returns
value corresponding to bin

◆ getCandidatePeakGroups_()

void getCandidatePeakGroups_ ( const Matrix< int > &  per_mass_abs_charge_ranges)
private

For selected masses in binned_log_masses_, select the peaks from the original spectrum. Also isotopic peaks are clustered in this function.

Parameters
per_mass_abs_charge_rangescharge range per mass

◆ getCosine()

static float getCosine ( const std::vector< float > &  a,
int  a_start,
int  a_end,
const IsotopeDistribution b,
int  offset,
int  min_iso_len 
)
static

calculate cosine between two vectors a and b with additional parameters for fast calculation

Parameters
avector a
a_startnon zero start index of a
a_endnon zero end index of a (exclusive)
bvector b
offsetelement index offset between a and b
min_iso_lenminimum isotope size. If isotope size is less than this, return 0

◆ getDeconvolvedSpectrum()

DeconvolvedSpectrum& getDeconvolvedSpectrum ( )

return deconvolved spectrum

◆ getIsotopeCosineAndIsoOffset()

static float getIsotopeCosineAndIsoOffset ( double  mono_mass,
const std::vector< float > &  per_isotope_intensities,
int &  offset,
const PrecalculatedAveragine avg,
int  iso_int_shift,
int  window_width,
const std::vector< double > &  excluded_masses 
)
static

Examine intensity distribution over isotope indices. Also determines the most plausible isotope index or, monoisotopic mono_mass.

Parameters
mono_massmonoisotopic mass
per_isotope_intensitiesvector of intensities associated with each isotope - aggregated through charges
offsetoutput offset between input monoisotopic mono_mass and determined monoisotopic mono_mass
avgprecalculated averagine
iso_int_shiftisotope shift in per_isotope_intensities.
window_widthisotope offset value range. If -1, set automatically.
excluded_massesmasses not considered in the calculation.
Returns
calculated cosine similarity score

◆ getMassFromMassBin_()

double getMassFromMassBin_ ( Size  mass_bin,
double  bin_mul_factor 
) const
private

get mass value for input mass bin: used for debugging

◆ getMzFromMzBin_()

double getMzFromMzBin_ ( Size  mass_bin,
double  bin_mul_factor 
) const
private

get mz value for input mz bin: used for debugging

◆ getNominalMass()

static int getNominalMass ( double  mass)
static

convert double to nominal mass

◆ isPeakGroupInExcludedMassForDecoyRuns_()

bool isPeakGroupInExcludedMassForDecoyRuns_ ( const PeakGroup peak_group,
double  tol,
int  offset = 0 
) const
private

◆ operator=() [1/2]

SpectralDeconvolution& operator= ( const SpectralDeconvolution fd)
default

assignment operator

◆ operator=() [2/2]

SpectralDeconvolution& operator= ( SpectralDeconvolution &&  fd)
default

move assignment operator

◆ performSpectrumDeconvolution()

void performSpectrumDeconvolution ( const MSSpectrum spec,
int  scan_number,
const PeakGroup precursor_peak_group 
)

main deconvolution function that generates the deconvolved target and dummy spectrum based on the original spectrum.

Parameters
specthe original spectrum
scan_numberscan number from input spectrum.
precursor_peak_groupprecursor peak group

◆ prepareNoiseDecoySpectrum_()

void prepareNoiseDecoySpectrum_ ( const MSSpectrum spec)
private

prepare noise decoy spectrum by filtering signal peaks

◆ prepareSignalDecoyExclusions_()

void prepareSignalDecoyExclusions_ ( )
private

prepare signal decoy exclusions for decoy runs

◆ registerPrecursorForMSn_()

void registerPrecursorForMSn_ ( const PeakGroup precursor_peak_group)
private

register precursor information for MSn spectra (n > 1)

◆ removeChargeErrorPeakGroups_()

void removeChargeErrorPeakGroups_ ( DeconvolvedSpectrum dspec,
const PeakGroup::TargetDecoyType target_decoy_type 
) const
private

filter out charge error masses

◆ removeExcludedMasses_()

void removeExcludedMasses_ ( DeconvolvedSpectrum dspec,
std::vector< double >  excluded_masses,
double  tol 
) const
private

filter out excluded masses

◆ removeOverlappingPeakGroups_()

static void removeOverlappingPeakGroups_ ( DeconvolvedSpectrum dspec,
double  tol,
PeakGroup::TargetDecoyType  target_decoy_type = PeakGroup::TargetDecoyType::target 
)
staticprivate

filter out overlapping masses

◆ scoreAndFilterPeakGroups_()

void scoreAndFilterPeakGroups_ ( )
private

function for peak group scoring and filtering

◆ setAveragine()

void setAveragine ( const PrecalculatedAveragine avg)

set calculated averagine

◆ setFilters_()

void setFilters_ ( )
private

Make the universal pattern.

◆ setTargetDecoyType()

void setTargetDecoyType ( PeakGroup::TargetDecoyType  target_decoy_type,
const DeconvolvedSpectrum target_dspec_for_decoy_calcualtion 
)

set target dummy type for the SpectralDeconvolution run. All masses from the target SpectralDeconvolution run will have the target_decoy_type_.

Parameters
target_decoy_typeThis target_decoy_type_ specifies if a PeakGroup is a target (0), charge dummy (1), noise dummy (2), or isotope dummy (3)
target_dspec_for_decoy_calcualtiontarget masses from normal deconvolution

◆ setTargetMasses()

void setTargetMasses ( const std::vector< double > &  masses,
bool  exclude = false 
)

set targeted or excluded masses for targeted deconvolution. Masses are targeted or excluded in all ms levels.

Parameters
massestarget masses to set
excludeif set, masses are excluded.

◆ setTargetPrecursorCharge_()

void setTargetPrecursorCharge_ ( )
private

◆ setToleranceEstimation()

void setToleranceEstimation ( )
inline

when estimating tolerance, max_mass_dalton_tolerance_ should be large

◆ updateCandidateMassBins_()

void updateCandidateMassBins_ ( std::vector< float > &  mass_intensities,
const std::vector< float > &  mz_intensities 
)
private

Subfunction of updateMassBins_. It select candidate masses and update binned_log_masses_ using the universal pattern, eliminate possible harmonic masses.

Parameters
mass_intensitiesmass bin intensities which are updated in this function
mz_intensitiesmz bin intensities

◆ updateLogMzPeaks_()

void updateLogMzPeaks_ ( )
private

generate log mz peaks from the input spectrum

◆ updateMassBins_()

Matrix<int> updateMassBins_ ( const std::vector< float > &  mz_intensities)
private

Update binned_log_masses_. It select candidate mass bins using the universal pattern, eliminate possible harmonic masses. This function does not perform deisotoping.

Parameters
mz_intensitiesper mz bin intensity
Returns
a matrix containing charge ranges for all found masses

◆ updateMembers_()

void updateMembers_ ( )
overrideprotectedvirtual

This method is used to update extra member variables at the end of the setParameters() method.

Also call it at the end of the derived classes' copy constructor and assignment operator.

The default implementation is empty.

Reimplemented from DefaultParamHandler.

Member Data Documentation

◆ allowed_iso_error_

int allowed_iso_error_ = -1
private

FLASHDeconv parameters.

allowed isotope error in deconvolved mass to calculate qvalue

◆ avg_

precalculated averagine distributions for fast averagine generation

◆ bin_mul_factors_

DoubleList bin_mul_factors_
private

bin multiplication factor (log mz * bin_mul_factors_ = bin number) - for fast convolution, binning is used

◆ binned_harmonic_patterns

Matrix<int> binned_harmonic_patterns
private

This stores the patterns for harmonic reduction in binned dimension.

◆ binned_log_masses_

boost::dynamic_bitset binned_log_masses_
private

binned_log_masses_ stores the selected bins for this spectrum + overlapped spectrum (previous a few spectra).

◆ binned_log_mz_peaks_

boost::dynamic_bitset binned_log_mz_peaks_
private

binned_log_mz_peaks_ stores the binned log mz peaks

◆ binned_universal_pattern_

std::vector<int> binned_universal_pattern_
private

This stores the "universal pattern" in binned dimension.

◆ current_max_charge_

int current_max_charge_
private

current_max_charge_: controlled by precursor charge for MSn n>1; otherwise just max_abs_charge_

◆ current_max_mass_

double current_max_mass_
private

max mass is controlled by precursor mass for MSn n>1; otherwise just max_mass

◆ current_min_mass_

double current_min_mass_
private

max mass is max_mass for MS1 and 50 for MS2

◆ deconvolved_spectrum_

DeconvolvedSpectrum deconvolved_spectrum_
private

selected_peak_groups_ stores the deconvolved mass peak groups

◆ excluded_mass_bins_for_decoy_runs_

boost::dynamic_bitset excluded_mass_bins_for_decoy_runs_
private

mass bins that are previously deconvolved and excluded for decoy mass generation

◆ excluded_masses_

std::vector<double> excluded_masses_
private

mass bins that are excluded for FLASHIda global targeting mode

◆ excluded_masses_for_decoy_runs_

std::vector<double> excluded_masses_for_decoy_runs_
private

◆ excluded_peak_masses_for_decoy_runs_

std::vector<double> excluded_peak_masses_for_decoy_runs_
private

◆ harmonic_pattern_matrix_

Matrix<double> harmonic_pattern_matrix_
private

This stores the patterns for harmonic reduction.

◆ is_positive_

bool is_positive_
private

is positive mode

◆ iso_da_distance_

double iso_da_distance_
private

isotope dalton distance

◆ log_mz_peaks_

std::vector<LogMzPeak> log_mz_peaks_
private

Stores log mz peaks.

◆ mass_bin_min_value_

double mass_bin_min_value_
private

minimum mass and mz values representing the first bin of massBin and mzBin, respectively: to save memory space

◆ max_abs_charge_

int max_abs_charge_
private

◆ max_mass_

double max_mass_
private

◆ max_mass_dalton_tolerance_

double max_mass_dalton_tolerance_ = .16
private

this is additional mass tolerance in Da to get more high signal-to-ratio peaks in this candidate peakgroup finding

◆ min_abs_charge_

int min_abs_charge_
private

min charge and max charge subject to analysis, set by users

◆ min_iso_size

const int min_iso_size = 2
static

minimum isotopologue count in a peak group

◆ min_isotope_cosine_

DoubleList min_isotope_cosine_
private

Isotope cosine threshold for each MS level.

◆ min_mass_

double min_mass_
private

mass ranges of deconvolution, set by users

◆ min_qscore_

double min_qscore_ = .2
private

minimum Qscore of a deconvolved mass

◆ min_snr_

DoubleList min_snr_
private

snr threshold for each MS level

◆ min_support_peak_count_

const int min_support_peak_count_ = 2
staticprivate

minimum number of peaks supporting a mass

◆ ms_level_

uint ms_level_
private

current ms Level

◆ mz_bin_min_value_

double mz_bin_min_value_
private

◆ noise_iso_delta_

double noise_iso_delta_ = .9444
private

isotope distance for noise decoy. This distance has been determined so i) it is close to the mass of 13C - 12C yet is not close to any correct isotope distances of different charges (from 1 to 10).

◆ target_decoy_type_

PeakGroup::TargetDecoyType target_decoy_type_ = PeakGroup::TargetDecoyType::target
private

◆ target_dspec_for_decoy_calculation_

const DeconvolvedSpectrum* target_dspec_for_decoy_calculation_
private

the peak group vector from normal run. This is used when dummy masses are generated.

◆ target_mass_bins_

boost::dynamic_bitset target_mass_bins_
private

mass bins that are targeted for FLASHIda global targeting mode

◆ target_mono_masses_

std::vector<double> target_mono_masses_
private

◆ target_precursor_charge_

int target_precursor_charge_ = 0
private

for precursor targetting.

◆ target_precursor_mz_

double target_precursor_mz_ = 0
private

◆ tolerance_

DoubleList tolerance_
private

tolerance in ppm for each MS level

◆ universal_pattern_

std::vector<double> universal_pattern_
private

This stores the "universal pattern".