OpenMS
FLASHDeconvAlgorithm.h
Go to the documentation of this file.
1 // Copyright (c) 2002-present, 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 FLASHDeconvAlgorithm : public DefaultParamHandler
34  {
35  public:
38 
41 
44 
47 
50 
53 
55  ~FLASHDeconvAlgorithm() = default;
56 
64  void performSpectrumDeconvolution(const MSSpectrum& spec, const std::vector<DeconvolvedSpectrum>& survey_scans, int scan_number,
65  const std::map<int, std::vector<std::vector<float>>>& precursor_map_for_FLASHIda);
66 
69 
72 
75 
81  void setTargetMasses(const std::vector<double>& masses, bool exclude = false);
82 
86  void calculateAveragine(bool use_RNA_averagine);
87 
89  static int getNominalMass(double mass);
90 
100  static float getCosine(const std::vector<float>& a, int a_start, int a_end, const IsotopeDistribution& b, int b_size, int offset, int min_iso_size);
101 
102 
115  static float getIsotopeCosineAndDetermineIsotopeIndex(double mono_mass, const std::vector<float>& per_isotope_intensities, int& offset, const PrecalculatedAveragine& avg, int iso_int_shift = 0,
116  int window_width = -1, int allowed_iso_error_for_second_best_cos = 0,
117  PeakGroup::TargetDummyType target_dummy_type = PeakGroup::TargetDummyType::target);
118 
124  static void addMZsToExcludsionList(const DeconvolvedSpectrum& dspec, std::unordered_set<double>& excluded_mzs);
125 
131  void setTargetDummyType(PeakGroup::TargetDummyType target_dummy_type, DeconvolvedSpectrum& target_dspec_for_dummy_calcualtion);
132 
133  protected:
134  void updateMembers_() override;
135 
136  private:
138 
140  const static int min_iso_size_ = 2;
141 
143  int allowed_iso_error_ = 1;
144 
146  double min_rt_, max_rt_;
148  double min_mz_, max_mz_;
150  int min_abs_charge_, max_abs_charge_;
154  double min_mass_, max_mass_;
164  const static int min_support_peak_count_ = 2;
173 
175  PeakGroup::TargetDummyType target_dummy_type_ = PeakGroup::TargetDummyType::target;
176 
179 
181  boost::dynamic_bitset<> target_mass_bins_;
182  std::vector<double> target_mono_masses_;
183 
185  std::vector<double> excluded_masses_;
186 
190  std::unordered_set<double> excluded_peak_mzs_;
191 
193  std::vector<LogMzPeak> log_mz_peaks_;
197  boost::dynamic_bitset<> mass_bins_;
199  boost::dynamic_bitset<> mz_bins_;
200 
202  std::vector<double> filter_;
205 
208 
210  std::vector<int> bin_offsets_;
213 
217 
219  uint ms_level_;
220 
223 
230  static double getBinValue_(Size bin, double min_value, double bin_mul_factor);
231 
238  static Size getBinNumber_(double value, double min_value, double bin_mul_factor);
239 
242 
247  void updateMzBins_(Size bin_number, std::vector<float>& mz_bin_intensities);
248 
249 
251  double getMassFromMassBin_(Size mass_bin, double bin_mul_factor) const;
252 
254  double getMzFromMzBin_(Size mass_bin, double bin_mul_factor) const;
255 
258 
263  Matrix<int> updateMassBins_(const std::vector<float>& mz_intensities);
264 
269  Matrix<int> filterMassBins_(const std::vector<float>& mass_intensities);
270 
275  void updateCandidateMassBins_(std::vector<float>& mass_intensities, const std::vector<float>& mz_intensities);
276 
280  void getCandidatePeakGroups_(const Matrix<int>& per_mass_abs_charge_ranges);
281 
283  void setFilters_();
284 
287 
290 
293 
302  bool registerPrecursor_(const std::vector<DeconvolvedSpectrum>& survey_scans, const std::map<int, std::vector<std::vector<float>>>& precursor_map_for_real_time_acquisition);
303  };
304 } // 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
FLASHDeconv algorithm: ultrafast mass deconvolution algorithm for top down mass spectrometry dataset ...
Definition: FLASHDeconvAlgorithm.h:34
std::vector< double > excluded_masses_
mass bins that are excluded for FLASHIda global targeting mode
Definition: FLASHDeconvAlgorithm.h:185
DeconvolvedSpectrum & getDeconvolvedSpectrum()
return deconvolved spectrum
FLASHDeconvAlgorithm()
default constructor
const PrecalculatedAveragine & getAveragine()
get calculated averagine. Call after calculateAveragine is called.
double getMassFromMassBin_(Size mass_bin, double bin_mul_factor) const
get mass value for input mass bin
std::vector< double > filter_
This stores the "universal pattern".
Definition: FLASHDeconvAlgorithm.h:202
Matrix< int > updateMassBins_(const std::vector< float > &mz_intensities)
Update mass_bins_. It select candidate mass bins using the universal pattern, eliminate possible harm...
void updateLogMzPeaks_()
generate log mz peaks from the input spectrum
DeconvolvedSpectrum deconvolved_spectrum_
selected_peak_groups_ stores the deconvolved mass peak groups
Definition: FLASHDeconvAlgorithm.h:195
void removeOverlappingPeakGroups_(DeconvolvedSpectrum &dspec, double tol)
filter out overlapping masses
~FLASHDeconvAlgorithm()=default
destructor
double max_mz_
Definition: FLASHDeconvAlgorithm.h:148
double current_max_mass_
max mass is controlled by precursor mass for MSn n>1; otherwise just max_mass
Definition: FLASHDeconvAlgorithm.h:158
static float getIsotopeCosineAndDetermineIsotopeIndex(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_iso_error_for_second_best_cos=0, PeakGroup::TargetDummyType target_dummy_type=PeakGroup::TargetDummyType::target)
Examine intensity distribution over isotope indices. Also determines the most plausible isotope index...
static double getBinValue_(Size bin, double min_value, double bin_mul_factor)
static function that converts bin to value
FLASHDeconvHelperStructs::LogMzPeak LogMzPeak
Definition: FLASHDeconvAlgorithm.h:37
bool is_positive_
is positive mode
Definition: FLASHDeconvAlgorithm.h:152
void updateCandidateMassBins_(std::vector< float > &mass_intensities, const std::vector< float > &mz_intensities)
Subfunction of updateMassBins_. It select candidate masses and update mass_bins_ using the universal ...
void updateMzBins_(Size bin_number, std::vector< float > &mz_bin_intensities)
generate mz bins and intensity per mz bin from log mz peaks
static int getNominalMass(double mass)
convert double to nominal mass
FLASHDeconvAlgorithm(FLASHDeconvAlgorithm &&other)=default
move constructor
DoubleList bin_mul_factors_
bin multiplication factor (log mz * bin_mul_factors_ = bin number) - for fast convolution,...
Definition: FLASHDeconvAlgorithm.h:168
boost::dynamic_bitset target_mass_bins_
mass bins that are targeted for FLASHIda global targeting mode
Definition: FLASHDeconvAlgorithm.h:181
boost::dynamic_bitset previously_deconved_mass_bins_for_dummy_
mass bins that are previsouly deconvolved and excluded for dummy mass generation
Definition: FLASHDeconvAlgorithm.h:188
FLASHDeconvHelperStructs::PrecalculatedAveragine avg_
precalculated averagine distributions for fast averagine generation
Definition: FLASHDeconvAlgorithm.h:178
Matrix< int > filterMassBins_(const std::vector< float > &mass_intensities)
Subfunction of updateMassBins_.
static void addMZsToExcludsionList(const DeconvolvedSpectrum &dspec, std::unordered_set< double > &excluded_mzs)
uint ms_level_
current ms Level
Definition: FLASHDeconvAlgorithm.h:219
double getMzFromMzBin_(Size mass_bin, double bin_mul_factor) const
get mz value for input mz bin
void setTargetDummyType(PeakGroup::TargetDummyType target_dummy_type, DeconvolvedSpectrum &target_dspec_for_dummy_calcualtion)
double mass_bin_min_value_
minimum mass and mz values representing the first bin of massBin and mzBin, respectively: to save mem...
Definition: FLASHDeconvAlgorithm.h:215
int current_max_charge_
current_max_charge_: controlled by precursor charge for MSn n>1; otherwise just max_abs_charge_
Definition: FLASHDeconvAlgorithm.h:156
double max_mass_
Definition: FLASHDeconvAlgorithm.h:154
void setAveragine(const PrecalculatedAveragine &avg)
set calculated averagine
FLASHDeconvAlgorithm(const FLASHDeconvAlgorithm &)=default
copy constructor
double max_rt_
Definition: FLASHDeconvAlgorithm.h:146
FLASHDeconvAlgorithm & operator=(const FLASHDeconvAlgorithm &fd)=default
assignment operator
double isolation_window_size_
default precursor isolation window size.
Definition: FLASHDeconvAlgorithm.h:222
DoubleList min_isotope_cosine_
cosine threshold between observed and theoretical isotope patterns for each MS level
Definition: FLASHDeconvAlgorithm.h:170
Matrix< int > harmonic_bin_offset_matrix_
This stores the patterns for harmonic reduction in binned dimension.
Definition: FLASHDeconvAlgorithm.h:212
std::unordered_set< double > excluded_peak_mzs_
Definition: FLASHDeconvAlgorithm.h:190
void calculateAveragine(bool use_RNA_averagine)
precalculate averagine (for predefined mass bins) to speed up averagine generation
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...
static float getCosine(const std::vector< float > &a, int a_start, int a_end, const IsotopeDistribution &b, int b_size, int offset, int min_iso_size)
double intensity_threshold_
peak intensity threshold subject to analysis
Definition: FLASHDeconvAlgorithm.h:162
boost::dynamic_bitset mass_bins_
mass_bins_ stores the selected bins for this spectrum + overlapped spectrum (previous a few spectra).
Definition: FLASHDeconvAlgorithm.h:197
std::vector< int > bin_offsets_
This stores the "universal pattern" in binned dimension.
Definition: FLASHDeconvAlgorithm.h:210
bool registerPrecursor_(const std::vector< DeconvolvedSpectrum > &survey_scans, const std::map< int, std::vector< std::vector< float >>> &precursor_map_for_real_time_acquisition)
register the precursor peak as well as the precursor peak group (or mass) if possible for MSn (n>1) s...
std::vector< double > previously_deconved_mono_masses_for_dummy_
Definition: FLASHDeconvAlgorithm.h:189
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
void performSpectrumDeconvolution(const MSSpectrum &spec, const std::vector< DeconvolvedSpectrum > &survey_scans, int scan_number, const std::map< int, std::vector< std::vector< float >>> &precursor_map_for_FLASHIda)
main deconvolution function that generates the deconvolved target and dummy spectrum based on the ori...
Matrix< double > harmonic_filter_matrix_
This stores the patterns for harmonic reduction.
Definition: FLASHDeconvAlgorithm.h:204
double iso_da_distance_
isotope dalton distance
Definition: FLASHDeconvAlgorithm.h:207
boost::dynamic_bitset mz_bins_
mz_bins_ stores the binned log mz peaks
Definition: FLASHDeconvAlgorithm.h:199
std::vector< LogMzPeak > log_mz_peaks_
Stores log mz peaks.
Definition: FLASHDeconvAlgorithm.h:193
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: FLASHDeconvAlgorithm.h:160
std::vector< double > target_mono_masses_
Definition: FLASHDeconvAlgorithm.h:182
void scoreAndFilterPeakGroups_()
function for peak group scoring and filtering
double mz_bin_min_value_
Definition: FLASHDeconvAlgorithm.h:216
FLASHDeconvHelperStructs::PrecalculatedAveragine PrecalculatedAveragine
Definition: FLASHDeconvAlgorithm.h:36
void removeChargeErrorPeakGroups_(DeconvolvedSpectrum &dspec) const
filter out charge error masses
DoubleList tolerance_
tolerance in ppm for each MS level
Definition: FLASHDeconvAlgorithm.h:166
void setFilters_()
Make the universal pattern.
FLASHDeconvAlgorithm & operator=(FLASHDeconvAlgorithm &&fd)=default
move assignment operator
void generatePeakGroupsFromSpectrum_()
Generate peak groups from the input spectrum.
void getCandidatePeakGroups_(const Matrix< int > &per_mass_abs_charge_ranges)
For selected masses in mass_bins_, select the peaks from the original spectrum. Also isotopic peaks a...
int max_abs_charge_
Definition: FLASHDeconvAlgorithm.h:150
DeconvolvedSpectrum * target_dspec_for_dummy_calcualtion_
the deconvolved spectrum from normal run. This is used when dummy masses are generated.
Definition: FLASHDeconvAlgorithm.h:172
log transformed peak. After deconvolution, all necessary information from deconvolution such as charg...
Definition: FLASHDeconvHelperStructs.h:139
Averagine patterns pre-calculated for speed up. Other variables are also calculated for fast cosine c...
Definition: FLASHDeconvHelperStructs.h:34
Definition: IsotopeDistribution.h:39
The representation of a 1D spectrum.
Definition: MSSpectrum.h:44
TargetDummyType
target dummy type of PeakGroup. This specifies if a PeakGroup is a target (0), charge dummy (1),...
Definition: PeakGroup.h:33
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
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19