OpenMS
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 
14 // forward decl
15 namespace Eigen
16 {
17  template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
18  class Matrix;
19  using MatrixXd = Matrix<double, -1, -1, 0, -1, -1>;
20  using VectorXd = Matrix<double, -1, 1, 0, -1, 1>;
21 }
22 
23 namespace OpenMS
24 {
25  class IsobaricQuantitationMethod;
26  class IsobaricQuantifierStatistics;
27  class ConsensusMap;
28  class ConsensusFeature;
29 
50  class OPENMS_DLLAPI IsobaricIsotopeCorrector
51  {
52 public:
64  ConsensusMap& consensus_map_out,
65  const IsobaricQuantitationMethod* quant_method);
66 
81  static void
82  correctIsotopicImpurities(std::vector<double> & intensities,
83  const IsobaricQuantitationMethod* quant_method);
84 
85  // No instance methods as this is a purely static class
86 
87 private:
107  Matrix<double>& m_b,
108  const ConsensusFeature& cf,
109  const ConsensusMap& cm);
110 
121  static std::vector<double> getIntensities_(const IsobaricQuantitationMethod* quant_method, const ConsensusFeature& cf, const ConsensusMap& cm);
122 
135  static void solveNNLS_(const Matrix<double>& correction_matrix,
136  const Matrix<double>& m_b, Matrix<double>& m_x);
137 
148  static void solveNNLS_(Matrix<double>::EigenMatrixType& correction_matrix, std::vector<double>& b, std::vector<double>& x);
149 
154  static void computeStats_(const std::vector<double>& m_x,
155  const Eigen::MatrixXd& x, const float cf_intensity,
156  const IsobaricQuantitationMethod* quant_method, IsobaricQuantifierStatistics& stats);
157 
158 
171  static void computeStats_(const Matrix<double>& m_x,
172  const Eigen::MatrixXd& x, const float cf_intensity,
173  const IsobaricQuantitationMethod* quant_method, IsobaricQuantifierStatistics& stats);
174 
187  static float updateOutputMap_(const ConsensusMap& consensus_map_in,
188  ConsensusMap& consensus_map_out,
189  Size current_cf,
190  const std::vector<double>& m_x);
191 
204  static float updateOutputMap_(const ConsensusMap& consensus_map_in,
205  ConsensusMap& consensus_map_out,
206  Size current_cf,
207  const Matrix<double>& m_x);
208  };
209 } // namespace
Definition: IsobaricIsotopeCorrector.h:18
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:51
static void solveNNLS_(Matrix< double >::EigenMatrixType &correction_matrix, std::vector< double > &b, std::vector< double > &x)
Solve the non-negative least squares problem using Eigen matrices.
static void computeStats_(const std::vector< double > &m_x, const Eigen::MatrixXd &x, const float cf_intensity, const IsobaricQuantitationMethod *quant_method, IsobaricQuantifierStatistics &stats)
static void fillInputVector_(Eigen::VectorXd &b, Matrix< double > &m_b, const ConsensusFeature &cf, const ConsensusMap &cm)
Fills the input vector for the Eigen/NNLS step given the ConsensusFeature.
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 std::vector< double > getIntensities_(const IsobaricQuantitationMethod *quant_method, const ConsensusFeature &cf, const ConsensusMap &cm)
Extract channel intensities from a ConsensusFeature.
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 computeStats_(const Matrix< double > &m_x, const Eigen::MatrixXd &x, 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 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.
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
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:97
Definition: IsobaricIsotopeCorrector.h:16
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19