OpenMS
SpectralDeconvolution.h
Go to the documentation of this file.
1 // Copyright (c) 2002-2023, The OpenMS Team -- EKU Tuebingen, ETH Zurich, and FU Berlin
2 // SPDX-License-Identifier: BSD-3-Clause
3 //
4 // --------------------------------------------------------------------------
5 // $Maintainer: Kyowon Jeong, Jihyung Kim $
6 // $Authors: Kyowon Jeong, Jihyung Kim $
7 // --------------------------------------------------------------------------
8 
9 #pragma once
10 
18 #include <boost/dynamic_bitset.hpp>
19 #include <iostream>
20 
21 namespace OpenMS
22 {
33  class OPENMS_DLLAPI SpectralDeconvolution : public DefaultParamHandler
34  {
35  public:
38 
41 
44 
47 
50 
53 
56 
63  void performSpectrumDeconvolution(const MSSpectrum& spec, int scan_number, const PeakGroup& precursor_peak_group);
64 
67 
70 
73 
78  void setTargetMasses(const std::vector<double>& masses, bool exclude = false);
79 
84  void calculateAveragine(bool use_RNA_averagine, const bool is_centroid = true);
85 
88  {
89  max_mass_dalton_tolerance_ = 100;
90  }
91 
93  static int getNominalMass(double mass);
94 
103  static float getCosine(const std::vector<float>& a, int a_start, int a_end, const IsotopeDistribution& b, int offset, int min_iso_len);
104 
105 
118  static float getIsotopeCosineAndIsoOffset(double mono_mass, const std::vector<float>& per_isotope_intensities, int& offset, const PrecalculatedAveragine& avg, int iso_int_shift = 0,
119  int window_width = -1, int allowed_isotope_error = 0,
120  PeakGroup::TargetDecoyType target_decoy_type = PeakGroup::TargetDecoyType::target);
121 
127  void setTargetDecoyType(PeakGroup::TargetDecoyType target_decoy_type, const DeconvolvedSpectrum& target_dspec_for_decoy_calcualtion);
128 
130  const static int min_iso_size = 2;
131 
132  protected:
133  void updateMembers_() override;
134 
135  private:
137 
139  int allowed_iso_error_ = 1;
140 
142  int min_abs_charge_, max_abs_charge_;
146  double min_mass_, max_mass_;
154  double noise_iso_delta_ = .9444; // should be described in the FDR paper
156  const static int min_support_peak_count_ = 2;
165 
168 
170  PeakGroup::TargetDecoyType target_decoy_type_ = PeakGroup::TargetDecoyType::target;
171 
174 
176  boost::dynamic_bitset<> target_mass_bins_;
177  std::vector<double> target_mono_masses_;
178 
180  std::vector<double> excluded_masses_;
181 
185 
187  std::vector<LogMzPeak> log_mz_peaks_;
191  boost::dynamic_bitset<> binned_log_masses_;
193  boost::dynamic_bitset<> binned_log_mz_peaks_;
194 
196  std::vector<double> universal_pattern_;
199 
201  std::vector<int> binned_universal_pattern_;
204 
208 
210  uint ms_level_;
211 
214 
215  int target_precursor_charge_ = 0;
216  double target_precursor_mz_ = 0;
217 
219  double max_mass_dalton_tolerance_ = .16;
220 
227  static double getBinValue_(Size bin, double min_value, double bin_mul_factor);
228 
235  static Size getBinNumber_(double value, double min_value, double bin_mul_factor);
236 
239 
244  void binLogMzPeaks_(const Size bin_number, std::vector<float>& binned_log_mz_peak_intensities);
245 
247  double getMassFromMassBin_(Size mass_bin, double bin_mul_factor) const;
248 
250  double getMzFromMzBin_(Size mass_bin, double bin_mul_factor) const;
251 
254 
256  void removeOverlappingPeakGroups_(DeconvolvedSpectrum& dspec, double tol, PeakGroup::TargetDecoyType target_decoy_type = PeakGroup::TargetDecoyType::target);
257 
262  Matrix<int> updateMassBins_(const std::vector<float>& mz_intensities);
263 
268  Matrix<int> filterMassBins_(const std::vector<float>& mass_intensities);
269 
274  void updateCandidateMassBins_(std::vector<float>& mass_intensities, const std::vector<float>& mz_intensities);
275 
279  void getCandidatePeakGroups_(const Matrix<int>& per_mass_abs_charge_ranges);
280 
282  void setFilters_();
283 
286 
289 
291  void removeExcludedMasses_(DeconvolvedSpectrum& dspec, std::vector<double> excluded_masses) const;
292 
294 
295  };
296 } // namespace OpenMS
A class representing a deconvolved spectrum. DeconvolvedSpectrum consists of PeakGroups representing ...
Definition: DeconvolvedSpectrum.h:30
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:66
log transformed peak. After deconvolution, all necessary information from deconvolution such as charg...
Definition: FLASHHelperClasses.h:177
Averagine patterns pre-calculated for speed up. Other variables are also calculated for fast cosine c...
Definition: FLASHHelperClasses.h:35
Definition: IsotopeDistribution.h:39
The representation of a 1D spectrum.
Definition: MSSpectrum.h:44
Class describing a deconvolved mass. A mass contains multiple (LogMz) peaks of different charges and ...
Definition: PeakGroup.h:26
TargetDecoyType
target decoy type of PeakGroup. This specifies if a PeakGroup is a target (0), charge decoy (1),...
Definition: PeakGroup.h:33
FLASHDeconv algorithm: ultrafast mass deconvolution algorithm for top down mass spectrometry dataset ...
Definition: SpectralDeconvolution.h:34
const DeconvolvedSpectrum * target_dspec_for_decoy_calculation_
the peak group vector from normal run. This is used when dummy masses are generated.
Definition: SpectralDeconvolution.h:167
std::vector< double > excluded_masses_
mass bins that are excluded for FLASHIda global targeting mode
Definition: SpectralDeconvolution.h:180
DeconvolvedSpectrum & getDeconvolvedSpectrum()
return deconvolved spectrum
const PrecalculatedAveragine & getAveragine()
get calculated averagine. Call after calculateAveragine is called.
FLASHHelperClasses::PrecalculatedAveragine PrecalculatedAveragine
Definition: SpectralDeconvolution.h:36
FLASHHelperClasses::PrecalculatedAveragine avg_
precalculated averagine distributions for fast averagine generation
Definition: SpectralDeconvolution.h:173
double getMassFromMassBin_(Size mass_bin, double bin_mul_factor) const
get mass value for input mass bin
Matrix< int > updateMassBins_(const std::vector< float > &mz_intensities)
Update binned_log_masses_. It select candidate mass bins using the universal pattern,...
void calculateAveragine(bool use_RNA_averagine, const bool is_centroid=true)
precalculate averagine (for predefined mass bins) to speed up averagine generation
void updateLogMzPeaks_()
generate log mz peaks from the input spectrum
DeconvolvedSpectrum deconvolved_spectrum_
selected_peak_groups_ stores the deconvolved mass peak groups
Definition: SpectralDeconvolution.h:189
double current_max_mass_
max mass is controlled by precursor mass for MSn n>1; otherwise just max_mass
Definition: SpectralDeconvolution.h:150
static double getBinValue_(Size bin, double min_value, double bin_mul_factor)
static function that converts bin to value
bool is_positive_
is positive mode
Definition: SpectralDeconvolution.h:144
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 un...
void removeChargeErrorPeakGroups_(DeconvolvedSpectrum &dspec, const PeakGroup::TargetDecoyType &target_decoy_type) const
filter out charge error masses
void removeExcludedMasses_(DeconvolvedSpectrum &dspec, std::vector< double > excluded_masses) const
filter out excluded masses
Matrix< int > binned_harmonic_patterns
This stores the patterns for harmonic reduction in binned dimension.
Definition: SpectralDeconvolution.h:203
static int getNominalMass(double mass)
convert double to nominal mass
DoubleList bin_mul_factors_
bin multiplication factor (log mz * bin_mul_factors_ = bin number) - for fast convolution,...
Definition: SpectralDeconvolution.h:160
static float getIsotopeCosineAndIsoOffset(double mono_mass, const std::vector< float > &per_isotope_intensities, int &offset, const PrecalculatedAveragine &avg, int iso_int_shift=0, int window_width=-1, int allowed_isotope_error=0, PeakGroup::TargetDecoyType target_decoy_type=PeakGroup::TargetDecoyType::target)
Examine intensity distribution over isotope indices. Also determines the most plausible isotope index...
SpectralDeconvolution & operator=(SpectralDeconvolution &&fd)=default
move assignment operator
~SpectralDeconvolution()=default
destructor
boost::dynamic_bitset target_mass_bins_
mass bins that are targeted for FLASHIda global targeting mode
Definition: SpectralDeconvolution.h:176
void removeOverlappingPeakGroups_(DeconvolvedSpectrum &dspec, double tol, PeakGroup::TargetDecoyType target_decoy_type=PeakGroup::TargetDecoyType::target)
filter out overlapping masses
Matrix< int > filterMassBins_(const std::vector< float > &mass_intensities)
Subfunction of updateMassBins_.
uint ms_level_
current ms Level
Definition: SpectralDeconvolution.h:210
double getMzFromMzBin_(Size mass_bin, double bin_mul_factor) const
get mz value for input mz bin
std::vector< double > previously_deconved_peak_masses_for_decoy_
Definition: SpectralDeconvolution.h:184
void setToleranceEstimation()
when estimating tolerance, max_mass_dalton_tolerance_ should be large
Definition: SpectralDeconvolution.h:87
SpectralDeconvolution & operator=(const SpectralDeconvolution &fd)=default
assignment operator
void binLogMzPeaks_(const Size bin_number, std::vector< float > &binned_log_mz_peak_intensities)
generate mz bins and intensity per mz bin from log mz peaks
double mass_bin_min_value_
minimum mass and mz values representing the first bin of massBin and mzBin, respectively: to save mem...
Definition: SpectralDeconvolution.h:206
boost::dynamic_bitset previously_deconved_mass_bins_for_decoy_
mass bins that are previously deconvolved and excluded for dummy mass generation
Definition: SpectralDeconvolution.h:183
int current_max_charge_
current_max_charge_: controlled by precursor charge for MSn n>1; otherwise just max_abs_charge_
Definition: SpectralDeconvolution.h:148
double max_mass_
Definition: SpectralDeconvolution.h:146
void setAveragine(const PrecalculatedAveragine &avg)
set calculated averagine
boost::dynamic_bitset binned_log_masses_
binned_log_masses_ stores the selected bins for this spectrum + overlapped spectrum (previous a few s...
Definition: SpectralDeconvolution.h:191
std::vector< double > universal_pattern_
This stores the "universal pattern".
Definition: SpectralDeconvolution.h:196
std::vector< int > binned_universal_pattern_
This stores the "universal pattern" in binned dimension.
Definition: SpectralDeconvolution.h:201
DoubleList min_isotope_cosine_
Isotope cosine threshold for each MS level.
Definition: SpectralDeconvolution.h:162
FLASHHelperClasses::LogMzPeak LogMzPeak
Definition: SpectralDeconvolution.h:37
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...
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 ori...
DoubleList min_snr_
SNR threshold for each MS level.
Definition: SpectralDeconvolution.h:164
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
void setTargetDecoyType(PeakGroup::TargetDecoyType target_decoy_type, const DeconvolvedSpectrum &target_dspec_for_decoy_calcualtion)
double iso_da_distance_
isotope dalton distance
Definition: SpectralDeconvolution.h:213
static float getCosine(const std::vector< float > &a, int a_start, int a_end, const IsotopeDistribution &b, int offset, int min_iso_len)
boost::dynamic_bitset binned_log_mz_peaks_
binned_log_mz_peaks_ stores the binned log mz peaks
Definition: SpectralDeconvolution.h:193
SpectralDeconvolution(SpectralDeconvolution &&other)=default
move constructor
SpectralDeconvolution(const SpectralDeconvolution &)=default
copy constructor
std::vector< LogMzPeak > log_mz_peaks_
Stores log mz peaks.
Definition: SpectralDeconvolution.h:187
static Size getBinNumber_(double value, double min_value, double bin_mul_factor)
static function that converts value to bin
double current_min_mass_
max mass is max_mass for MS1 and 50 for MS2
Definition: SpectralDeconvolution.h:152
std::vector< double > target_mono_masses_
Definition: SpectralDeconvolution.h:177
void scoreAndFilterPeakGroups_()
function for peak group scoring and filtering
double mz_bin_min_value_
Definition: SpectralDeconvolution.h:207
DoubleList tolerance_
tolerance in ppm for each MS level
Definition: SpectralDeconvolution.h:158
void setFilters_()
Make the universal pattern.
void generatePeakGroupsFromSpectrum_()
Generate peak groups from the input spectrum.
SpectralDeconvolution()
default constructor
void getCandidatePeakGroups_(const Matrix< int > &per_mass_abs_charge_ranges)
For selected masses in binned_log_masses_, select the peaks from the original spectrum....
int max_abs_charge_
Definition: SpectralDeconvolution.h:142
Matrix< double > harmonic_pattern_matrix_
This stores the patterns for harmonic reduction.
Definition: SpectralDeconvolution.h:198
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:97
std::vector< double > DoubleList
Vector of double precision real types.
Definition: ListUtils.h:36
FLASHIda C++ to C# (or vice versa) bridge functions The functions here are called in C# to invoke fun...
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19
static FLASHHelperClasses::PrecalculatedAveragine avg
keeps the precalculated averagine to calculate average masses from monoisotopic masses
Definition: FLASHIdaBridgeFunctions.h:92