OpenMS
Loading...
Searching...
No Matches
IsobaricIsotopeCorrector Class Reference

Performs isotope impurity correction on intensities extracted from isobaric labeling experiments. More...

#include <OpenMS/ANALYSIS/QUANTITATION/IsobaricIsotopeCorrector.h>

Static Public Member Functions

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 Private Member Functions

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 std::vector< double > getIntensities_ (const IsobaricQuantitationMethod *quant_method, const ConsensusFeature &cf, const ConsensusMap &cm)
 Extract channel intensities from a ConsensusFeature.
 
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 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 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.
 
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 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.
 

Detailed Description

Performs isotope impurity correction on intensities extracted from isobaric labeling experiments.

This class implements algorithms for correcting isotope impurities in quantitative proteomics data obtained from isobaric labeling experiments such as iTRAQ or TMT. Isotope impurities arise from the fact that the reagents used for labeling are not 100% pure and contain isotopic variants that can contribute to the signal in neighboring channels.

The correction is performed using a non-negative least squares (NNLS) approach, which solves the linear system Ax = b, where:

  • A is the correction matrix derived from the isotope impurity information provided by the reagent manufacturer
  • b is the vector of observed intensities in each channel
  • x is the vector of corrected intensities

The NNLS approach ensures that the corrected intensities remain non-negative, which is physically meaningful for mass spectrometry data.

See also
IsobaricQuantitationMethod
IsobaricQuantifierStatistics

Member Function Documentation

◆ computeStats_() [1/2]

static void computeStats_ ( const Matrix< double > &  m_x,
const std::vector< double > &  x_naive,
const float  cf_intensity,
const IsobaricQuantitationMethod quant_method,
IsobaricQuantifierStatistics stats 
)
staticprivate

Compute statistics for the correction process using OpenMS matrices.

Calculates various statistics about the correction process by comparing the NNLS solution (guaranteed non-negative) with the naive LU-decomposition solution (which may have negative values).

Parameters
[in]m_xThe NNLS-corrected intensities as OpenMS matrix (non-negative)
[in]x_naiveThe naive LU-decomposition solution (may contain negative values)
[in]cf_intensityThe original intensity of the consensus feature
[in]quant_methodThe isobaric quantitation method
[in,out]statsThe statistics object to update

◆ computeStats_() [2/2]

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 
)
staticprivate

Compute statistics for the correction process.

Calculates various statistics about the correction process by comparing the NNLS solution (guaranteed non-negative) with the naive LU-decomposition solution (which may have negative values).

Parameters
[in]m_xThe NNLS-corrected intensities (non-negative)
[in]x_naiveThe naive LU-decomposition solution (may contain negative values)
[in]cf_intensityThe original intensity of the consensus feature
[in]quant_methodThe isobaric quantitation method
[in,out]statsThe statistics object to update

◆ correctIsotopicImpurities() [1/2]

static IsobaricQuantifierStatistics correctIsotopicImpurities ( const ConsensusMap consensus_map_in,
ConsensusMap consensus_map_out,
const IsobaricQuantitationMethod quant_method 
)
static

Apply isotope correction to the given input map and store the corrected values in the output map.

Parameters
[in]consensus_map_inThe map containing the values that should be corrected.
[in,out]consensus_map_outThe map where the corrected values should be stored.
[in]quant_methodIsobaricQuantitationMethod (e.g., iTRAQ 4 plex)
Exceptions
Exception::FailedAPICallIf the least-squares fit fails.
Exception::InvalidParameterIf the given correction matrix is invalid.

◆ correctIsotopicImpurities() [2/2]

static void correctIsotopicImpurities ( std::vector< double > &  intensities,
const IsobaricQuantitationMethod quant_method 
)
static

Apply isotope correction to a vector of channel intensities.

This method applies the isotope correction directly to a vector of intensities representing the different isobaric channels. The vector is modified in-place to contain the corrected values.

Parameters
[in,out]intensitiesVector of channel intensities to be corrected (modified in-place)
[in]quant_methodIsobaricQuantitationMethod providing the correction matrix (e.g., iTRAQ 4 plex)
Exceptions
Exception::FailedAPICallIf the least-squares fit fails.
Exception::InvalidParameterIf the given correction matrix is invalid.
Note
The size of the intensities vector must match the number of channels in the quantitation method.

◆ fillInputVector_()

static void fillInputVector_ ( std::vector< double > &  b,
Matrix< double > &  m_b,
const ConsensusFeature cf,
const ConsensusMap cm 
)
staticprivate

Fills the input vector for the NNLS step given the ConsensusFeature.

Warning, assumes that the consensusMap and its ConsensusFeatures have exactly the same cardinality as the number of channels as in the quantitation method and are in the same order as the channels.

I.e. for a TMT16plex, although the whole ConsensusMap has 160 potential map_index values because we had 10 files, every ConsensusFeature is only allowed to have exactly 16 map_index values (one for each channel) and they are in the same order as the channels in the quantitation method.

Parameters
[out]bVector to be filled with intensities
[out]m_bOpenMS matrix to be filled with intensities (alternative representation)
[in]cfConsensusFeature containing the channel intensities
[in]cmConsensusMap containing the feature
Precondition
The ConsensusFeature must contain exactly the same number of elements as there are channels in the quantitation method

◆ getIntensities_()

static std::vector< double > getIntensities_ ( const IsobaricQuantitationMethod quant_method,
const ConsensusFeature cf,
const ConsensusMap cm 
)
staticprivate

Extract channel intensities from a ConsensusFeature.

Extracts the intensities for each channel from the given ConsensusFeature.

Parameters
[in]quant_methodThe isobaric quantitation method defining the channels
[in]cfConsensusFeature containing the channel intensities
[in]cmConsensusMap containing the feature
Returns
Vector of intensities, one for each channel

◆ solveNNLS_() [1/2]

static void solveNNLS_ ( const Matrix< double > &  correction_matrix,
const Matrix< double > &  m_b,
Matrix< double > &  m_x 
)
staticprivate

Solve the non-negative least squares problem using OpenMS matrices.

Solves the NNLS problem Ax = b, where A is the correction matrix and b is the vector of observed intensities. This version uses OpenMS Matrix objects for the computation.

Parameters
[in]correction_matrixThe isotope correction matrix (A)
[in]m_bThe vector of observed intensities (b)
[in]m_xThe output vector of corrected intensities (x)
Exceptions
Exception::FailedAPICallIf the NNLS solver fails

◆ solveNNLS_() [2/2]

static void solveNNLS_ ( Matrix< double > &  correction_matrix,
std::vector< double > &  b,
std::vector< double > &  x 
)
staticprivate

Solve the non-negative least squares problem using Matrix and vectors.

Solves the NNLS problem Ax = b, where A is the correction matrix and b is the vector of observed intensities.

Note
This overload mutates the correction matrix in-place for efficiency. Use the const-preserving overload (taking const Matrix<double>&) if the original matrix must remain unchanged.
Parameters
[in,out]correction_matrixThe isotope correction matrix (A). Modified in-place by the NNLS solver.
[in]bThe vector of observed intensities (b)
[out]xThe output vector of corrected intensities (x)

◆ updateOutputMap_() [1/2]

static float updateOutputMap_ ( const ConsensusMap consensus_map_in,
ConsensusMap consensus_map_out,
Size  current_cf,
const Matrix< double > &  m_x 
)
staticprivate

Update the output consensus map with corrected intensities using OpenMS Matrix.

Copies the consensus feature from the input map to the output map and updates its intensities with the corrected values.

Parameters
[in]consensus_map_inThe input consensus map
[out]consensus_map_outThe output consensus map to be updated
[in]current_cfIndex of the current consensus feature
[out]m_xMatrix of corrected intensities
Returns
The sum of corrected intensities

◆ updateOutputMap_() [2/2]

static float updateOutputMap_ ( const ConsensusMap consensus_map_in,
ConsensusMap consensus_map_out,
Size  current_cf,
const std::vector< double > &  m_x 
)
staticprivate

Update the output consensus map with corrected intensities using std::vector.

Copies the consensus feature from the input map to the output map and updates its intensities with the corrected values.

Parameters
[in]consensus_map_inThe input consensus map
[in,out]consensus_map_outThe output consensus map to be updated
[in]current_cfIndex of the current consensus feature
[in]m_xVector of corrected intensities
Returns
The sum of corrected intensities