OpenMS
Loading...
Searching...
No Matches
IsobaricIsotopeCorrector.h
Go to the documentation of this file.
1// Copyright (c) 2002-present, OpenMS Inc. -- EKU Tuebingen, ETH Zurich, and FU Berlin
2// SPDX-License-Identifier: BSD-3-Clause
3//
4// --------------------------------------------------------------------------
5// $Maintainer: Timo Sachsenberg $
6// $Authors: Stephan Aiche, Chris Bielow $
7// --------------------------------------------------------------------------
8
9#pragma once
10
12#include <memory>
13#include <vector>
14
15namespace OpenMS
16{
17 class IsobaricQuantitationMethod;
18 class IsobaricQuantifierStatistics;
19 class ConsensusMap;
20 class ConsensusFeature;
21
42 class OPENMS_DLLAPI IsobaricIsotopeCorrector
43 {
44public:
56 ConsensusMap& consensus_map_out,
57 const IsobaricQuantitationMethod* quant_method);
58
73 static void
74 correctIsotopicImpurities(std::vector<double>& intensities,
75 const IsobaricQuantitationMethod* quant_method);
76
77 // No instance methods as this is a purely static class
78
79private:
98 static void fillInputVector_(std::vector<double>& b,
99 Matrix<double>& m_b,
100 const ConsensusFeature& cf,
101 const ConsensusMap& cm);
102
113 static std::vector<double> getIntensities_(const IsobaricQuantitationMethod* quant_method,
114 const ConsensusFeature& cf,
115 const ConsensusMap& cm);
116
129 static void solveNNLS_(const Matrix<double>& correction_matrix,
130 const Matrix<double>& m_b, Matrix<double>& m_x);
131
144 static void solveNNLS_(Matrix<double>& correction_matrix,
145 std::vector<double>& b,
146 std::vector<double>& x);
147
160 static void computeStats_(const std::vector<double>& m_x,
161 const std::vector<double>& x_naive,
162 const float cf_intensity,
163 const IsobaricQuantitationMethod* quant_method,
165
178 static void computeStats_(const Matrix<double>& m_x,
179 const std::vector<double>& x_naive,
180 const float cf_intensity,
181 const IsobaricQuantitationMethod* quant_method,
183
196 static float updateOutputMap_(const ConsensusMap& consensus_map_in,
197 ConsensusMap& consensus_map_out,
198 Size current_cf,
199 const std::vector<double>& m_x);
200
213 static float updateOutputMap_(const ConsensusMap& consensus_map_in,
214 ConsensusMap& consensus_map_out,
215 Size current_cf,
216 const Matrix<double>& m_x);
217 };
218} // namespace
A consensus feature spanning multiple LC-MS/MS experiments.
Definition ConsensusFeature.h:45
A container for consensus elements.
Definition ConsensusMap.h:68
Performs isotope impurity correction on intensities extracted from isobaric labeling experiments.
Definition IsobaricIsotopeCorrector.h:43
static void solveNNLS_(Matrix< double > &correction_matrix, std::vector< double > &b, std::vector< double > &x)
Solve the non-negative least squares problem using Matrix and vectors.
static float updateOutputMap_(const ConsensusMap &consensus_map_in, ConsensusMap &consensus_map_out, Size current_cf, const std::vector< double > &m_x)
Update the output consensus map with corrected intensities using std::vector.
static float updateOutputMap_(const ConsensusMap &consensus_map_in, ConsensusMap &consensus_map_out, Size current_cf, const Matrix< double > &m_x)
Update the output consensus map with corrected intensities using OpenMS Matrix.
static void fillInputVector_(std::vector< double > &b, Matrix< double > &m_b, const ConsensusFeature &cf, const ConsensusMap &cm)
Fills the input vector for the NNLS step given the ConsensusFeature.
static void computeStats_(const Matrix< double > &m_x, const std::vector< double > &x_naive, const float cf_intensity, const IsobaricQuantitationMethod *quant_method, IsobaricQuantifierStatistics &stats)
Compute statistics for the correction process using OpenMS matrices.
static IsobaricQuantifierStatistics correctIsotopicImpurities(const ConsensusMap &consensus_map_in, ConsensusMap &consensus_map_out, const IsobaricQuantitationMethod *quant_method)
Apply isotope correction to the given input map and store the corrected values in the output map.
static std::vector< double > getIntensities_(const IsobaricQuantitationMethod *quant_method, const ConsensusFeature &cf, const ConsensusMap &cm)
Extract channel intensities from a ConsensusFeature.
static void correctIsotopicImpurities(std::vector< double > &intensities, const IsobaricQuantitationMethod *quant_method)
Apply isotope correction to a vector of channel intensities.
static void solveNNLS_(const Matrix< double > &correction_matrix, const Matrix< double > &m_b, Matrix< double > &m_x)
Solve the non-negative least squares problem using OpenMS matrices.
static void computeStats_(const std::vector< double > &m_x, const std::vector< double > &x_naive, const float cf_intensity, const IsobaricQuantitationMethod *quant_method, IsobaricQuantifierStatistics &stats)
Compute statistics for the correction process.
Statistics for quantitation performance and comparison of NNLS vs. naive method (aka matrix inversion...
Definition IsobaricQuantifierStatistics.h:23
Abstract base class describing an isobaric quantitation method in terms of the used channels and an i...
Definition IsobaricQuantitationMethod.h:32
A 2D matrix class with efficient buffer access for NumPy interoperability.
Definition Matrix.h:35
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition Types.h:97
Main OpenMS namespace.
Definition openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19