OpenMS
IsobaricWorkflow

Extracts and normalizes isobaric labeling information from an LC-MS/MS experiment.

pot. predecessor tools → IsobaricAnalyzer → pot. successor tools
PeakPickerHiRes IDMapper
FileFilter

The input MSn spectra have to be in centroid mode for the tool to work properly. Use e.g. PeakPickerHiRes to perform centroiding of profile data, if necessary.

This tool currently supports iTRAQ 4-plex and 8-plex, and TMT 6-plex, 10-plex, 11-plex, 16-plex, and 18-plex and higher labeling methods. It extracts the isobaric reporter ion intensities from centroided MS2 or MS3 data (MSn), then performs isotope correction and stores the resulting quantitation in a consensus map, in which each consensus feature represents one relevant MSn scan (e.g. HCD; see parameters select_activation and min_precursor_intensity). The MS level for quantification is chosen automatically, i.e. if MS3 is present, MS2 will be ignored. For intensity, the closest non-zero m/z signal to the theoretical position is taken as reporter ion abundance. The position (RT, m/z) of the consensus centroid is the precursor position in MS1 (from the MS2 spectrum); the consensus sub-elements correspond to the theoretical channel m/z (with m/z values of 113-121 Th for iTRAQ and 126-131 Th for TMT, respectively).

For all labeling techniques, the search radius (reporter_mass_shift) should be set as small as possible, to avoid picking up false-positive ions as reporters. Usually, Orbitraps deliver precision of about 0.0001 Th at this low mass range. Low intensity reporters might have a slightly higher deviation. By default, the mass range is set to ~0.002 Th, which should be sufficient for all instruments (~15 ppm). The tool will throw an Exception if you set it below 0.0001 Th (~0.7ppm). The tool will also throw an Exception if you set reporter_mass_shift > 0.003 Th for TMT-10plex and TMT-11plex, since this could lead to ambiguities with neighbouring channels (which are ~0.006 Th apart in most cases).

For quality control purposes, the tool reports the median distance between the theoretical vs. observed reporter ion peaks in each channel. The search radius is fixed to 0.5 Th (regardless of the user defined search radius). This allows to track calibration issues. For TMT-10plex, these results are automatically omitted if they could be confused with a neighbouring channel, i.e. exceed the tolerance to a neighbouring channel with the same nominal mass (C/N channels). If the distance is too large, you might have a m/z calibration problem (see InternalCalibration).

Note
If none of the reporter ions can be detected in an MSn scan, a consensus feature will still be generated, but the intensities of the overall feature and of all its sub-elements will be zero. (If desired, such features can be removed by applying an intensity filter in FileFilter.) However, if the spectrum is completely empty (no ions whatsoever), no consensus feature will be generated.

Isotope correction is done using non-negative least squares (NNLS), i.e.:
Minimize ||Ax - b||, subject to x >= 0, where b is the vector of observed reporter intensities (with "contaminating" isotope species), A is a correction matrix (as supplied by the manufacturer of the labeling kit) and x is the desired vector of corrected (real) reporter intensities. Other software tools solve this problem by using an inverse matrix multiplication, but this can yield entries in x which are negative. In a real sample, this solution cannot possibly be true, so usually negative values (= negative reporter intensities) are set to zero. However, a negative result usually means that noise was not properly accounted for in the calculation. We thus use NNLS to get a non-negative solution, without the need to truncate negative values. In the (usual) case that inverse matrix multiplication yields only positive values, our NNLS will give the exact same optimal solution.

The correction matrices can be found (and changed) in the INI file (parameter correction_matrix of the corresponding labeling method). However, these matrices for both 4-plex and 8-plex iTRAQ are now stable, and every kit delivered should have the same isotope correction values. Thus, there should be no need to change them, but feel free to compare the values in the INI file with your kit's certificate. For TMT (6-plex and 10-plex) the values have to be adapted for each kit: Modify the correction matrix according to the data in the product data sheet of your charge:

  Data sheet:
  Mass Tag  Repoter Ion -2      -1      Monoisotopic    +1     +2
  126       126.12776   0.0%    0.0%        100%        5.0%   0.0%
  127N      127.124761  0.0%    0.2%        100%        4.6%   0.0%
  ...
  

Corresponding correction matrix:

  [0.0/0.0/5.0/0.0,
  0.0/0.2/4.6/0.0,
  ...
  

The command line parameters of this tool are:

IsobaricWorkflow -- Calculates isobaric quantitative values for peptides
Full documentation: http://www.openms.de/doxygen/nightly/html/TOPP_IsobaricWorkflow.html
Version: 3.5.0-pre-nightly-2025-11-12 Nov 13 2025, 02:53:18, Revision: 9b55e34
To cite OpenMS:
 + Pfeuffer, J., Bielow, C., Wein, S. et al.. OpenMS 3 enables reproducible analysis of large-scale mass spec
   trometry data. Nat Methods (2024). doi:10.1038/s41592-024-02197-7.

Usage:
  IsobaricWorkflow <options>

This tool has algorithm parameters that are not shown here! Please check the ini file for a detailed descript
ion or use the --helphelp option

Options (mandatory options marked with '*'):
  -type <mode>                                                            Isobaric Quantitation method used 
                                                                          in the experiment. (default: 'itraq
                                                                          4plex') (valid: 'itraq4plex', 'itra
                                                                          q8plex', 'tmt10plex', 'tmt11plex', 
                                                                          'tmt16plex', 'tmt18plex', 'tmt6plex
                                                                          ')
  -in <file>*                                                             Input centroided spectrum files 
                                                                          (valid formats: 'mzML')
  -in_id <file>*                                                          Corresponding input PSMs (valid 
                                                                          formats: 'idXML')
  -exp_design <file>                                                      Experimental design file (optional)
                                                                          . If not given, the design is assum
                                                                          ed to be unfractionated. (valid 
                                                                          formats: 'tsv')
  -out <file>*                                                            Output consensusXML file (valid 
                                                                          formats: 'consensusXML')
  -out_mzTab <file>*                                                      Output mzTab file with quantitative
                                                                           information (valid formats: 'mzTab
                                                                          ')
  -calculate_id_purity                                                    Calculate the purity of the precurs
                                                                          or ion based on the MS1 spectrum. 
                                                                          Only used for MS3, otherwise it is 
                                                                          the same as the quant. precursor 
                                                                          purity.
  -psm_score <score>                                                      The score which should be reached 
                                                                          by a peptide hit to be kept.  (use 
                                                                          'NAN' to disable this filter) (defa
                                                                          ult: 'NaN')
  -protein_score <score>                                                  The score which should be reached 
                                                                          by a protein hit to be kept. All 
                                                                          proteins are filtered based on thei
                                                                          r singleton scores irrespective of 
                                                                          grouping. Use in combination with 
                                                                          'delete_unreferenced_peptide_hits' 
                                                                          to remove affected peptides. (use 
                                                                          'NAN' to disable this filter) (defa
                                                                          ult: 'NaN')
  -delete_unreferenced_peptide_hits                                       Peptides not referenced by any prot
                                                                          ein are deleted in the IDs.
  -inference_method <option>                                              Methods used for protein inference 
                                                                          (default: 'aggregation') (valid: 
                                                                          'aggregation', 'bayesian')
  -proteinFDR <threshold>                                                 Protein FDR threshold (0.05=5%). 
                                                                          (default: '1.0') (min: '0.0' max: 
                                                                          '1.0')
  -psmFDR <threshold>                                                     FDR threshold for sub-protein level
                                                                           (e.g. 0.05=5%). (default: '1.0') 
                                                                          (min: '0.0' max: '1.0')

ProteinQuantification:
  -ProteinQuantification:method <choice>                                  - top - quantify based on three 
                                                                          most abundant peptides (number can 
                                                                          be changed in 'top').
                                                                          - iBAQ (intensity based absolute 
                                                                          quantification), calculate the sum 
                                                                          of all peptide peak intensities 
                                                                          divided by the number of theoretica
                                                                          lly observable tryptic peptides 
                                                                          ...
                                                                          input is allowed! (default: 'top') 
                                                                          (valid: 'top', 'iBAQ')
  -ProteinQuantification:best_charge_and_fraction                         Distinguish between fraction and 
                                                                          charge states of a peptide. For 
                                                                          peptides, abundances will be report
                                                                          ed separately for each fraction 
                                                                          and charge;
                                                                          for proteins, abundances will be 
                                                                          computed based only on the most 
                                                                          prevalent charge observed of each 
                                                                          ...
                                                                          over all charge states.

Additional options for custom quantification using top N peptides.:
  -ProteinQuantification:top:N <number>                                   Calculate protein abundance from 
                                                                          this number of proteotypic peptides
                                                                           (most abundant first; '0' for all)
                                                                           (default: '3') (min: '0')
  -ProteinQuantification:top:aggregate <choice>                           Aggregation method used to compute 
                                                                          protein abundances from peptide 
                                                                          abundances (default: 'median') (val
                                                                          id: 'median', 'mean', 'weighted_mea
                                                                          n', 'sum')

Additional options for consensus maps (and identification results comprising multiple runs):
  -ProteinQuantification:consensus:normalize                              Scale peptide abundances so that 
                                                                          medians of all samples are equal
  -ProteinQuantification:consensus:fix_peptides                           Use the same peptides for protein 
                                                                          quantification across all samples.
                                                                          With 'N 0',all peptides that occur 
                                                                          in every sample are considered.
                                                                          Otherwise ('N'), the N peptides 
                                                                          that occur in the most samples (ind
                                                                          ependently of each other) are selec
                                                                          ted,
                                                                          ...
                                                                          n!).

BasicProteinInference:
  -BasicProteinInference:min_peptides_per_protein <number>                Minimal number of peptides needed 
                                                                          for a protein identification. If 
                                                                          set to zero, unmatched proteins 
                                                                          get a score of -Infinity. If bigger
                                                                           than zero, proteins with less pept
                                                                          ides are filtered and evidences 
                                                                          removed from the PSMs. PSMs that 
                                                                          do not reference any proteins anymo
                                                                          re are removed but the spectrum 
                                                                          info is kept. (default: '1') (min: 
                                                                          '0')
  -BasicProteinInference:score_aggregation_method <choice>                How to aggregate scores of peptides
                                                                           matching to the same protein? (def
                                                                          ault: 'best') (valid: 'best', 'prod
                                                                          uct', 'sum', 'maximum')
  -BasicProteinInference:treat_charge_variants_separately <choice>        If this is true, different charge 
                                                                          variants of the same peptide sequen
                                                                          ce count as individual evidences. 
                                                                          (default: 'true') (valid: 'true', 
                                                                          'false')
  -BasicProteinInference:treat_modification_variants_separately <choice>  If this is true, different modifica
                                                                          tion variants of the same peptide 
                                                                          sequence count as individual eviden
                                                                          ces. (default: 'true') (valid: 'tru
                                                                          e', 'false')
  -BasicProteinInference:use_shared_peptides <choice>                     If this is true, shared peptides 
                                                                          are used as evidences. Note: shared
                                                                          _peptides are not deleted and poten
                                                                          tially resolved in postprocessing 
                                                                          as well. (default: 'true') (valid: 
                                                                          'true', 'false')
  -BasicProteinInference:skip_count_annotation                            If this is set, peptide counts won'
                                                                          t be annotated at the proteins.
  -BasicProteinInference:score_type <choice>                              PSM score type to use for inference
                                                                          . (default: empty = main score) 
                                                                          (valid: '', 'PEP', 'q-value', 'RAW'
                                                                          )

BayesianProteinInference:
  -BayesianProteinInference:psm_probability_cutoff <value>                Remove PSMs with probabilities less
                                                                           than this cutoff (default: '1.0e-0
                                                                          3') (min: '0.0' max: '1.0')
  -BayesianProteinInference:top_PSMs <number>                             Consider only top X PSMs per spectr
                                                                          um. 0 considers all. (default: '1')
                                                                           (min: '0')

                                                                          
Common TOPP options:
  -ini <file>                                                             Use the given TOPP INI file
  -threads <n>                                                            Sets the number of threads allowed 
                                                                          to be used by the TOPP tool (defaul
                                                                          t: '1')
  -write_ini <file>                                                       Writes the default configuration 
                                                                          file
  --help                                                                  Shows options
  --helphelp                                                              Shows all options (including advanc
                                                                          ed)

The following configuration subsections are valid:
 - extraction       Parameters for the channel extraction.
 - itraq4plex       Algorithm parameters for iTRAQ 4-plex
 - itraq8plex       Algorithm parameters for iTRAQ 8-plex
 - quantification   Parameters for the peptide quantification.
 - tmt10plex        Algorithm parameters for TMT 10-plex
 - tmt11plex        Algorithm parameters for TMT 11-plex
 - tmt16plex        Algorithm parameters for TMT 16-plex
 - tmt18plex        Algorithm parameters for TMT 18-plex
 - tmt6plex         Algorithm parameters for TMT 6-plex

You can write an example INI file using the '-write_ini' option.
Documentation of subsection parameters can be found in the doxygen documentation or the INIFileEditor.
For more information, please consult the online documentation for this tool:
  - http://www.openms.de/doxygen/nightly/html/TOPP_IsobaricWorkflow.html

INI file documentation of this tool:

Legend:
required parameter
advanced parameter
+IsobaricWorkflowCalculates isobaric quantitative values for peptides
version3.5.0-pre-nightly-2025-11-12 Version of the tool that generated this parameters file.
++1Instance '1' section for 'IsobaricWorkflow'
typeitraq4plex Isobaric Quantitation method used in the experiment.itraq4plex, itraq8plex, tmt10plex, tmt11plex, tmt16plex, tmt18plex, tmt6plex
in[] input centroided spectrum filesinput file*.mzML
in_id[] corresponding input PSMsinput file*.idXML
exp_design experimental design file (optional). If not given, the design is assumed to be unfractionated.input file*.tsv
out output consensusXML fileoutput file*.consensusXML
out_mzTab output mzTab file with quantitative informationoutput file*.mzTab
calculate_id_purityfalse Calculate the purity of the precursor ion based on the MS1 spectrum. Only used for MS3, otherwise it is the same as the quant. precursor purity.true, false
psm_scoreNaN The score which should be reached by a peptide hit to be kept. (use 'NAN' to disable this filter)
protein_scoreNaN The score which should be reached by a protein hit to be kept. All proteins are filtered based on their singleton scores irrespective of grouping. Use in combination with 'delete_unreferenced_peptide_hits' to remove affected peptides. (use 'NAN' to disable this filter)
delete_unreferenced_peptide_hitsfalse Peptides not referenced by any protein are deleted in the IDs.true, false
inference_methodaggregation Methods used for protein inferenceaggregation, bayesian
picked_fdrfalse Use a picked protein FDRtrue, false
picked_decoy_string If using picked protein FDRs, which decoy string was used? Leave blank for auto-detection.
picked_decoy_prefixprefix If using picked protein FDRs, was the decoy string a prefix or suffix? Ignored during auto-detection.prefix, suffix
FDR_typePSM Sub-protein FDR level. PSM, PSM+peptide (best PSM q-value).PSM, PSM+peptide
proteinFDR1.0 Protein FDR threshold (0.05=5%).0.0:1.0
psmFDR1.0 FDR threshold for sub-protein level (e.g. 0.05=5%).0.0:1.0
protein_quantificationunique_peptides Quantify proteins based on:
unique_peptides = use peptides mapping to single proteins or a group of indistinguishable proteins(according to the set of experimentally identified peptides).
strictly_unique_peptides = use peptides mapping to a unique single protein only.
shared_peptides = use shared peptides only for its best group (by inference score)
unique_peptides, strictly_unique_peptides, shared_peptides
log Name of log file (created only when specified)
debug0 Sets the debug level
threads1 Sets the number of threads allowed to be used by the TOPP tool
no_progressfalse Disables progress logging to command linetrue, false
forcefalse Overrides tool-specific checkstrue, false
testfalse Enables the test mode (needed for internal use only)true, false
+++ProteinQuantification
methodtop - top - quantify based on three most abundant peptides (number can be changed in 'top').
- iBAQ (intensity based absolute quantification), calculate the sum of all peptide peak intensities divided by the number of theoretically observable tryptic peptides (https://rdcu.be/cND1J). Warning: only consensusXML or featureXML input is allowed!
top, iBAQ
best_charge_and_fractionfalse Distinguish between fraction and charge states of a peptide. For peptides, abundances will be reported separately for each fraction and charge;
for proteins, abundances will be computed based only on the most prevalent charge observed of each peptide (over all fractions).
By default, abundances are summed over all charge states.
true, false
++++topAdditional options for custom quantification using top N peptides.
N3 Calculate protein abundance from this number of proteotypic peptides (most abundant first; '0' for all)0:∞
aggregatemedian Aggregation method used to compute protein abundances from peptide abundancesmedian, mean, weighted_mean, sum
include_alltrue Include results for proteins with fewer proteotypic peptides than indicated by 'N' (no effect if 'N' is 0 or 1)true, false
++++consensusAdditional options for consensus maps (and identification results comprising multiple runs)
normalizefalse Scale peptide abundances so that medians of all samples are equaltrue, false
fix_peptidesfalse Use the same peptides for protein quantification across all samples.
With 'N 0',all peptides that occur in every sample are considered.
Otherwise ('N'), the N peptides that occur in the most samples (independently of each other) are selected,
breaking ties by total abundance (there is no guarantee that the best co-ocurring peptides are chosen!).
true, false
+++BasicProteinInference
min_peptides_per_protein1 Minimal number of peptides needed for a protein identification. If set to zero, unmatched proteins get a score of -Infinity. If bigger than zero, proteins with less peptides are filtered and evidences removed from the PSMs. PSMs that do not reference any proteins anymore are removed but the spectrum info is kept.0:∞
score_aggregation_methodbest How to aggregate scores of peptides matching to the same protein?best, product, sum, maximum
treat_charge_variants_separatelytrue If this is true, different charge variants of the same peptide sequence count as individual evidences.true, false
treat_modification_variants_separatelytrue If this is true, different modification variants of the same peptide sequence count as individual evidences.true, false
use_shared_peptidestrue If this is true, shared peptides are used as evidences. Note: shared_peptides are not deleted and potentially resolved in postprocessing as well.true, false
skip_count_annotationfalse If this is set, peptide counts won't be annotated at the proteins.true, false
score_type PSM score type to use for inference. (default: empty = main score), PEP, q-value, RAW
+++BayesianProteinInference
psm_probability_cutoff1.0e-03 Remove PSMs with probabilities less than this cutoff0.0:1.0
top_PSMs1 Consider only top X PSMs per spectrum. 0 considers all.0:∞
keep_best_PSM_onlytrue Epifany uses the best PSM per peptide for inference. Discard the rest (true) or keepe.g. for quantification/reporting?true, false
++++model_parametersModel parameters for the Bayesian network
prot_prior-1.0 Protein prior probability ('gamma' parameter). Negative values enable grid search for this param.-1.0:1.0
pep_emission-1.0 Peptide emission probability ('alpha' parameter). Negative values enable grid search for this param.-1.0:1.0
pep_spurious_emission-1.0 Spurious peptide identification probability ('beta' parameter). Usually much smaller than emission from proteins. Negative values enable grid search for this param.-1.0:1.0
pep_prior0.1 Peptide prior probability (experimental, should be covered by combinations of the other params).0.0:1.0
regularizefalse Regularize the number of proteins that produce a peptide together (experimental, should be activated when using higher p-norms).true, false
extended_modelfalse Uses information from different peptidoforms also across runs (automatically activated if an experimental design is given!)true, false
++++loopy_belief_propagationSettings for the loopy belief propagation algorithm.
scheduling_typepriority (Not used yet) How to pick the next message: priority = based on difference to last message (higher = more important). fifo = first in first out. subtree = message passing follows a random spanning tree in each iterationpriority, fifo, subtree
convergence_threshold1.0e-05 Initial threshold under which MSE difference a message is considered to be converged.1.0e-09:1.0
dampening_lambda1.0e-03 Initial value for how strongly should messages be updated in each step. 0 = new message overwrites old completely (no dampening; only recommended for trees),0.5 = equal contribution of old and new message (stay below that),In-between it will be a convex combination of both. Prevents oscillations but hinders convergence.0.0:0.49999
max_nr_iterations2147483647 (Usually auto-determined by estimated but you can set a hard limit here). If not all messages converge, how many iterations should be done at max per connected component?
p_norm_inference1.0 P-norm used for marginalization of multidimensional factors. 1 == sum-product inference (all configurations vote equally) (default),<= 0 == infinity = max-product inference (only best configurations propagate)The higher the value the more important high probability configurations get.
++++param_optimizeSettings for the parameter optimization.
aucweight0.3 How important is target decoy AUC vs calibration of the posteriors? 0 = maximize calibration only, 1 = maximize AUC only, between = convex combination.0.0:1.0
conservative_fdrtrue Use (D+1)/(T) instead of (D+1)/(T+D) for parameter estimation.true, false
regularized_fdrtrue Use a regularized FDR for proteins without unique peptides.true, false
+++extractionParameters for the channel extraction.
select_activationauto Operate only on MSn scans where any of its precursors features a certain activation method. Setting to "auto" uses HCD and HCID spectra. Set to empty string if you want to disable filtering.auto, Collision-induced dissociation, Post-source decay, Plasma desorption, Surface-induced dissociation, Blackbody infrared radiative dissociation, Electron capture dissociation, Infrared multiphoton dissociation, Sustained off-resonance irradiation, High-energy collision-induced dissociation, Low-energy collision-induced dissociation, Photodissociation, Electron transfer dissociation, Electron transfer and collision-induced dissociation, Electron transfer and higher-energy collision dissociation, Pulsed q dissociation, trap-type collision-induced dissociation, beam-type collision-induced dissociation, in-source collision-induced dissociation, any
reporter_mass_shift2.0e-03 Allowed shift (left to right) in Th from the expected position.1.0e-04:0.5
min_precursor_intensity1.0 Minimum intensity of the precursor to be extracted. MS/MS scans having a precursor with a lower intensity will not be considered for quantitation.0.0:∞
keep_unannotated_precursortrue Flag if precursor with missing intensity value or missing precursor spectrum should be included or not.true, false
min_reporter_intensity0.0 Minimum intensity of the individual reporter ions to be extracted.0.0:∞
discard_low_intensity_quantificationsfalse Remove all reporter intensities if a single reporter is below the threshold given in 'min_reporter_intensity'.true, false
min_precursor_purity0.0 Minimum fraction of the total intensity in the isolation window of the precursor spectrum attributable to the selected precursor.0.0:1.0
precursor_isotope_deviation10.0 Maximum allowed deviation (in ppm) between theoretical and observed isotopic peaks of the precursor peak in the isolation window to be counted as part of the precursor.0.0:∞
purity_interpolationtrue If set to true the algorithm will try to compute the purity as a time weighted linear combination of the precursor scan and the following scan. If set to false, only the precursor scan will be used.true, false
+++itraq4plexAlgorithm parameters for iTRAQ 4-plex
channel_114_description Description for the content of the 114 channel.
channel_115_description Description for the content of the 115 channel.
channel_116_description Description for the content of the 116 channel.
channel_117_description Description for the content of the 117 channel.
reference_channel114 Number of the reference channel (114-117).114:117
correction_matrix[0.0/1.0/5.9/0.2, 0.0/2.0/5.6/0.1, 0.0/3.0/4.5/0.1, 0.1/4.0/3.5/0.1] Correction matrix for isotope distributions (see documentation); use the following format: <-2Da>/<-1Da>/<+1Da>/<+2Da>; e.g. '0/0.3/4/0', '0.1/0.3/3/0.2'
+++itraq8plexAlgorithm parameters for iTRAQ 8-plex
channel_113_description Description for the content of the 113 channel.
channel_114_description Description for the content of the 114 channel.
channel_115_description Description for the content of the 115 channel.
channel_116_description Description for the content of the 116 channel.
channel_117_description Description for the content of the 117 channel.
channel_118_description Description for the content of the 118 channel.
channel_119_description Description for the content of the 119 channel.
channel_121_description Description for the content of the 121 channel.
reference_channel113 Number of the reference channel (113-121). Please note that 120 is not valid.113:121
correction_matrix[0.00/0.00/6.89/0.22, 0.00/0.94/5.90/0.16, 0.00/1.88/4.90/0.10, 0.00/2.82/3.90/0.07, 0.06/3.77/2.99/0.00, 0.09/4.71/1.88/0.00, 0.14/5.66/0.87/0.00, 0.27/7.44/0.18/0.00] Correction matrix for isotope distributions (see documentation); use the following format: <-2Da>/<-1Da>/<+1Da>/<+2Da>; e.g. '0/0.3/4/0', '0.1/0.3/3/0.2'
+++quantificationParameters for the peptide quantification.
isotope_correctiontrue Enable isotope correction (highly recommended). Note that you need to provide a correct isotope correction matrix otherwise the tool will fail or produce invalid results.true, false
normalizationfalse Enable normalization of channel intensities with respect to the reference channel. The normalization is done by using the Median of Ratios (every channel / Reference). Also the ratio of medians (from any channel and reference) is provided as control measure!true, false
+++tmt10plexAlgorithm parameters for TMT 10-plex
channel_126_description Description for the content of the 126 channel.
channel_127N_description Description for the content of the 127N channel.
channel_127C_description Description for the content of the 127C channel.
channel_128N_description Description for the content of the 128N channel.
channel_128C_description Description for the content of the 128C channel.
channel_129N_description Description for the content of the 129N channel.
channel_129C_description Description for the content of the 129C channel.
channel_130N_description Description for the content of the 130N channel.
channel_130C_description Description for the content of the 130C channel.
channel_131_description Description for the content of the 131 channel.
reference_channel126 The reference channel (126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131).126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131
correction_matrix[0.0/0.0/5.09/0.0, 0.0/0.25/5.27/0.0, 0.0/0.37/5.36/0.15, 0.0/0.65/4.17/0.1, 0.08/0.49/3.06/0.0, 0.01/0.71/3.07/0.0, 0.0/1.32/2.62/0.0, 0.02/1.28/2.75/2.53, 0.03/2.08/2.23/0.0, 0.08/1.99/1.65/0.0] Correction matrix for isotope distributions (see documentation); use the following format: <-2Da>/<-1Da>/<+1Da>/<+2Da>; e.g. '0/0.3/4/0', '0.1/0.3/3/0.2'
+++tmt11plexAlgorithm parameters for TMT 11-plex
channel_126_description Description for the content of the 126 channel.
channel_127N_description Description for the content of the 127N channel.
channel_127C_description Description for the content of the 127C channel.
channel_128N_description Description for the content of the 128N channel.
channel_128C_description Description for the content of the 128C channel.
channel_129N_description Description for the content of the 129N channel.
channel_129C_description Description for the content of the 129C channel.
channel_130N_description Description for the content of the 130N channel.
channel_130C_description Description for the content of the 130C channel.
channel_131N_description Description for the content of the 131N channel.
channel_131C_description Description for the content of the 131C channel.
reference_channel126 The reference channel (126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C).126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C
correction_matrix[0.0/0.0/8.6/0.3, 0.0/0.1/7.8/0.1, 0.0/0.8/6.9/0.1, 0.0/7.4/7.4/0.0, 0.0/1.5/6.2/0.2, 0.0/1.5/5.7/0.1, 0.0/2.6/4.8/0.0, 0.0/2.2/4.6/0.0, 0.0/2.8/4.5/0.1, 0.1/2.9/3.8/0.0, 0.0/3.9/2.8/0.0] Correction matrix for isotope distributions (see documentation); use the following format: <-2Da>/<-1Da>/<+1Da>/<+2Da>; e.g. '0/0.3/4/0', '0.1/0.3/3/0.2'
+++tmt16plexAlgorithm parameters for TMT 16-plex
channel_126_description Description for the content of the 126 channel.
channel_127N_description Description for the content of the 127N channel.
channel_127C_description Description for the content of the 127C channel.
channel_128N_description Description for the content of the 128N channel.
channel_128C_description Description for the content of the 128C channel.
channel_129N_description Description for the content of the 129N channel.
channel_129C_description Description for the content of the 129C channel.
channel_130N_description Description for the content of the 130N channel.
channel_130C_description Description for the content of the 130C channel.
channel_131N_description Description for the content of the 131N channel.
channel_131C_description Description for the content of the 131C channel.
channel_132N_description Description for the content of the 132N channel.
channel_132C_description Description for the content of the 132C channel.
channel_133N_description Description for the content of the 133N channel.
channel_133C_description Description for the content of the 133C channel.
channel_134N_description Description for the content of the 134N channel.
reference_channel126 The reference channel (126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C, 132N, 132C, 133N, 133C, 134N).126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C, 132N, 132C, 133N, 133C, 134N
correction_matrix[NA/NA / NA/NA / 0.31/9.09 / 0.02/0.32, NA/NA / NA/0.78 / NA/9.41 / NA/0.33, NA/NA / 0.93/NA / 0.35/8.63 / 0.01/0.27, NA/0.00 / 0.82/0.65 / NA/8.13 / NA/0.26, 0.00/NA / 1.47/NA / 0.34/6.91 / 0.00/0.15, 0.00/0.00 / 1.46/1.28 / NA/6.86 / NA/0.15, 0.13/NA / 2.59/NA / 0.32/6.07 / 0.1/0.09, 0.13/0.00 / 2.41/0.27 / NA/5.58 / NA/0.10, 0.04/NA / 3.10/NA / 0.42/4.82 / 0.02/0.06, 0.03/0.00 / 2.78/0.63 / NA/4.57 / NA/0.12, 0.08/NA / 3.90/NA / 0.47/3.57 / 0.00/0.04, 0.15/0.01 / 3.58/0.72 / NA/1.80 / NA/0.00, 0.11/NA / 4.55/NA / 0.43/1.86 / 0.00/0.00, 0.07/0.01 / 3.14/0.73 / NA/3.40 / NA/0.03, 0.22/NA / 4.96/NA / 0.34/1.03 / 0.00/NA, 0.30/0.03 / 5.49/0.62 / NA/1.14 / NA/NA] Correction matrix for isotope distributions in percent from the Thermo data sheet (see documentation); Please provide 16 entries (rows), separated by comma, where each entry contains 8 values in the following format: <-2C13>/<-N15-C13>/<-C13>/<-N15>/<+N15>/<+C13>/<+N15+C13>/<+2C13> e.g. one row may look like this: 'NA/0.00 / 0.82/0.65 / NA/8.13 / NA/0.26'. You may use whitespaces at your leisure to ease reading.
+++tmt18plexAlgorithm parameters for TMT 18-plex
channel_126_description Description for the content of the 126 channel.
channel_127N_description Description for the content of the 127N channel.
channel_127C_description Description for the content of the 127C channel.
channel_128N_description Description for the content of the 128N channel.
channel_128C_description Description for the content of the 128C channel.
channel_129N_description Description for the content of the 129N channel.
channel_129C_description Description for the content of the 129C channel.
channel_130N_description Description for the content of the 130N channel.
channel_130C_description Description for the content of the 130C channel.
channel_131N_description Description for the content of the 131N channel.
channel_131C_description Description for the content of the 131C channel.
channel_132N_description Description for the content of the 132N channel.
channel_132C_description Description for the content of the 132C channel.
channel_133N_description Description for the content of the 133N channel.
channel_133C_description Description for the content of the 133C channel.
channel_134N_description Description for the content of the 134N channel.
channel_134C_description Description for the content of the 134C channel.
channel_135N_description Description for the content of the 135N channel.
reference_channel126 The reference channel (126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C, 132N, 132C, 133N, 133C, 134N, 134C, 135N).126, 127N, 127C, 128N, 128C, 129N, 129C, 130N, 130C, 131N, 131C, 132N, 132C, 133N, 133C, 134N, 134C, 135N
correction_matrix[NA/NA /NA/NA /0.31/9.09 /0.02/0.32, NA/NA /NA/0.78 /NA/9.41 /NA/0.33, NA/NA /0.93/NA /0.35/8.63 /0.01/0.27, NA/0.00 /0.82/0.65 /NA/8.13 /NA/0.26, 0.00/NA /1.47/NA /0.34/6.91 /0.00/0.15, 0.00/0.00 /1.46/1.28 /NA/6.86 /NA/0.15, 0.13/NA /2.59/NA /0.32/6.07 /0.1/0.09, 0.13/0.00 /2.41/0.27 /NA/5.58 /NA/0.10, 0.04/NA /3.10/NA /0.42/4.82 /0.02/0.06, 0.03/0.00 /2.78/0.63 /NA/4.57 /NA/0.12, 0.08/NA /3.90/NA /0.47/3.57 /0.00/0.04, 0.15/0.01 /3.58/0.72 /NA/1.80 /NA/0.00, 0.11/NA /4.55/NA /0.43/1.86 /0.00/0.00, 0.07/0.01 /3.14/0.73 /NA/3.40 /NA/0.03, 0.22/NA /4.96/NA /0.34/1.03 /0.00/NA, 0.30/0.03 /5.49/0.62 /NA/1.14 /NA/NA, 0.14/NA /5.81/NA /0.31/NA /NA/NA, 0.19/0.02 /5.42/0.36 /NA/NA /NA/NA] Correction matrix for isotope distributions in percent from the Thermo data sheet (see documentation); Please provide 18 entries (rows), separated by comma, where each entry contains 8 values in the following format: <-2C13>/<-N15-C13>/<-C13>/<-N15>/<+N15>/<+C13>/<+N15+C13>/<+2C13> e.g. one row may look like this: 'NA/0.00 / 0.82/0.65 / NA/8.13 / NA/0.26'. You may use whitespaces at your leisure to ease reading.
+++tmt6plexAlgorithm parameters for TMT 6-plex
channel_126_description Description for the content of the 126 channel.
channel_127_description Description for the content of the 127 channel.
channel_128_description Description for the content of the 128 channel.
channel_129_description Description for the content of the 129 channel.
channel_130_description Description for the content of the 130 channel.
channel_131_description Description for the content of the 131 channel.
reference_channel126 Number of the reference channel (126-131).126:131
correction_matrix[0.0/0.0/8.6/0.3, 0.0/0.1/7.8/0.1, 0.0/1.5/6.2/0.2, 0.0/1.5/5.7/0.1, 0.0/3.1/3.6/0.0, 0.1/2.9/3.8/0.0] Correction matrix for isotope distributions (see documentation); use the following format: <-2Da>/<-1Da>/<+1Da>/<+2Da>; e.g. '0/0.3/4/0', '0.1/0.3/3/0.2'