OpenMS  2.5.0
MRMTransitionGroupPicker.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2020.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Hannes Roest $
32 // $Authors: Hannes Roest $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
43 
46 
49 
50 // Cross-correlation
53 
54 #include <numeric>
55 
56 //#define DEBUG_TRANSITIONGROUPPICKER
57 
58 namespace OpenMS
59 {
60 
79  class OPENMS_DLLAPI MRMTransitionGroupPicker :
80  public DefaultParamHandler
81  {
82 
83 public:
84 
88 
90  ~MRMTransitionGroupPicker() override;
92 
112  template <typename SpectrumT, typename TransitionT>
114  {
115  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
116  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
117 
118  std::vector<MSChromatogram > picked_chroms;
119  std::vector<MSChromatogram > smoothed_chroms;
120 
121  // Pick fragment ion chromatograms
122  for (Size k = 0; k < transition_group.getChromatograms().size(); k++)
123  {
124  MSChromatogram& chromatogram = transition_group.getChromatograms()[k];
125  String native_id = chromatogram.getNativeID();
126 
127  // only pick detecting transitions (skip all others)
128  if (transition_group.getTransitions().size() > 0 &&
129  transition_group.hasTransition(native_id) &&
130  !transition_group.getTransition(native_id).isDetectingTransition() )
131  {
132  continue;
133  }
134 
135  MSChromatogram picked_chrom, smoothed_chrom;
136  picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
137  picked_chrom.sortByIntensity();
138  picked_chroms.push_back(std::move(picked_chrom));
139  smoothed_chroms.push_back(std::move(smoothed_chrom));
140  }
141 
142  // Pick precursor chromatograms
143  if (use_precursors_)
144  {
145  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
146  {
147  SpectrumT picked_chrom, smoothed_chrom;
148  SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
149 
150  picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
151  picked_chrom.sortByIntensity();
152  picked_chroms.push_back(picked_chrom);
153  smoothed_chroms.push_back(smoothed_chrom);
154  }
155  }
156 
157  // Find features (peak groups) in this group of transitions.
158  // While there are still peaks left, one will be picked and used to create
159  // a feature. Whenever we run out of peaks, we will get -1 back as index
160  // and terminate.
161  int chr_idx, peak_idx, cnt = 0;
162  std::vector<MRMFeature> features;
163  while (true)
164  {
165  chr_idx = -1; peak_idx = -1;
166 
167  if (boundary_selection_method_ == "largest")
168  {
169  findLargestPeak(picked_chroms, chr_idx, peak_idx);
170  }
171  else if (boundary_selection_method_ == "widest")
172  {
173  findWidestPeakIndices(picked_chroms, chr_idx, peak_idx);
174  }
175 
176  if (chr_idx == -1 && peak_idx == -1) break;
177 
178  // Compute a feature from the individual chromatograms and add non-zero features
179  MRMFeature mrm_feature = createMRMFeature(transition_group, picked_chroms, smoothed_chroms, chr_idx, peak_idx);
180  double total_xic = 0;
181  double intensity = mrm_feature.getIntensity();
182  if (intensity > 0)
183  {
184  total_xic = mrm_feature.getMetaValue("total_xic");
185  features.push_back(std::move(mrm_feature));
186  }
187 
188  cnt++;
189  if (stop_after_feature_ > 0 && cnt > stop_after_feature_) {break;}
190  if (intensity > 0 && intensity / total_xic < stop_after_intensity_ratio_)
191  {
192  break;
193  }
194  }
195 
196  // Check for completely overlapping features
197  for (Size i = 0; i < features.size(); i++)
198  {
199  MRMFeature& mrm_feature = features[i];
200  bool skip = false;
201  for (Size j = 0; j < i; j++)
202  {
203  if ((double)mrm_feature.getMetaValue("leftWidth") >= (double)features[j].getMetaValue("leftWidth") &&
204  (double)mrm_feature.getMetaValue("rightWidth") <= (double)features[j].getMetaValue("rightWidth") )
205  { skip = true; }
206  }
207  if (mrm_feature.getIntensity() > 0 && !skip)
208  {
209  transition_group.addFeature(mrm_feature);
210  }
211  }
212 
213  }
214 
216  template <typename SpectrumT, typename TransitionT>
218  std::vector<SpectrumT>& picked_chroms,
219  const std::vector<SpectrumT>& smoothed_chroms,
220  const int chr_idx,
221  const int peak_idx)
222  {
223  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
224  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
225 
226  MRMFeature mrmFeature;
227  mrmFeature.setIntensity(0.0);
228  double best_left = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][peak_idx];
229  double best_right = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][peak_idx];
230  double peak_apex = picked_chroms[chr_idx][peak_idx].getRT();
231  OPENMS_LOG_DEBUG << "**** Creating MRMFeature for peak " << chr_idx << " " << peak_idx << " " <<
232  picked_chroms[chr_idx][peak_idx] << " with borders " << best_left << " " <<
233  best_right << " (" << best_right - best_left << ")" << std::endl;
234 
235  if (use_consensus_ && recalculate_peaks_)
236  {
237  // This may change best_left / best_right
238  recalculatePeakBorders_(picked_chroms, best_left, best_right, recalculate_peaks_max_z_);
239  if (peak_apex < best_left || peak_apex > best_right)
240  {
241  // apex fell out of range, lets correct it
242  peak_apex = (best_left + best_right) / 2.0;
243  }
244  }
245 
246  std::vector< double > left_edges;
247  std::vector< double > right_edges;
248  double min_left = best_left;
249  double max_right = best_right;
250  if (use_consensus_)
251  {
252  // Remove other, overlapping, picked peaks (in this and other
253  // chromatograms) and then ensure that at least one peak is set to zero
254  // (the currently best peak).
255  remove_overlapping_features(picked_chroms, best_left, best_right);
256  }
257  else
258  {
259  pickApex(picked_chroms, best_left, best_right, peak_apex,
260  min_left, max_right, left_edges, right_edges);
261 
262  } // end !use_consensus_
263  picked_chroms[chr_idx][peak_idx].setIntensity(0.0); // ensure that we set at least one peak to zero
264 
265  // Check for minimal peak width -> return empty feature (Intensity zero)
266  if (use_consensus_)
267  {
268  if (min_peak_width_ > 0.0 && std::fabs(best_right - best_left) < min_peak_width_)
269  {
270  return mrmFeature;
271  }
272 
273  if (compute_peak_quality_)
274  {
275  String outlier = "none";
276  double qual = computeQuality_(transition_group, picked_chroms, chr_idx, best_left, best_right, outlier);
277  if (qual < min_qual_)
278  {
279  return mrmFeature;
280  }
281  mrmFeature.setMetaValue("potentialOutlier", outlier);
282  mrmFeature.setMetaValue("initialPeakQuality", qual);
283  mrmFeature.setOverallQuality(qual);
284  }
285  }
286 
287  // Prepare linear resampling of all the chromatograms, here creating the
288  // empty master_peak_container with the same RT (m/z) values as the
289  // reference chromatogram. We use the overall minimal left boundary and
290  // maximal right boundary to prepare the container.
291  SpectrumT master_peak_container;
292  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
293  prepareMasterContainer_(ref_chromatogram, master_peak_container, min_left, max_right);
294 
295  // Iterate over initial transitions / chromatograms (note that we may
296  // have a different number of picked chromatograms than total transitions
297  // as not all are detecting transitions).
298  double total_intensity = 0; double total_peak_apices = 0; double total_xic = 0; double total_mi = 0;
299  pickFragmentChromatograms(transition_group, picked_chroms, mrmFeature, smoothed_chroms,
300  best_left, best_right, use_consensus_,
301  total_intensity, total_xic, total_mi, total_peak_apices,
302  master_peak_container, left_edges, right_edges,
303  chr_idx, peak_idx);
304 
305  // Also pick the precursor chromatogram(s); note total_xic is not
306  // extracted here, only for fragment traces
307  pickPrecursorChromatograms(transition_group,
308  picked_chroms, mrmFeature, smoothed_chroms,
309  best_left, best_right, use_consensus_,
310  total_intensity, master_peak_container, left_edges, right_edges,
311  chr_idx, peak_idx);
312 
313  mrmFeature.setRT(peak_apex);
314  mrmFeature.setIntensity(total_intensity);
315  mrmFeature.setMetaValue("PeptideRef", transition_group.getTransitionGroupID());
316  mrmFeature.setMetaValue("leftWidth", best_left);
317  mrmFeature.setMetaValue("rightWidth", best_right);
318  mrmFeature.setMetaValue("total_xic", total_xic);
319  if (compute_total_mi_)
320  {
321  mrmFeature.setMetaValue("total_mi", total_mi);
322  }
323  mrmFeature.setMetaValue("peak_apices_sum", total_peak_apices);
324 
325  mrmFeature.ensureUniqueId();
326  return mrmFeature;
327  }
328 
341  template <typename SpectrumT>
342  void pickApex(std::vector<SpectrumT>& picked_chroms,
343  const double best_left, const double best_right, const double peak_apex,
344  double &min_left, double &max_right,
345  std::vector< double > & left_edges, std::vector< double > & right_edges)
346  {
347  for (Size k = 0; k < picked_chroms.size(); k++)
348  {
349  double peak_apex_dist_min = std::numeric_limits<double>::max();
350  int min_dist = -1;
351  for (Size i = 0; i < picked_chroms[k].size(); i++)
352  {
353  PeakIntegrator::PeakArea pa_tmp = pi_.integratePeak( // get the peak apex
354  picked_chroms[k],
355  picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i],
356  picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i]);
357  if (pa_tmp.apex_pos > 0.0 && std::fabs(pa_tmp.apex_pos - peak_apex) < peak_apex_dist_min)
358  { // update best candidate
359  peak_apex_dist_min = std::fabs(pa_tmp.apex_pos - peak_apex);
360  min_dist = (int)i;
361  }
362  }
363 
364  // Select master peak boundaries, or in the case we found at least one peak, the local peak boundaries
365  double l = best_left;
366  double r = best_right;
367  if (min_dist >= 0)
368  {
369  l = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][min_dist];
370  r = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][min_dist];
371  picked_chroms[k][min_dist].setIntensity(0.0); // only remove one peak per transition
372  }
373 
374  left_edges.push_back(l);
375  right_edges.push_back(r);
376  // ensure we remember the overall maxima / minima
377  if (l < min_left) {min_left = l;}
378  if (r > max_right) {max_right = r;}
379  }
380  }
381 
382  template <typename SpectrumT, typename TransitionT>
384  const std::vector<SpectrumT>& picked_chroms,
385  MRMFeature& mrmFeature,
386  const std::vector<SpectrumT>& smoothed_chroms,
387  const double best_left, const double best_right,
388  const bool use_consensus_,
389  double & total_intensity,
390  double & total_xic,
391  double & total_mi,
392  double & total_peak_apices,
393  const SpectrumT & master_peak_container,
394  const std::vector< double > & left_edges,
395  const std::vector< double > & right_edges,
396  const int chr_idx,
397  const int peak_idx)
398  {
399  for (Size k = 0; k < transition_group.getTransitions().size(); k++)
400  {
401 
402  double local_left = best_left;
403  double local_right = best_right;
404  if (!use_consensus_)
405  {
406  // We cannot have any non-detecting transitions (otherwise we have
407  // too few left / right edges) as we skipped those when doing peak
408  // picking and smoothing.
409  if (!transition_group.getTransitions()[k].isDetectingTransition())
410  {
411  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
412  "When using non-censensus peak picker, all transitions need to be detecting transitions.");
413  }
414  local_left = left_edges[k];
415  local_right = right_edges[k];
416  }
417 
418  const SpectrumT& chromatogram = selectChromHelper_(transition_group, transition_group.getTransitions()[k].getNativeID());
419  if (transition_group.getTransitions()[k].isDetectingTransition())
420  {
421  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
422  {
423  total_xic += it->getIntensity();
424  }
425  }
426 
427  // Compute total intensity on transition-level
428  double transition_total_xic = 0;
429 
430  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
431  {
432  transition_total_xic += it->getIntensity();
433  }
434 
435  // Compute total mutual information on transition-level.
436  double transition_total_mi = 0;
437  if (compute_total_mi_)
438  {
439  std::vector<double> chrom_vect_id, chrom_vect_det;
440  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
441  {
442  chrom_vect_id.push_back(it->getIntensity());
443  }
444 
445  // compute baseline mutual information
446  int transition_total_mi_norm = 0;
447  for (Size m = 0; m < transition_group.getTransitions().size(); m++)
448  {
449  if (transition_group.getTransitions()[m].isDetectingTransition())
450  {
451  const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.getTransitions()[m].getNativeID());
452  chrom_vect_det.clear();
453  for (typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
454  {
455  chrom_vect_det.push_back(it->getIntensity());
456  }
457  transition_total_mi += OpenSwath::Scoring::rankedMutualInformation(chrom_vect_det, chrom_vect_id);
458  transition_total_mi_norm++;
459  }
460  }
461  if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
462 
463  if (transition_group.getTransitions()[k].isDetectingTransition())
464  {
465  // sum up all transition-level total MI and divide by the number of detection transitions to have peak group level total MI
466  total_mi += transition_total_mi / transition_total_mi_norm;
467  }
468  }
469 
470  SpectrumT used_chromatogram;
471  // resample the current chromatogram
472  if (peak_integration_ == "original")
473  {
474  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
475  }
476  else if (peak_integration_ == "smoothed")
477  {
478  if (smoothed_chroms.size() <= k)
479  {
480  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
481  "Tried to calculate peak area and height without any smoothed chromatograms");
482  }
483  used_chromatogram = resampleChromatogram_(smoothed_chroms[k], master_peak_container, local_left, local_right);
484  }
485  else
486  {
487  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
488  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
489  }
490 
491  Feature f;
492  double quality = 0;
493  f.setQuality(0, quality);
494  f.setOverallQuality(quality);
495 
496  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
497  double peak_integral = pa.area;
498  double peak_apex_int = pa.height;
499  f.setMetaValue("peak_apex_position", pa.apex_pos);
500  if (background_subtraction_ != "none")
501  {
502  double background{0};
503  double avg_noise_level{0};
504  if (background_subtraction_ == "original")
505  {
506  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
507  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
508  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
509  avg_noise_level = (intensity_right + intensity_left) / 2;
510  background = avg_noise_level * n_points;
511  }
512  else if (background_subtraction_ == "exact")
513  {
514  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
515  background = pb.area;
516  avg_noise_level = pb.height;
517  }
518  peak_integral -= background;
519  peak_apex_int -= avg_noise_level;
520  if (peak_integral < 0) {peak_integral = 0;}
521  if (peak_apex_int < 0) {peak_apex_int = 0;}
522 
523  f.setMetaValue("area_background_level", background);
524  f.setMetaValue("noise_background_level", avg_noise_level);
525  } // end background
526 
527  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
528  f.setIntensity(peak_integral);
529  ConvexHull2D hull;
530  hull.setHullPoints(pa.hull_points);
531  f.getConvexHulls().push_back(hull);
532 
533  f.setMZ(chromatogram.getProduct().getMZ());
534  mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());
535 
536  if (chromatogram.metaValueExists("product_mz")) // legacy code (ensures that old tests still work)
537  {
538  f.setMetaValue("MZ", chromatogram.getMetaValue("product_mz"));
539  f.setMZ(chromatogram.getMetaValue("product_mz"));
540  }
541 
542  f.setMetaValue("native_id", chromatogram.getNativeID());
543  f.setMetaValue("peak_apex_int", peak_apex_int);
544  f.setMetaValue("total_xic", transition_total_xic);
545  if (compute_total_mi_)
546  {
547  f.setMetaValue("total_mi", transition_total_mi);
548  }
549 
550  if (transition_group.getTransitions()[k].isQuantifyingTransition())
551  {
552  total_intensity += peak_integral;
553  total_peak_apices += peak_apex_int;
554  }
555 
556  // for backwards compatibility with TOPP tests
557  // Calculate peak shape metrics that will be used for later QC
558  PeakIntegrator::PeakShapeMetrics psm = pi_.calculatePeakShapeMetrics(used_chromatogram, local_left, local_right, peak_apex_int, pa.apex_pos);
559  f.setMetaValue("width_at_50", psm.width_at_50);
560  if (compute_peak_shape_metrics_)
561  {
562  f.setMetaValue("width_at_5", psm.width_at_5);
563  f.setMetaValue("width_at_10", psm.width_at_10);
564  f.setMetaValue("start_position_at_5", psm.start_position_at_5);
565  f.setMetaValue("start_position_at_10", psm.start_position_at_10);
566  f.setMetaValue("start_position_at_50", psm.start_position_at_50);
567  f.setMetaValue("end_position_at_5", psm.end_position_at_5);
568  f.setMetaValue("end_position_at_10", psm.end_position_at_10);
569  f.setMetaValue("end_position_at_50", psm.end_position_at_50);
570  f.setMetaValue("total_width", psm.total_width);
571  f.setMetaValue("tailing_factor", psm.tailing_factor);
572  f.setMetaValue("asymmetry_factor", psm.asymmetry_factor);
573  f.setMetaValue("slope_of_baseline", psm.slope_of_baseline);
574  f.setMetaValue("baseline_delta_2_height", psm.baseline_delta_2_height);
575  f.setMetaValue("points_across_baseline", psm.points_across_baseline);
576  f.setMetaValue("points_across_half_height", psm.points_across_half_height);
577  }
578 
579  mrmFeature.addFeature(f, chromatogram.getNativeID()); //map index and feature
580  }
581  }
582 
583  template <typename SpectrumT, typename TransitionT>
585  const std::vector<SpectrumT>& picked_chroms,
586  MRMFeature& mrmFeature,
587  const std::vector<SpectrumT>& smoothed_chroms,
588  const double best_left, const double best_right,
589  const bool use_consensus_,
590  double & total_intensity,
591  const SpectrumT & master_peak_container,
592  const std::vector< double > & left_edges,
593  const std::vector< double > & right_edges,
594  const int chr_idx,
595  const int peak_idx)
596  {
597  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
598  {
599  const SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
600 
601  // Identify precursor index
602  // note: this is only valid if all transitions are detecting transitions
603  Size prec_idx = transition_group.getChromatograms().size() + k;
604 
605  double local_left = best_left;
606  double local_right = best_right;
607  if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
608  {
609  local_left = left_edges[prec_idx];
610  local_right = right_edges[prec_idx];
611  }
612 
613  SpectrumT used_chromatogram;
614  // resample the current chromatogram
615  if (peak_integration_ == "original")
616  {
617  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
618  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
619  }
620  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= prec_idx)
621  {
622  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
623  "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
624  }
625  else if (peak_integration_ == "smoothed")
626  {
627  used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
628  }
629  else
630  {
631  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
632  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
633  }
634 
635  Feature f;
636  double quality = 0;
637  f.setQuality(0, quality);
638  f.setOverallQuality(quality);
639 
640  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
641  double peak_integral = pa.area;
642  double peak_apex_int = pa.height;
643 
644  if (background_subtraction_ != "none")
645  {
646  double background{0};
647  double avg_noise_level{0};
648  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= prec_idx)
649  {
650  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
651  "Tried to calculate background estimation without any smoothed chromatograms");
652  }
653  else if (background_subtraction_ == "original")
654  {
655  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
656  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
657  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
658  avg_noise_level = (intensity_right + intensity_left) / 2;
659  background = avg_noise_level * n_points;
660  }
661  else if (background_subtraction_ == "exact")
662  {
663  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
664  background = pb.area;
665  avg_noise_level = pb.height;
666  }
667  peak_integral -= background;
668  peak_apex_int -= avg_noise_level;
669  if (peak_integral < 0) {peak_integral = 0;}
670  if (peak_apex_int < 0) {peak_apex_int = 0;}
671 
672  f.setMetaValue("area_background_level", background);
673  f.setMetaValue("noise_background_level", avg_noise_level);
674  }
675 
676  f.setMZ(chromatogram.getPrecursor().getMZ());
677  if (k == 0) {mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());} // only use m/z if first (monoisotopic) isotope
678 
679  if (chromatogram.metaValueExists("precursor_mz")) // legacy code (ensures that old tests still work)
680  {
681  f.setMZ(chromatogram.getMetaValue("precursor_mz"));
682  if (k == 0) {mrmFeature.setMZ(chromatogram.getMetaValue("precursor_mz"));} // only use m/z if first (monoisotopic) isotope
683  }
684 
685  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
686  f.setIntensity(peak_integral);
687  ConvexHull2D hull;
688  hull.setHullPoints(pa.hull_points);
689  f.getConvexHulls().push_back(hull);
690  f.setMetaValue("native_id", chromatogram.getNativeID());
691  f.setMetaValue("peak_apex_int", peak_apex_int);
692 
693  if (use_precursors_ && transition_group.getTransitions().empty())
694  {
695  total_intensity += peak_integral;
696  }
697 
698  mrmFeature.addPrecursorFeature(f, chromatogram.getNativeID());
699  }
700  }
701 
702  // maybe private, but we have tests
703 
715  template <typename SpectrumT>
716  void remove_overlapping_features(std::vector<SpectrumT>& picked_chroms, double best_left, double best_right)
717  {
718  // delete all seeds that lie within the current seed
719  for (Size k = 0; k < picked_chroms.size(); k++)
720  {
721  for (Size i = 0; i < picked_chroms[k].size(); i++)
722  {
723  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
724  {
725  picked_chroms[k][i].setIntensity(0.0);
726  }
727  }
728  }
729 
730  // delete all seeds that overlap within the current seed
731  for (Size k = 0; k < picked_chroms.size(); k++)
732  {
733  for (Size i = 0; i < picked_chroms[k].size(); i++)
734  {
735  if (picked_chroms[k][i].getIntensity() <= 0.0) {continue; }
736 
737  double left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
738  double right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
739  if ((left > best_left && left < best_right)
740  || (right > best_left && right < best_right))
741  {
742  picked_chroms[k][i].setIntensity(0.0);
743  }
744  }
745  }
746  }
747 
749  void findLargestPeak(const std::vector<MSChromatogram >& picked_chroms, int& chr_idx, int& peak_idx);
750 
759  void findWidestPeakIndices(const std::vector<MSChromatogram>& picked_chroms, Int& chrom_idx, Int& point_idx) const;
760 
761 protected:
762 
764  void updateMembers_() override;
765 
767  MRMTransitionGroupPicker& operator=(const MRMTransitionGroupPicker& rhs);
768 
772  template <typename SpectrumT, typename TransitionT>
773  const SpectrumT& selectChromHelper_(const MRMTransitionGroup<SpectrumT, TransitionT>& transition_group, const String& native_id)
774  {
775  if (transition_group.hasChromatogram(native_id))
776  {
777  return transition_group.getChromatogram(native_id);
778  }
779  else if (transition_group.hasPrecursorChromatogram(native_id))
780  {
781  return transition_group.getPrecursorChromatogram(native_id);
782  }
783  else
784  {
785  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION, "Did not find chromatogram for id '" + native_id + "'.");
786  }
787  }
788 
805  template <typename SpectrumT, typename TransitionT>
807  const std::vector<SpectrumT>& picked_chroms,
808  const int chr_idx,
809  const double best_left,
810  const double best_right,
811  String& outlier)
812  {
813  // Resample all chromatograms around the current estimated peak and
814  // collect the raw intensities. For resampling, use a bit more on either
815  // side to correctly identify shoulders etc.
816  double resample_boundary = resample_boundary_; // sample 15 seconds more on each side
817  SpectrumT master_peak_container;
818  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
819  prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
820  std::vector<std::vector<double> > all_ints;
821  for (Size k = 0; k < picked_chroms.size(); k++)
822  {
823  const SpectrumT& chromatogram = selectChromHelper_(transition_group, picked_chroms[k].getNativeID());
824  const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
825  master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
826 
827  std::vector<double> int_here;
828  for (const auto& peak : used_chromatogram) int_here.push_back(peak.getIntensity());
829  // Remove chromatograms without a single peak
830  double tic = std::accumulate(int_here.begin(), int_here.end(), 0.0);
831  if (tic > 0.0) all_ints.push_back(int_here);
832  }
833 
834  // Compute the cross-correlation for the collected intensities
835  std::vector<double> mean_shapes;
836  std::vector<double> mean_coel;
837  for (Size k = 0; k < all_ints.size(); k++)
838  {
839  std::vector<double> shapes;
840  std::vector<double> coel;
841  for (Size i = 0; i < all_ints.size(); i++)
842  {
843  if (i == k) {continue;}
845  all_ints[k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
846 
847  // the first value is the x-axis (retention time) and should be an int -> it show the lag between the two
848  double res_coelution = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->first);
849  double res_shape = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->second);
850 
851  shapes.push_back(res_shape);
852  coel.push_back(res_coelution);
853  }
854 
855  // We have computed the cross-correlation of chromatogram k against
856  // all others. Use the mean of these computations as the value for k.
858  msc = std::for_each(shapes.begin(), shapes.end(), msc);
859  double shapes_mean = msc.mean();
860  msc = std::for_each(coel.begin(), coel.end(), msc);
861  double coel_mean = msc.mean();
862 
863  // mean shape scores below 0.5-0.6 should be a real sign of trouble ... !
864  // mean coel scores above 3.5 should be a real sign of trouble ... !
865  mean_shapes.push_back(shapes_mean);
866  mean_coel.push_back(coel_mean);
867  }
868 
869  // find the chromatogram with the minimal shape score and the maximal
870  // coelution score -> if it is the same chromatogram, the chance is
871  // pretty good that it is different from the others...
872  int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
873  int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
874 
875  // Look at the picked peaks that are within the current left/right borders
876  int missing_peaks = 0;
877  int multiple_peaks = 0;
878 
879  // collect all seeds that lie within the current seed
880  std::vector<double> left_borders;
881  std::vector<double> right_borders;
882  for (Size k = 0; k < picked_chroms.size(); k++)
883  {
884  double l_tmp;
885  double r_tmp;
886  double max_int = -1;
887 
888  int pfound = 0;
889  l_tmp = -1;
890  r_tmp = -1;
891  for (Size i = 0; i < picked_chroms[k].size(); i++)
892  {
893  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
894  {
895  pfound++;
896  if (picked_chroms[k][i].getIntensity() > max_int)
897  {
898  max_int = picked_chroms[k][i].getIntensity() > max_int;
899  l_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
900  r_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
901  }
902  }
903  }
904 
905  if (l_tmp > 0.0) left_borders.push_back(l_tmp);
906  if (r_tmp > 0.0) right_borders.push_back(r_tmp);
907 
908  if (pfound == 0) missing_peaks++;
909  if (pfound > 1) multiple_peaks++;
910  }
911 
912  // Check how many chromatograms had exactly one peak picked between our
913  // current left/right borders -> this would be a sign of consistency.
914  OPENMS_LOG_DEBUG << " Overall found missing : " << missing_peaks << " and multiple : " << multiple_peaks << std::endl;
915 
917 
918  // Is there one transitions that is very different from the rest (e.g.
919  // the same element has a bad shape and a bad coelution score) -> potential outlier
920  if (min_index_shape == max_index_coel)
921  {
922  OPENMS_LOG_DEBUG << " element " << min_index_shape << " is a candidate for removal ... " << std::endl;
923  outlier = String(picked_chroms[min_index_shape].getNativeID());
924  }
925  else
926  {
927  outlier = "none";
928  }
929 
930  // For the final score (larger is better), consider these scores:
931  // - missing_peaks (the more peaks are missing, the worse)
932  // - multiple_peaks
933  // - mean of the shapes (1 is very good, 0 is bad)
934  // - mean of the co-elution scores (0 is good, 1 is ok, above 1 is pretty bad)
935  double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
936  double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
937  coel_score = (coel_score - 1.0) / 2.0;
938 
939  double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
940 
941  OPENMS_LOG_DEBUG << " computed score " << score << " (from " << shape_score <<
942  " - " << coel_score << " - " << 1.0 * missing_peaks / picked_chroms.size() << ")" << std::endl;
943 
944  return score;
945  }
946 
956  template <typename SpectrumT>
957  void recalculatePeakBorders_(const std::vector<SpectrumT>& picked_chroms, double& best_left, double& best_right, double max_z)
958  {
959  // 1. Collect all seeds that lie within the current seed
960  // - Per chromatogram only the most intense one counts, otherwise very
961  // - low intense peaks can contribute disproportionally to the voting
962  // - procedure.
963  std::vector<double> left_borders;
964  std::vector<double> right_borders;
965  for (Size k = 0; k < picked_chroms.size(); k++)
966  {
967  double max_int = -1;
968  double left = -1;
969  double right = -1;
970  for (Size i = 0; i < picked_chroms[k].size(); i++)
971  {
972  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
973  {
974  if (picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i] > max_int)
975  {
976  max_int = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i];
977  left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
978  right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
979  }
980  }
981  }
982  if (max_int > -1 )
983  {
984  left_borders.push_back(left);
985  right_borders.push_back(right);
986  OPENMS_LOG_DEBUG << " * " << k << " left boundary " << left_borders.back() << " with int " << max_int << std::endl;
987  OPENMS_LOG_DEBUG << " * " << k << " right boundary " << right_borders.back() << " with int " << max_int << std::endl;
988  }
989  }
990 
991  // Return for empty peak list
992  if (right_borders.empty())
993  {
994  return;
995  }
996 
997  // FEATURE IDEA: instead of Z-score use modified Z-score for small data sets
998  // http://d-scholarship.pitt.edu/7948/1/Seo.pdf
999  // http://www.itl.nist.gov/div898/handbook/eda/section3/eda35h.htm
1000  // 1. calculate median
1001  // 2. MAD = calculate difference to median for each value -> take median of that
1002  // 3. Mi = 0.6745*(xi - median) / MAD
1003 
1004  // 2. Calculate mean and standard deviation
1005  // If the coefficient of variation is too large for one border, we use a
1006  // "pseudo-median" instead of the border of the most intense peak.
1007  double mean, stdev;
1008 
1009  // Right borders
1010  mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
1011  stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
1012  / right_borders.size() - mean * mean);
1013  std::sort(right_borders.begin(), right_borders.end());
1014 
1015  OPENMS_LOG_DEBUG << " - Recalculating right peak boundaries " << mean << " mean / best "
1016  << best_right << " std " << stdev << " : " << std::fabs(best_right - mean) / stdev
1017  << " coefficient of variation" << std::endl;
1018 
1019  // Compare right borders of best transition with the mean
1020  if (std::fabs(best_right - mean) / stdev > max_z)
1021  {
1022  best_right = right_borders[right_borders.size() / 2]; // pseudo median
1023  OPENMS_LOG_DEBUG << " - Setting right boundary to " << best_right << std::endl;
1024  }
1025 
1026  // Left borders
1027  mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
1028  stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
1029  / left_borders.size() - mean * mean);
1030  std::sort(left_borders.begin(), left_borders.end());
1031 
1032  OPENMS_LOG_DEBUG << " - Recalculating left peak boundaries " << mean << " mean / best "
1033  << best_left << " std " << stdev << " : " << std::fabs(best_left - mean) / stdev
1034  << " coefficient of variation" << std::endl;
1035 
1036  // Compare left borders of best transition with the mean
1037  if (std::fabs(best_left - mean) / stdev > max_z)
1038  {
1039  best_left = left_borders[left_borders.size() / 2]; // pseudo median
1040  OPENMS_LOG_DEBUG << " - Setting left boundary to " << best_left << std::endl;
1041  }
1042 
1043  }
1044 
1046 
1047 
1062  template <typename SpectrumT>
1063  void prepareMasterContainer_(const SpectrumT& ref_chromatogram,
1064  SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1065  {
1066  OPENMS_PRECONDITION(master_peak_container.empty(), "Master peak container must be empty")
1067 
1068  // get the start / end point of this chromatogram => then add one more
1069  // point beyond the two boundaries to make the resampling accurate also
1070  // at the edge.
1071  typename SpectrumT::const_iterator begin = ref_chromatogram.begin();
1072  while (begin != ref_chromatogram.end() && begin->getMZ() < left_boundary) {begin++; }
1073  if (begin != ref_chromatogram.begin()) {begin--; }
1074 
1075  typename SpectrumT::const_iterator end = begin;
1076  while (end != ref_chromatogram.end() && end->getMZ() < right_boundary) {end++; }
1077  if (end != ref_chromatogram.end()) {end++; }
1078 
1079  // resize the master container and set the m/z values to the ones of the master container
1080  master_peak_container.resize(distance(begin, end)); // initialize to zero
1081  typename SpectrumT::iterator it = master_peak_container.begin();
1082  for (typename SpectrumT::const_iterator chrom_it = begin; chrom_it != end; chrom_it++, it++)
1083  {
1084  it->setMZ(chrom_it->getMZ());
1085  }
1086  }
1087 
1098  template <typename SpectrumT>
1099  SpectrumT resampleChromatogram_(const SpectrumT& chromatogram,
1100  const SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1101  {
1102  // get the start / end point of this chromatogram => then add one more
1103  // point beyond the two boundaries to make the resampling accurate also
1104  // at the edge.
1105  typename SpectrumT::const_iterator begin = chromatogram.begin();
1106  while (begin != chromatogram.end() && begin->getMZ() < left_boundary) {begin++;}
1107  if (begin != chromatogram.begin()) {begin--;}
1108 
1109  typename SpectrumT::const_iterator end = begin;
1110  while (end != chromatogram.end() && end->getMZ() < right_boundary) {end++;}
1111  if (end != chromatogram.end()) {end++;}
1112 
1113  SpectrumT resampled_peak_container = master_peak_container; // copy the master container, which contains the RT values
1114  LinearResamplerAlign lresampler;
1115  lresampler.raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1116 
1117  return resampled_peak_container;
1118  }
1119 
1121 
1122  // Members
1131  double min_qual_;
1132 
1138 
1145 
1148  };
1149 }
1150 
1151 
OpenMS::MRMTransitionGroup::hasPrecursorChromatogram
bool hasPrecursorChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:243
OpenMS::Param
Management and storage of parameters / INI files.
Definition: Param.h:73
OpenMS::MRMTransitionGroup::isInternallyConsistent
bool isInternallyConsistent() const
Check whether internal state is consistent, e.g. same number of chromatograms and transitions are pre...
Definition: MRMTransitionGroup.h:292
TargetedExperiment.h
OpenMS::MRMTransitionGroup::chromatogramIdsMatch
bool chromatogramIdsMatch() const
Ensure that chromatogram native ids match their keys in the map.
Definition: MRMTransitionGroup.h:301
OpenMS::FeatureXMLFile
This class provides Input/Output functionality for feature maps.
Definition: FeatureXMLFile.h:68
StatsHelpers.h
OpenMS::MRMTransitionGroup::addFeature
void addFeature(const MRMFeature &feature)
Definition: MRMTransitionGroup.h:276
OpenSwath::Scoring::normalizedCrossCorrelation
OPENSWATHALGO_DLLAPI XCorrArrayType normalizedCrossCorrelation(std::vector< double > &data1, std::vector< double > &data2, const int &maxdelay, const int &lag)
OpenMS::MRMTransitionGroup::hasTransition
bool hasTransition(String key) const
Definition: MRMTransitionGroup.h:158
OpenMS::ChromatogramSettings::getNativeID
const String & getNativeID() const
returns the native identifier for the spectrum, used by the acquisition software.
OpenMS::Size
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
OpenMS::PeakIntegrator::PeakShapeMetrics::width_at_50
double width_at_50
Definition: PeakIntegrator.h:141
OpenMS::PeakIntegrator::PeakArea
Definition: PeakIntegrator.h:85
OpenMS::MRMTransitionGroup
The representation of a group of transitions in a targeted proteomics experiment.
Definition: MRMTransitionGroup.h:67
OpenMS::PeakIntegrator::PeakBackground
Definition: PeakIntegrator.h:110
OpenMS::PeakIntegrator::PeakShapeMetrics
Definition: PeakIntegrator.h:128
OpenMS::Exception::IllegalArgument
A method or algorithm argument contains illegal values.
Definition: Exception.h:648
OpenMS::MRMTransitionGroupPicker::background_subtraction_
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1124
OpenMS::Constants::k
const double k
OpenMS::MRMFeature::addFeature
void addFeature(const Feature &feature, const String &key)
Adds an feature from a single chromatogram into the feature.
LogStream.h
OpenMS::MRMTransitionGroup::getTransitions
const std::vector< TransitionType > & getTransitions() const
Definition: MRMTransitionGroup.h:142
OpenMS::LinearResamplerAlign::raster
void raster(SpecT &container)
Applies the resampling algorithm to a container (MSSpectrum or MSChromatogram).
Definition: LinearResamplerAlign.h:80
OpenMS::MRMTransitionGroupPicker::compute_peak_quality_
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1128
File.h
SimpleOpenMSSpectraAccessFactory.h
OpenMS::MRMTransitionGroupPicker::remove_overlapping_features
void remove_overlapping_features(std::vector< SpectrumT > &picked_chroms, double best_left, double best_right)
Remove overlapping features.
Definition: MRMTransitionGroupPicker.h:716
OpenMS::Feature::setOverallQuality
void setOverallQuality(QualityType q)
Set the overall quality.
OpenMS::MRMTransitionGroup::getTransitionGroupID
const String & getTransitionGroupID() const
Definition: MRMTransitionGroup.h:130
OpenMS::Feature::getConvexHulls
const std::vector< ConvexHull2D > & getConvexHulls() const
Non-mutable access to the convex hulls.
OpenMS::PeakIntegrator
Compute the area, background and shape metrics of a peak.
Definition: PeakIntegrator.h:72
OpenMS::MRMTransitionGroup::getChromatogram
ChromatogramType & getChromatogram(const String &key)
Definition: MRMTransitionGroup.h:197
OpenMS::Peak2D::setMZ
void setMZ(CoordinateType coordinate)
Mutable access to the m/z coordinate (index 1)
Definition: Peak2D.h:202
OpenMS::Feature::setSubordinates
void setSubordinates(const std::vector< Feature > &rhs)
mutable access to subordinate features
OpenMS::MRMTransitionGroupPicker::resample_boundary_
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1137
OpenMS::Peak2D::setIntensity
void setIntensity(IntensityType intensity)
Non-mutable access to the data point intensity (height)
Definition: Peak2D.h:172
MRMFeature.h
OpenMS::MRMTransitionGroupPicker::min_peak_width_
double min_peak_width_
Definition: MRMTransitionGroupPicker.h:1135
Exception.h
OpenMS::FeatureXMLFile::store
void store(const String &filename, const FeatureMap &feature_map)
stores the map feature_map in file with name filename.
OpenMS::MetaInfoInterface::getMetaValue
const DataValue & getMetaValue(const String &name, const DataValue &default_value=DataValue::EMPTY) const
Returns the value corresponding to a string, or a default value (default: DataValue::EMPTY) if not fo...
OpenMS::MRMTransitionGroup::getChromatograms
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:174
OpenMS::TraMLFile::load
void load(const String &filename, TargetedExperiment &id)
Loads a map from a TraML file.
LinearResamplerAlign.h
OpenMS::PeakIntegrator::PeakArea::apex_pos
double apex_pos
Definition: PeakIntegrator.h:98
OpenMS::LinearResamplerAlign
Linear Resampling of raw data with alignment.
Definition: LinearResamplerAlign.h:57
MRMTransitionGroupPicker.h
OpenMS::FeatureMap
A container for features.
Definition: FeatureMap.h:95
OpenMS::MRMTransitionGroup::getPrecursorChromatograms
std::vector< ChromatogramType > & getPrecursorChromatograms()
Definition: MRMTransitionGroup.h:215
OpenMS::PeakPickerMRM::IDX_ABUNDANCE
Definition: PeakPickerMRM.h:83
OpenMS::MRMTransitionGroupPicker::use_precursors_
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1126
MSExperiment.h
double
OpenSwath::ChromatogramPtr
boost::shared_ptr< Chromatogram > ChromatogramPtr
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/DataStructures.h:172
FeatureMap.h
OpenMS::PeakPickerMRM
The PeakPickerMRM finds peaks a single chromatogram.
Definition: PeakPickerMRM.h:68
OpenMS::TOPPBase
Base class for TOPP applications.
Definition: TOPPBase.h:144
ChromatogramPeak.h
MzMLFile.h
int
OpenMS::PeakIntegrator::PeakShapeMetrics::end_position_at_10
double end_position_at_10
Definition: PeakIntegrator.h:161
OpenMS::Feature
An LC-MS feature.
Definition: Feature.h:70
OpenMS::TraMLFile
File adapter for HUPO PSI TraML files.
Definition: TraMLFile.h:63
OpenMS::MRMTransitionGroup::getTransition
const TransitionType & getTransition(String key)
Definition: MRMTransitionGroup.h:163
Scoring.h
OpenMS::MRMTransitionGroup::getPrecursorChromatogram
ChromatogramType & getPrecursorChromatogram(const String &key)
Definition: MRMTransitionGroup.h:248
PeakIntegrator.h
DataAccessHelper.h
OpenMS::MRMTransitionGroupPicker::min_qual_
double min_qual_
Definition: MRMTransitionGroupPicker.h:1131
OpenMS::MRMTransitionGroupPicker::pickTransitionGroup
void pickTransitionGroup(MRMTransitionGroup< SpectrumT, TransitionT > &transition_group)
Pick a group of chromatograms belonging to the same peptide.
Definition: MRMTransitionGroupPicker.h:113
OpenMS::MetaInfoInterface::setMetaValue
void setMetaValue(const String &name, const DataValue &value)
Sets the DataValue corresponding to a name.
OpenMS::PeakIntegrator::PeakArea::area
double area
Definition: PeakIntegrator.h:90
OpenMS::PeakIntegrator::PeakShapeMetrics::start_position_at_5
double start_position_at_5
Definition: PeakIntegrator.h:145
OpenMS::MRMTransitionGroupPicker::pickFragmentChromatograms
void pickFragmentChromatograms(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, MRMFeature &mrmFeature, const std::vector< SpectrumT > &smoothed_chroms, const double best_left, const double best_right, const bool use_consensus_, double &total_intensity, double &total_xic, double &total_mi, double &total_peak_apices, const SpectrumT &master_peak_container, const std::vector< double > &left_edges, const std::vector< double > &right_edges, const int chr_idx, const int peak_idx)
Definition: MRMTransitionGroupPicker.h:383
OpenMS::UniqueIdInterface::ensureUniqueId
Size ensureUniqueId()
Assigns a valid unique id, but only if the present one is invalid. Returns 1 if the unique id was cha...
Definition: UniqueIdInterface.h:154
OPENMS_PRECONDITION
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:136
OpenMS::MRMTransitionGroupPicker::stop_after_feature_
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1133
OpenMS::FeatureMap::setPrimaryMSRunPath
void setPrimaryMSRunPath(const StringList &s)
set the file path to the primary MS run (usually the mzML file obtained after data conversion from ra...
OpenMS::ProgressLogger::setLogType
void setLogType(LogType type) const
Sets the progress log that should be used. The default type is NONE!
OpenMS::Feature::setQuality
void setQuality(Size index, QualityType q)
Set the quality in dimension c.
OpenMS::PeakIntegrator::PeakShapeMetrics::start_position_at_10
double start_position_at_10
Definition: PeakIntegrator.h:149
MRMTransitionGroup.h
FeatureXMLFile.h
OpenMS::MSChromatogram
The representation of a chromatogram.
Definition: MSChromatogram.h:54
OpenMS::MRMTransitionGroupPicker::picker_
PeakPickerMRM picker_
Definition: MRMTransitionGroupPicker.h:1146
OpenMS::MSChromatogram::sortByIntensity
void sortByIntensity(bool reverse=false)
Lexicographically sorts the peaks by their intensity.
OpenMS::PeakIntegrator::PeakShapeMetrics::baseline_delta_2_height
double baseline_delta_2_height
Definition: PeakIntegrator.h:203
TraMLFile.h
OpenMS::Peak2D::setRT
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:214
OpenMS::MRMTransitionGroupPicker::stop_after_intensity_ratio_
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1134
OpenMS::DefaultParamHandler::setParameters
void setParameters(const Param &param)
Sets the parameters.
OpenMS::MRMFeature
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:50
OPENMS_LOG_DEBUG
#define OPENMS_LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:470
OpenMS::PeakIntegrator::PeakShapeMetrics::slope_of_baseline
double slope_of_baseline
Definition: PeakIntegrator.h:198
OpenMS::MRMTransitionGroupPicker::prepareMasterContainer_
void prepareMasterContainer_(const SpectrumT &ref_chromatogram, SpectrumT &master_peak_container, double left_boundary, double right_boundary)
Create an empty master peak container that has the correct mz / RT values set.
Definition: MRMTransitionGroupPicker.h:1063
OpenMS::MRMTransitionGroupPicker::computeQuality_
double computeQuality_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, const int chr_idx, const double best_left, const double best_right, String &outlier)
Compute transition group quality (higher score is better)
Definition: MRMTransitionGroupPicker.h:806
OpenMS::MSExperiment
In-Memory representation of a mass spectrometry experiment.
Definition: MSExperiment.h:77
OpenMS::MRMTransitionGroupPicker::resampleChromatogram_
SpectrumT resampleChromatogram_(const SpectrumT &chromatogram, const SpectrumT &master_peak_container, double left_boundary, double right_boundary)
Resample a container at the positions indicated by the master peak container.
Definition: MRMTransitionGroupPicker.h:1099
OpenMS::MRMFeature::getFeatures
const std::vector< Feature > & getFeatures() const
get a list of features
TOPPBase.h
OpenMS::Param::copy
Param copy(const String &prefix, bool remove_prefix=false) const
Returns a new Param object containing all entries that start with prefix.
OpenMS::MRMTransitionGroupPicker::recalculate_peaks_
bool recalculate_peaks_
Definition: MRMTransitionGroupPicker.h:1125
OpenSwath::mean_and_stddev
functor to compute the mean and stddev of sequence using the std::foreach algorithm
Definition: StatsHelpers.h:169
OpenSwath::Scoring::xcorrArrayGetMaxPeak
OPENSWATHALGO_DLLAPI XCorrArrayType::const_iterator xcorrArrayGetMaxPeak(const XCorrArrayType &array)
Find best peak in an cross-correlation (highest apex)
OpenMS::PeakIntegrator::PeakShapeMetrics::end_position_at_5
double end_position_at_5
Definition: PeakIntegrator.h:157
OpenMS::DefaultParamHandler::getDefaults
const Param & getDefaults() const
Non-mutable access to the default parameters.
OpenMS::PeakIntegrator::PeakShapeMetrics::end_position_at_50
double end_position_at_50
Definition: PeakIntegrator.h:165
MSChromatogram.h
ProgressLogger.h
OpenMS::Peak2D::getIntensity
IntensityType getIntensity() const
Definition: Peak2D.h:166
DefaultParamHandler.h
OpenMS::DefaultParamHandler
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:91
OpenMS::MRMTransitionGroupPicker::compute_total_mi_
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1130
OpenMS::MzMLFile::load
void load(const String &filename, PeakMap &map)
Loads a map from a MzML file. Spectra and chromatograms are sorted by default (this can be disabled u...
OpenMS::MRMTransitionGroupPicker::peak_integration_
String peak_integration_
Definition: MRMTransitionGroupPicker.h:1123
OpenMS::MRMTransitionGroupPicker::recalculatePeakBorders_
void recalculatePeakBorders_(const std::vector< SpectrumT > &picked_chroms, double &best_left, double &best_right, double max_z)
Recalculate the borders of the peak.
Definition: MRMTransitionGroupPicker.h:957
OpenMS::PeakIntegrator::PeakShapeMetrics::width_at_5
double width_at_5
Definition: PeakIntegrator.h:133
OpenSwath::SpectrumAccessPtr
boost::shared_ptr< ISpectrumAccess > SpectrumAccessPtr
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:89
ISpectrumAccess.h
OpenMS::MRMTransitionGroupPicker::pickPrecursorChromatograms
void pickPrecursorChromatograms(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, MRMFeature &mrmFeature, const std::vector< SpectrumT > &smoothed_chroms, const double best_left, const double best_right, const bool use_consensus_, double &total_intensity, const SpectrumT &master_peak_container, const std::vector< double > &left_edges, const std::vector< double > &right_edges, const int chr_idx, const int peak_idx)
Definition: MRMTransitionGroupPicker.h:584
OpenMS::ChromatogramSettings::setNativeID
void setNativeID(const String &native_id)
sets the native identifier for the spectrum, used by the acquisition software.
OpenMS::MzMLFile
File adapter for MzML files.
Definition: MzMLFile.h:55
OpenMS::MRMTransitionGroupPicker::pi_
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1147
OpenMS::MRMTransitionGroupPicker::boundary_selection_method_
String boundary_selection_method_
Which method to use for selecting peaks' boundaries.
Definition: MRMTransitionGroupPicker.h:1144
OpenMS::PeakIntegrator::PeakShapeMetrics::start_position_at_50
double start_position_at_50
Definition: PeakIntegrator.h:153
OpenMS::PeakPickerMRM::IDX_LEFTBORDER
Definition: PeakPickerMRM.h:83
OpenMS::PeakPickerMRM::IDX_RIGHTBORDER
Definition: PeakPickerMRM.h:83
OpenMS::MRMTransitionGroupPicker::compute_peak_shape_metrics_
bool compute_peak_shape_metrics_
Definition: MRMTransitionGroupPicker.h:1129
OpenMS::String
A more convenient string class.
Definition: String.h:58
OpenMS::MRMTransitionGroupPicker::recalculate_peaks_max_z_
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1136
OpenMS::MRMTransitionGroupPicker
The MRMTransitionGroupPicker finds peaks in chromatograms that belong to the same precursors.
Definition: MRMTransitionGroupPicker.h:79
OpenMS::MRMFeature::addPrecursorFeature
void addPrecursorFeature(const Feature &feature, const String &key)
Adds a precursor feature from a single chromatogram into the feature.
OpenMS::PeakIntegrator::PeakArea::height
double height
Definition: PeakIntegrator.h:94
OpenMS::Math::mean
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:133
MSSpectrum.h
OpenSwath::Scoring::rankedMutualInformation
OPENSWATHALGO_DLLAPI double rankedMutualInformation(std::vector< double > &data1, std::vector< double > &data2)
OpenSwath::Scoring::XCorrArrayType
Definition: Scoring.h:61
OpenMS::ConvexHull2D::setHullPoints
void setHullPoints(const PointArrayType &points)
accessor for the outer(!) points (no checking is performed if this is actually a convex hull)
OpenMS::MRMTransitionGroup::hasChromatogram
bool hasChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:192
OpenMS::ReactionMonitoringTransition
This class stores a SRM/MRM transition.
Definition: ReactionMonitoringTransition.h:56
OpenMS
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
OpenMS::ConvexHull2D
A 2-dimensional hull representation in [counter]clockwise direction - depending on axis labelling.
Definition: ConvexHull2D.h:72
OpenMS::PeakIntegrator::PeakShapeMetrics::asymmetry_factor
double asymmetry_factor
Definition: PeakIntegrator.h:193
OpenMS::UInt
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
PeakPickerMRM.h
OpenMS::MRMTransitionGroupPicker::selectChromHelper_
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:773
OpenMS::MRMTransitionGroupPicker::use_consensus_
bool use_consensus_
Definition: MRMTransitionGroupPicker.h:1127
OpenMS::PeakIntegrator::PeakShapeMetrics::width_at_10
double width_at_10
Definition: PeakIntegrator.h:137
OpenMS::PeakIntegrator::PeakShapeMetrics::points_across_baseline
Int points_across_baseline
Definition: PeakIntegrator.h:207
OpenMS::PeakIntegrator::PeakArea::hull_points
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:102
OpenMS::MRMTransitionGroupPicker::pickApex
void pickApex(std::vector< SpectrumT > &picked_chroms, const double best_left, const double best_right, const double peak_apex, double &min_left, double &max_right, std::vector< double > &left_edges, std::vector< double > &right_edges)
Apex-based peak picking.
Definition: MRMTransitionGroupPicker.h:342
OpenMS::PeakIntegrator::PeakShapeMetrics::points_across_half_height
Int points_across_half_height
Definition: PeakIntegrator.h:211
main
int main(int argc, const char **argv)
Definition: INIFileEditor.cpp:73
OpenMS::TargetedExperiment
A description of a targeted experiment containing precursor and production ions.
Definition: TargetedExperiment.h:64
OpenMS::PeakIntegrator::PeakShapeMetrics::total_width
double total_width
Definition: PeakIntegrator.h:169
OpenMS::MRMTransitionGroupPicker::createMRMFeature
MRMFeature createMRMFeature(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, std::vector< SpectrumT > &picked_chroms, const std::vector< SpectrumT > &smoothed_chroms, const int chr_idx, const int peak_idx)
Create feature from a vector of chromatograms and a specified peak.
Definition: MRMTransitionGroupPicker.h:217
OpenSwath::mean_and_stddev::mean
double mean() const
Definition: StatsHelpers.h:207
OpenMS::PeakIntegrator::PeakShapeMetrics::tailing_factor
double tailing_factor
Definition: PeakIntegrator.h:183