OpenMS
MRMTransitionGroupPicker.h
Go to the documentation of this file.
1 // Copyright (c) 2002-present, The OpenMS Team -- EKU Tuebingen, ETH Zurich, and FU Berlin
2 // SPDX-License-Identifier: BSD-3-Clause
3 //
4 // --------------------------------------------------------------------------
5 // $Maintainer: Hannes Roest $
6 // $Authors: Hannes Roest $
7 // --------------------------------------------------------------------------
8 
9 #pragma once
10 
17 
20 
23 
24 // Cross-correlation
27 
28 #include <numeric>
29 
30 //#define DEBUG_TRANSITIONGROUPPICKER
31 
32 namespace OpenMS
33 {
34 
53  class OPENMS_DLLAPI MRMTransitionGroupPicker :
54  public DefaultParamHandler
55  {
56 
57 public:
58 
62 
66 
86  template <typename SpectrumT, typename TransitionT>
88  {
89  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
90  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
91 
92  std::vector<MSChromatogram > picked_chroms;
93  std::vector<MSChromatogram > smoothed_chroms;
94 
95  // Pick fragment ion chromatograms
96  for (Size k = 0; k < transition_group.getChromatograms().size(); k++)
97  {
98  MSChromatogram& chromatogram = transition_group.getChromatograms()[k];
99  String native_id = chromatogram.getNativeID();
100 
101  // only pick detecting transitions (skip all others)
102  if (transition_group.getTransitions().size() > 0 &&
103  transition_group.hasTransition(native_id) &&
104  !transition_group.getTransition(native_id).isDetectingTransition() )
105  {
106  continue;
107  }
108 
109  MSChromatogram picked_chrom, smoothed_chrom;
110  smoothed_chrom.setNativeID(native_id);
111  picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
112  picked_chrom.sortByIntensity();
113  picked_chroms.push_back(std::move(picked_chrom));
114  smoothed_chroms.push_back(std::move(smoothed_chrom));
115  }
116 
117  // Pick precursor chromatograms
118  if (use_precursors_)
119  {
120  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
121  {
122  SpectrumT picked_chrom, smoothed_chrom;
123  SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
124 
125  picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
126  picked_chrom.sortByIntensity();
127  picked_chroms.push_back(picked_chrom);
128  smoothed_chroms.push_back(smoothed_chrom);
129  }
130  }
131 
132  // Find features (peak groups) in this group of transitions.
133  // While there are still peaks left, one will be picked and used to create
134  // a feature. Whenever we run out of peaks, we will get -1 back as index
135  // and terminate.
136  int chr_idx, peak_idx, cnt = 0;
137  std::vector<MRMFeature> features;
138  while (true)
139  {
140  chr_idx = -1; peak_idx = -1;
141 
142  if (boundary_selection_method_ == "largest")
143  {
144  findLargestPeak(picked_chroms, chr_idx, peak_idx);
145  }
146  else if (boundary_selection_method_ == "widest")
147  {
148  findWidestPeakIndices(picked_chroms, chr_idx, peak_idx);
149  }
150 
151  if (chr_idx == -1 && peak_idx == -1)
152  {
153  OPENMS_LOG_DEBUG << "**** No more peaks left. Nr. chroms: " << picked_chroms.size() << std::endl;
154  break;
155  }
156 
157  // Compute a feature from the individual chromatograms and add non-zero features
158  MRMFeature mrm_feature = createMRMFeature(transition_group, picked_chroms, smoothed_chroms, chr_idx, peak_idx);
159  double total_xic = 0;
160  double intensity = mrm_feature.getIntensity();
161  if (intensity > 0)
162  {
163  total_xic = mrm_feature.getMetaValue("total_xic");
164  features.push_back(std::move(mrm_feature));
165  cnt++;
166  }
167 
168  if (stop_after_feature_ > 0 && cnt > stop_after_feature_) {
169  // If you set this, you only expect one feature anyway. No logging necessary why it stopped.
170  break;
171  }
172  if (intensity > 0 && intensity / total_xic < stop_after_intensity_ratio_)
173  {
174  OPENMS_LOG_DEBUG << "**** Minimum intensity ratio reached. Nr. chroms: " << picked_chroms.size() << std::endl;
175  break;
176  }
177  }
178 
179  // Check for completely overlapping features
180  for (Size i = 0; i < features.size(); i++)
181  {
182  MRMFeature& mrm_feature = features[i];
183  bool skip = false;
184  for (Size j = 0; j < i; j++)
185  {
186  if ((double)mrm_feature.getMetaValue("leftWidth") >= (double)features[j].getMetaValue("leftWidth") &&
187  (double)mrm_feature.getMetaValue("rightWidth") <= (double)features[j].getMetaValue("rightWidth") )
188  { skip = true; }
189  }
190  if (mrm_feature.getIntensity() > 0 && !skip)
191  {
192  transition_group.addFeature(mrm_feature);
193  }
194  }
195 
196  }
197 
199  template <typename SpectrumT, typename TransitionT>
201  std::vector<SpectrumT>& picked_chroms,
202  const std::vector<SpectrumT>& smoothed_chroms,
203  const int chr_idx,
204  const int peak_idx)
205  {
206  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
207  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
208 
209  MRMFeature mrmFeature;
210  mrmFeature.setIntensity(0.0);
211  double best_left = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][peak_idx];
212  double best_right = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][peak_idx];
213  double peak_apex = picked_chroms[chr_idx][peak_idx].getRT();
214  OPENMS_LOG_DEBUG << "**** Creating MRMFeature for peak " << peak_idx << " in chrom. " << chr_idx << " with " <<
215  picked_chroms[chr_idx][peak_idx] << " and borders " << best_left << " " <<
216  best_right << " (" << best_right - best_left << ")" << std::endl;
217 
218  if (use_consensus_ && recalculate_peaks_)
219  {
220  // This may change best_left / best_right
221  recalculatePeakBorders_(picked_chroms, best_left, best_right, recalculate_peaks_max_z_);
222  if (peak_apex < best_left || peak_apex > best_right)
223  {
224  // apex fell out of range, lets correct it
225  peak_apex = (best_left + best_right) / 2.0;
226  }
227  }
228 
229  std::vector< double > left_edges;
230  std::vector< double > right_edges;
231  double min_left = best_left;
232  double max_right = best_right;
233  if (use_consensus_)
234  {
235  // Remove other, overlapping, picked peaks (in this and other
236  // chromatograms) and then ensure that at least one peak is set to zero
237  // (the currently best peak).
238  remove_overlapping_features(picked_chroms, best_left, best_right);
239  }
240  else
241  {
242  pickApex(picked_chroms, best_left, best_right, peak_apex,
243  min_left, max_right, left_edges, right_edges);
244 
245  } // end !use_consensus_
246  picked_chroms[chr_idx][peak_idx].setIntensity(0.0); // ensure that we set at least one peak to zero
247 
248  // Check for minimal peak width -> return empty feature (Intensity zero)
249  if (use_consensus_)
250  {
251  if (min_peak_width_ > 0.0 && std::fabs(best_right - best_left) < min_peak_width_)
252  {
253  return mrmFeature;
254  }
255 
256  if (compute_peak_quality_)
257  {
258  String outlier = "none";
259  double qual = computeQuality_(transition_group, picked_chroms, chr_idx, best_left, best_right, outlier);
260  if (qual < min_qual_)
261  {
262  return mrmFeature;
263  }
264  mrmFeature.setMetaValue("potentialOutlier", outlier);
265  mrmFeature.setMetaValue("initialPeakQuality", qual);
266  mrmFeature.setOverallQuality(qual);
267  }
268  }
269 
270  // Prepare linear resampling of all the chromatograms, here creating the
271  // empty master_peak_container with the same RT (m/z) values as the
272  // reference chromatogram. We use the overall minimal left boundary and
273  // maximal right boundary to prepare the container.
274  SpectrumT master_peak_container;
275  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
276  prepareMasterContainer_(ref_chromatogram, master_peak_container, min_left, max_right);
277 
278  // Iterate over initial transitions / chromatograms (note that we may
279  // have a different number of picked chromatograms than total transitions
280  // as not all are detecting transitions).
281  double total_intensity = 0; double total_peak_apices = 0; double total_xic = 0; double total_mi = 0;
282  pickFragmentChromatograms(transition_group, picked_chroms, mrmFeature, smoothed_chroms,
283  best_left, best_right, use_consensus_,
284  total_intensity, total_xic, total_mi, total_peak_apices,
285  master_peak_container, left_edges, right_edges,
286  chr_idx, peak_idx);
287 
288  // Also pick the precursor chromatogram(s); note total_xic is not
289  // extracted here, only for fragment traces
290  pickPrecursorChromatograms(transition_group,
291  picked_chroms, mrmFeature, smoothed_chroms,
292  best_left, best_right, use_consensus_,
293  total_intensity, master_peak_container, left_edges, right_edges,
294  chr_idx, peak_idx);
295 
296  mrmFeature.setRT(peak_apex);
297  mrmFeature.setIntensity(total_intensity);
298  mrmFeature.setMetaValue("PeptideRef", transition_group.getTransitionGroupID());
299  mrmFeature.setMetaValue("leftWidth", best_left);
300  mrmFeature.setMetaValue("rightWidth", best_right);
301  mrmFeature.setMetaValue("total_xic", total_xic);
302  if (compute_total_mi_)
303  {
304  mrmFeature.setMetaValue("total_mi", total_mi);
305  }
306  mrmFeature.setMetaValue("peak_apices_sum", total_peak_apices);
307 
308  mrmFeature.ensureUniqueId();
309  return mrmFeature;
310  }
311 
324  template <typename SpectrumT>
325  void pickApex(std::vector<SpectrumT>& picked_chroms,
326  const double best_left, const double best_right, const double peak_apex,
327  double &min_left, double &max_right,
328  std::vector< double > & left_edges, std::vector< double > & right_edges)
329  {
330  for (Size k = 0; k < picked_chroms.size(); k++)
331  {
332  double peak_apex_dist_min = std::numeric_limits<double>::max();
333  int min_dist = -1;
334  for (Size i = 0; i < picked_chroms[k].size(); i++)
335  {
336  PeakIntegrator::PeakArea pa_tmp = pi_.integratePeak( // get the peak apex
337  picked_chroms[k],
338  picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][i],
339  picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][i]);
340  if (pa_tmp.apex_pos > 1e-11 && std::fabs(pa_tmp.apex_pos - peak_apex) < peak_apex_dist_min)
341  { // update best candidate
342  peak_apex_dist_min = std::fabs(pa_tmp.apex_pos - peak_apex);
343  min_dist = (int)i;
344  }
345  }
346 
347  // Select master peak boundaries, or in the case we found at least one peak, the local peak boundaries
348  double l = best_left;
349  double r = best_right;
350  if (min_dist >= 0)
351  {
352  l = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][min_dist];
353  r = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][min_dist];
354  picked_chroms[k][min_dist].setIntensity(0.0); // only remove one peak per transition
355  }
356 
357  left_edges.push_back(l);
358  right_edges.push_back(r);
359  // ensure we remember the overall maxima / minima
360  if (l < min_left) {min_left = l;}
361  if (r > max_right) {max_right = r;}
362  }
363  }
364 
365  template <typename SpectrumT, typename TransitionT>
367  const std::vector<SpectrumT>& picked_chroms,
368  MRMFeature& mrmFeature,
369  const std::vector<SpectrumT>& smoothed_chroms,
370  const double best_left, const double best_right,
371  const bool use_consensus_,
372  double & total_intensity,
373  double & total_xic,
374  double & total_mi,
375  double & total_peak_apices,
376  const SpectrumT & master_peak_container,
377  const std::vector< double > & left_edges,
378  const std::vector< double > & right_edges,
379  const int chr_idx,
380  const int peak_idx)
381  {
382  for (Size k = 0; k < transition_group.getTransitions().size(); k++)
383  {
384 
385  double local_left = best_left;
386  double local_right = best_right;
387  if (!use_consensus_)
388  {
389  // We cannot have any non-detecting transitions (otherwise we have
390  // too few left / right edges) as we skipped those when doing peak
391  // picking and smoothing.
392  if (!transition_group.getTransitions()[k].isDetectingTransition())
393  {
394  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
395  "When using non-consensus peak picker, all transitions need to be detecting transitions.");
396  }
397  local_left = left_edges[k];
398  local_right = right_edges[k];
399  }
400 
401  const SpectrumT& chromatogram = selectChromHelper_(transition_group, transition_group.getTransitions()[k].getNativeID());
402  if (transition_group.getTransitions()[k].isDetectingTransition())
403  {
404  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
405  {
406  total_xic += it->getIntensity();
407  }
408  }
409 
410  // Compute total intensity on transition-level
411  double transition_total_xic = 0;
412 
413  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
414  {
415  transition_total_xic += it->getIntensity();
416  }
417 
418  // Compute total mutual information on transition-level.
419  double transition_total_mi = 0;
420  if (compute_total_mi_)
421  {
422  std::vector<unsigned int> chrom_vect_id_ranked, chrom_vect_det_ranked;
423  std::vector<double> chrom_vect_id, chrom_vect_det;
424  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
425  {
426  chrom_vect_id.push_back(it->getIntensity());
427  }
428  unsigned int max_rank_det = OpenSwath::Scoring::computeAndAppendRank(chrom_vect_id, chrom_vect_det_ranked);
429  // compute baseline mutual information
430  int transition_total_mi_norm = 0;
431  for (Size m = 0; m < transition_group.getTransitions().size(); m++)
432  {
433  if (transition_group.getTransitions()[m].isDetectingTransition())
434  {
435  const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.getTransitions()[m].getNativeID());
436  chrom_vect_det.clear();
437  for (typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
438  {
439  chrom_vect_det.push_back(it->getIntensity());
440  }
441  unsigned int max_rank_id = OpenSwath::Scoring::computeAndAppendRank(chrom_vect_det, chrom_vect_id_ranked);
442  transition_total_mi += OpenSwath::Scoring::rankedMutualInformation(chrom_vect_id_ranked, chrom_vect_det_ranked, max_rank_id, max_rank_det);
443  transition_total_mi_norm++;
444  }
445  }
446  if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
447 
448  if (transition_group.getTransitions()[k].isDetectingTransition())
449  {
450  // sum up all transition-level total MI and divide by the number of detection transitions to have peak group level total MI
451  total_mi += transition_total_mi / transition_total_mi_norm;
452  }
453  }
454 
455  SpectrumT used_chromatogram;
456  // resample the current chromatogram
457  if (peak_integration_ == "original")
458  {
459  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
460  }
461  else if (peak_integration_ == "smoothed")
462  {
463  if (smoothed_chroms.size() <= k)
464  {
465  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
466  "Tried to calculate peak area and height without any smoothed chromatograms");
467  }
468  used_chromatogram = resampleChromatogram_(smoothed_chroms[k], master_peak_container, local_left, local_right);
469  }
470  else
471  {
472  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
473  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
474  }
475 
476  Feature f;
477  double quality = 0;
478  f.setQuality(0, quality);
479  f.setOverallQuality(quality);
480 
481  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
482  double peak_integral = pa.area;
483  double peak_apex_int = pa.height;
484  f.setMetaValue("peak_apex_position", pa.apex_pos);
485  if (background_subtraction_ != "none")
486  {
487  double background{0};
488  double avg_noise_level{0};
489  if (background_subtraction_ == "original")
490  {
491  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
492  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
493  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
494  avg_noise_level = (intensity_right + intensity_left) / 2;
495  background = avg_noise_level * n_points;
496  }
497  else if (background_subtraction_ == "exact")
498  {
499  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
500  background = pb.area;
501  avg_noise_level = pb.height;
502  }
503  peak_integral -= background;
504  peak_apex_int -= avg_noise_level;
505  if (peak_integral < 0) {peak_integral = 0;}
506  if (peak_apex_int < 0) {peak_apex_int = 0;}
507 
508  f.setMetaValue("area_background_level", background);
509  f.setMetaValue("noise_background_level", avg_noise_level);
510  } // end background
511 
512  f.setRT(picked_chroms[chr_idx][peak_idx].getPos());
513  f.setIntensity(peak_integral);
514  ConvexHull2D hull;
515  hull.setHullPoints(pa.hull_points);
516  f.getConvexHulls().push_back(hull);
517 
518  f.setMZ(chromatogram.getProduct().getMZ());
519  mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());
520 
521  if (chromatogram.metaValueExists("product_mz")) // legacy code (ensures that old tests still work)
522  {
523  f.setMetaValue("MZ", chromatogram.getMetaValue("product_mz"));
524  f.setMZ(chromatogram.getMetaValue("product_mz"));
525  }
526 
527  f.setMetaValue("native_id", chromatogram.getNativeID());
528  f.setMetaValue("peak_apex_int", peak_apex_int);
529  f.setMetaValue("total_xic", transition_total_xic);
530  if (compute_total_mi_)
531  {
532  f.setMetaValue("total_mi", transition_total_mi);
533  }
534 
535  if (transition_group.getTransitions()[k].isQuantifyingTransition())
536  {
537  total_intensity += peak_integral;
538  total_peak_apices += peak_apex_int;
539  }
540 
541  // for backwards compatibility with TOPP tests
542  // Calculate peak shape metrics that will be used for later QC
543  PeakIntegrator::PeakShapeMetrics psm = pi_.calculatePeakShapeMetrics(used_chromatogram, local_left, local_right, peak_apex_int, pa.apex_pos);
544  f.setMetaValue("width_at_50", psm.width_at_50);
545  if (compute_peak_shape_metrics_)
546  {
547  f.setMetaValue("width_at_5", psm.width_at_5);
548  f.setMetaValue("width_at_10", psm.width_at_10);
549  f.setMetaValue("start_position_at_5", psm.start_position_at_5);
550  f.setMetaValue("start_position_at_10", psm.start_position_at_10);
551  f.setMetaValue("start_position_at_50", psm.start_position_at_50);
552  f.setMetaValue("end_position_at_5", psm.end_position_at_5);
553  f.setMetaValue("end_position_at_10", psm.end_position_at_10);
554  f.setMetaValue("end_position_at_50", psm.end_position_at_50);
555  f.setMetaValue("total_width", psm.total_width);
556  f.setMetaValue("tailing_factor", psm.tailing_factor);
557  f.setMetaValue("asymmetry_factor", psm.asymmetry_factor);
558  f.setMetaValue("slope_of_baseline", psm.slope_of_baseline);
559  f.setMetaValue("baseline_delta_2_height", psm.baseline_delta_2_height);
560  f.setMetaValue("points_across_baseline", psm.points_across_baseline);
561  f.setMetaValue("points_across_half_height", psm.points_across_half_height);
562  }
563 
564  mrmFeature.addFeature(f, chromatogram.getNativeID()); //map index and feature
565  }
566  }
567 
568  template <typename SpectrumT, typename TransitionT>
570  const std::vector<SpectrumT>& picked_chroms,
571  MRMFeature& mrmFeature,
572  const std::vector<SpectrumT>& smoothed_chroms,
573  const double best_left, const double best_right,
574  const bool use_consensus_,
575  double & total_intensity,
576  const SpectrumT & master_peak_container,
577  const std::vector< double > & left_edges,
578  const std::vector< double > & right_edges,
579  const int chr_idx,
580  const int peak_idx)
581  {
582  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
583  {
584  const SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
585 
586  // Identify precursor index
587  // note: this is only valid if all transitions are detecting transitions
588  Size prec_idx = transition_group.getChromatograms().size() + k;
589 
590  double local_left = best_left;
591  double local_right = best_right;
592  if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
593  {
594  local_left = left_edges[prec_idx];
595  local_right = right_edges[prec_idx];
596  }
597 
598  SpectrumT used_chromatogram;
599  // resample the current chromatogram
600  if (peak_integration_ == "original")
601  {
602  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
603  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
604  }
605  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= prec_idx)
606  {
607  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
608  "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
609  }
610  else if (peak_integration_ == "smoothed")
611  {
612  used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
613  }
614  else
615  {
616  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
617  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
618  }
619 
620  Feature f;
621  double quality = 0;
622  f.setQuality(0, quality);
623  f.setOverallQuality(quality);
624 
625  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
626  double peak_integral = pa.area;
627  double peak_apex_int = pa.height;
628 
629  if (background_subtraction_ != "none")
630  {
631  double background{0};
632  double avg_noise_level{0};
633  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= prec_idx)
634  {
635  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
636  "Tried to calculate background estimation without any smoothed chromatograms");
637  }
638  else if (background_subtraction_ == "original")
639  {
640  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
641  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
642  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
643  avg_noise_level = (intensity_right + intensity_left) / 2;
644  background = avg_noise_level * n_points;
645  }
646  else if (background_subtraction_ == "exact")
647  {
648  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
649  background = pb.area;
650  avg_noise_level = pb.height;
651  }
652  peak_integral -= background;
653  peak_apex_int -= avg_noise_level;
654  if (peak_integral < 0) {peak_integral = 0;}
655  if (peak_apex_int < 0) {peak_apex_int = 0;}
656 
657  f.setMetaValue("area_background_level", background);
658  f.setMetaValue("noise_background_level", avg_noise_level);
659  }
660 
661  f.setMZ(chromatogram.getPrecursor().getMZ());
662  if (k == 0) {mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());} // only use m/z if first (monoisotopic) isotope
663 
664  if (chromatogram.metaValueExists("precursor_mz")) // legacy code (ensures that old tests still work)
665  {
666  f.setMZ(chromatogram.getMetaValue("precursor_mz"));
667  if (k == 0) {mrmFeature.setMZ(chromatogram.getMetaValue("precursor_mz"));} // only use m/z if first (monoisotopic) isotope
668  }
669 
670  f.setRT(picked_chroms[chr_idx][peak_idx].getPos());
671  f.setIntensity(peak_integral);
672  ConvexHull2D hull;
673  hull.setHullPoints(pa.hull_points);
674  f.getConvexHulls().push_back(hull);
675  f.setMetaValue("native_id", chromatogram.getNativeID());
676  f.setMetaValue("peak_apex_int", peak_apex_int);
677 
678  if (use_precursors_ && transition_group.getTransitions().empty())
679  {
680  total_intensity += peak_integral;
681  }
682 
683  mrmFeature.addPrecursorFeature(f, chromatogram.getNativeID());
684  }
685  }
686 
687  // maybe private, but we have tests
688 
700  template <typename SpectrumT>
701  void remove_overlapping_features(std::vector<SpectrumT>& picked_chroms, double best_left, double best_right)
702  {
703  // delete all seeds that lie within the current seed
704  Size count_inside = 0;
705  for (Size k = 0; k < picked_chroms.size(); k++)
706  {
707  for (Size i = 0; i < picked_chroms[k].size(); i++)
708  {
709  if (picked_chroms[k][i].getPos() >= best_left && picked_chroms[k][i].getPos() <= best_right)
710  {
711  picked_chroms[k][i].setIntensity(0.0);
712  count_inside++;
713  }
714  }
715  }
716 
717  Size count_overlap = 0;
718  // delete all seeds that overlap 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].getIntensity() <= 1e-11) {continue; }
724 
725  double left = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][i];
726  double right = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][i];
727  if ((left > best_left && left < best_right)
728  || (right > best_left && right < best_right))
729  {
730  picked_chroms[k][i].setIntensity(0.0);
731  count_overlap++;
732  }
733  }
734  }
735  OPENMS_LOG_DEBUG << " ** Removed " << count_inside << " peaks enclosed in and " <<
736  count_overlap << " peaks overlapping with added feature." << std::endl;
737  }
738 
740  void findLargestPeak(const std::vector<MSChromatogram >& picked_chroms, int& chr_idx, int& peak_idx);
741 
750  void findWidestPeakIndices(const std::vector<MSChromatogram>& picked_chroms, Int& chrom_idx, Int& point_idx) const;
751 
752 protected:
753 
755  void updateMembers_() override;
756 
759 
763  template <typename SpectrumT, typename TransitionT>
764  const SpectrumT& selectChromHelper_(const MRMTransitionGroup<SpectrumT, TransitionT>& transition_group, const String& native_id)
765  {
766  if (transition_group.hasChromatogram(native_id))
767  {
768  return transition_group.getChromatogram(native_id);
769  }
770  else if (transition_group.hasPrecursorChromatogram(native_id))
771  {
772  return transition_group.getPrecursorChromatogram(native_id);
773  }
774  else
775  {
776  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION, "Did not find chromatogram for id '" + native_id + "'.");
777  }
778  }
779 
796  template <typename SpectrumT, typename TransitionT>
798  const std::vector<SpectrumT>& picked_chroms,
799  const int chr_idx,
800  const double best_left,
801  const double best_right,
802  String& outlier)
803  {
804  // Resample all chromatograms around the current estimated peak and
805  // collect the raw intensities. For resampling, use a bit more on either
806  // side to correctly identify shoulders etc.
807  double resample_boundary = resample_boundary_; // sample 15 seconds more on each side
808  SpectrumT master_peak_container;
809  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
810  prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
811  std::vector<std::vector<double> > all_ints;
812  for (Size k = 0; k < picked_chroms.size(); k++)
813  {
814  const SpectrumT& chromatogram = selectChromHelper_(transition_group, picked_chroms[k].getNativeID());
815  const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
816  master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
817 
818  std::vector<double> int_here;
819  for (const auto& peak : used_chromatogram) int_here.push_back(peak.getIntensity());
820  // Remove chromatograms without a single peak
821  double tic = std::accumulate(int_here.begin(), int_here.end(), 0.0);
822  if (tic > 1e-11) all_ints.push_back(int_here);
823  }
824 
825  // Compute the cross-correlation for the collected intensities
826  std::vector<double> mean_shapes;
827  std::vector<double> mean_coel;
828  for (Size k = 0; k < all_ints.size(); k++)
829  {
830  std::vector<double> shapes;
831  std::vector<double> coel;
832  for (Size i = 0; i < all_ints.size(); i++)
833  {
834  if (i == k) {continue;}
836  all_ints[k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
837 
838  // the first value is the x-axis (retention time) and should be an int -> it show the lag between the two
839  double res_coelution = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->first);
840  double res_shape = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->second);
841 
842  shapes.push_back(res_shape);
843  coel.push_back(res_coelution);
844  }
845 
846  // We have computed the cross-correlation of chromatogram k against
847  // all others. Use the mean of these computations as the value for k.
849  msc = std::for_each(shapes.begin(), shapes.end(), msc);
850  double shapes_mean = msc.mean();
851  msc = std::for_each(coel.begin(), coel.end(), msc);
852  double coel_mean = msc.mean();
853 
854  // mean shape scores below 0.5-0.6 should be a real sign of trouble ... !
855  // mean coel scores above 3.5 should be a real sign of trouble ... !
856  mean_shapes.push_back(shapes_mean);
857  mean_coel.push_back(coel_mean);
858  }
859 
860  // find the chromatogram with the minimal shape score and the maximal
861  // coelution score -> if it is the same chromatogram, the chance is
862  // pretty good that it is different from the others...
863  int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
864  int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
865 
866  // Look at the picked peaks that are within the current left/right borders
867  int missing_peaks = 0;
868  int multiple_peaks = 0;
869 
870  // collect all seeds that lie within the current seed
871  std::vector<double> left_borders;
872  std::vector<double> right_borders;
873  for (Size k = 0; k < picked_chroms.size(); k++)
874  {
875  double l_tmp;
876  double r_tmp;
877  double max_int = -1;
878 
879  int pfound = 0;
880  l_tmp = -1;
881  r_tmp = -1;
882  for (Size i = 0; i < picked_chroms[k].size(); i++)
883  {
884  if (picked_chroms[k][i].getPos() >= best_left && picked_chroms[k][i].getPos() <= best_right)
885  {
886  pfound++;
887  if (picked_chroms[k][i].getIntensity() > max_int)
888  {
889  max_int = picked_chroms[k][i].getIntensity() > max_int;
890  l_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][i];
891  r_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][i];
892  }
893  }
894  }
895 
896  if (l_tmp > 1e-11) left_borders.push_back(l_tmp);
897  if (r_tmp > 1e-11) right_borders.push_back(r_tmp);
898 
899  if (pfound == 0) missing_peaks++;
900  if (pfound > 1) multiple_peaks++;
901  }
902 
903  // Check how many chromatograms had exactly one peak picked between our
904  // current left/right borders -> this would be a sign of consistency.
905  OPENMS_LOG_DEBUG << " Overall found missing : " << missing_peaks << " and multiple : " << multiple_peaks << std::endl;
906 
908 
909  // Is there one transitions that is very different from the rest (e.g.
910  // the same element has a bad shape and a bad coelution score) -> potential outlier
911  if (min_index_shape == max_index_coel)
912  {
913  OPENMS_LOG_DEBUG << " Element " << min_index_shape << " is a candidate for removal ... " << std::endl;
914  outlier = String(picked_chroms[min_index_shape].getNativeID());
915  }
916  else
917  {
918  outlier = "none";
919  }
920 
921  // For the final score (larger is better), consider these scores:
922  // - missing_peaks (the more peaks are missing, the worse)
923  // - multiple_peaks
924  // - mean of the shapes (1 is very good, 0 is bad)
925  // - mean of the co-elution scores (0 is good, 1 is ok, above 1 is pretty bad)
926  double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
927  double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
928  coel_score = (coel_score - 1.0) / 2.0;
929 
930  double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
931 
932  OPENMS_LOG_DEBUG << " Computed score " << score << " (from " << shape_score <<
933  " - " << coel_score << " - " << 1.0 * missing_peaks / picked_chroms.size() << ")" << std::endl;
934 
935  return score;
936  }
937 
947  template <typename SpectrumT>
948  void recalculatePeakBorders_(const std::vector<SpectrumT>& picked_chroms, double& best_left, double& best_right, double max_z)
949  {
950  // 1. Collect all seeds that lie within the current seed
951  // - Per chromatogram only the most intense one counts, otherwise very
952  // - low intense peaks can contribute disproportionally to the voting
953  // - procedure.
954  // TODO Especially with DDA MS1 "transitions" from FFID, you often have exactly the same
955  // borders and the logs get confusing and you get -Nan CVs. You also might save time if you use a set here.
956  std::vector<double> left_borders;
957  std::vector<double> right_borders;
958  for (Size k = 0; k < picked_chroms.size(); k++)
959  {
960  double max_int = -1;
961  double left = -1;
962  double right = -1;
963  for (Size i = 0; i < picked_chroms[k].size(); i++)
964  {
965  if (picked_chroms[k][i].getPos() >= best_left && picked_chroms[k][i].getPos() <= best_right)
966  {
967  if (picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_ABUNDANCE][i] > max_int)
968  {
969  max_int = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_ABUNDANCE][i];
970  left = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_LEFTBORDER][i];
971  right = picked_chroms[k].getFloatDataArrays()[PeakPickerChromatogram::IDX_RIGHTBORDER][i];
972  }
973  }
974  }
975  if (max_int > -1 )
976  {
977  left_borders.push_back(left);
978  right_borders.push_back(right);
979  OPENMS_LOG_DEBUG << " * peak " << k << " left boundary " << left_borders.back() << " with inty " << max_int << std::endl;
980  OPENMS_LOG_DEBUG << " * peak " << k << " right boundary " << right_borders.back() << " with inty " << max_int << std::endl;
981  }
982  }
983 
984  // Return for empty peak list
985  if (right_borders.empty())
986  {
987  return;
988  }
989 
990  // FEATURE IDEA: instead of Z-score use modified Z-score for small data sets
991  // http://d-scholarship.pitt.edu/7948/1/Seo.pdf
992  // http://www.itl.nist.gov/div898/handbook/eda/section3/eda35h.htm
993  // 1. calculate median
994  // 2. MAD = calculate difference to median for each value -> take median of that
995  // 3. Mi = 0.6745*(xi - median) / MAD
996 
997  // 2. Calculate mean and standard deviation
998  // If the coefficient of variation is too large for one border, we use a
999  // "pseudo-median" instead of the border of the most intense peak.
1000  double mean, stdev;
1001 
1002  // Right borders
1003  mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
1004  stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
1005  / right_borders.size() - mean * mean);
1006  std::sort(right_borders.begin(), right_borders.end());
1007 
1008  OPENMS_LOG_DEBUG << " - Recalculating right peak boundaries " << mean << " mean / best "
1009  << best_right << " std " << stdev << " : " << std::fabs(best_right - mean) / stdev
1010  << " coefficient of variation" << std::endl;
1011 
1012  // Compare right borders of best transition with the mean
1013  if (std::fabs(best_right - mean) / stdev > max_z)
1014  {
1015  best_right = right_borders[right_borders.size() / 2]; // pseudo median
1016  OPENMS_LOG_DEBUG << " - CV too large: correcting right boundary to " << best_right << std::endl;
1017  }
1018 
1019  // Left borders
1020  mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
1021  stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
1022  / left_borders.size() - mean * mean);
1023  std::sort(left_borders.begin(), left_borders.end());
1024 
1025  OPENMS_LOG_DEBUG << " - Recalculating left peak boundaries " << mean << " mean / best "
1026  << best_left << " std " << stdev << " : " << std::fabs(best_left - mean) / stdev
1027  << " coefficient of variation" << std::endl;
1028 
1029  // Compare left borders of best transition with the mean
1030  if (std::fabs(best_left - mean) / stdev > max_z)
1031  {
1032  best_left = left_borders[left_borders.size() / 2]; // pseudo median
1033  OPENMS_LOG_DEBUG << " - CV too large: correcting left boundary to " << best_left << std::endl;
1034  }
1035 
1036  }
1037 
1039 
1040 
1055  template <typename PeakContainerT>
1056  void prepareMasterContainer_(const PeakContainerT& ref_chromatogram,
1057  PeakContainerT& master_peak_container,
1058  double left_boundary,
1059  double right_boundary)
1060  {
1061  OPENMS_PRECONDITION(master_peak_container.empty(), "Master peak container must be empty")
1062 
1063  // get the start / end point of this chromatogram => then add one more
1064  // point beyond the two boundaries to make the resampling accurate also
1065  // at the edge.
1066  auto begin = ref_chromatogram.begin();
1067  while (begin != ref_chromatogram.end() && begin->getPos() < left_boundary) {begin++; }
1068  if (begin != ref_chromatogram.begin()) { begin--; }
1069 
1070  auto end = begin;
1071  while (end != ref_chromatogram.end() && end->getPos() < right_boundary) {end++; }
1072  if (end != ref_chromatogram.end()) { end++; }
1073 
1074  // resize the master container and set the m/z values to the ones of the master container
1075  master_peak_container.resize(distance(begin, end)); // initialize to zero
1076  auto it = master_peak_container.begin();
1077  for (auto chrom_it = begin; chrom_it != end; chrom_it++, it++)
1078  {
1079  it->setPos(chrom_it->getPos());
1080  }
1081  }
1082 
1093  template <typename PeakContainerT>
1094  PeakContainerT resampleChromatogram_(const PeakContainerT& chromatogram,
1095  const PeakContainerT& master_peak_container,
1096  double left_boundary,
1097  double right_boundary)
1098  {
1099  // get the start / end point of this chromatogram => then add one more
1100  // point beyond the two boundaries to make the resampling accurate also
1101  // at the edge.
1102  auto begin = chromatogram.begin();
1103  while (begin != chromatogram.end() && begin->getPos() < left_boundary) {begin++;}
1104  if (begin != chromatogram.begin()) { begin--; }
1105 
1106  auto end = begin;
1107  while (end != chromatogram.end() && end->getPos() < right_boundary) {end++;}
1108  if (end != chromatogram.end()) { end++; }
1109 
1110  auto resampled_peak_container = master_peak_container; // copy the master container, which contains the RT values
1111  LinearResamplerAlign lresampler;
1112  lresampler.raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1113 
1114  return resampled_peak_container;
1115  }
1116 
1118 
1119  // Members
1128  double min_qual_;
1129 
1135 
1142 
1145  };
1146 }
1147 
1148 
#define OPENMS_LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:454
const String & getNativeID() const
returns the native identifier for the spectrum, used by the acquisition software.
void setNativeID(const String &native_id)
sets the native identifier for the spectrum, used by the acquisition software.
Definition: ConvexHull2D.h:47
void setHullPoints(const PointArrayType &points)
accessor for the outer(!) points (no checking is performed if this is actually a convex hull)
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:66
A method or algorithm argument contains illegal values.
Definition: Exception.h:629
An LC-MS feature.
Definition: Feature.h:46
const std::vector< ConvexHull2D > & getConvexHulls() const
Non-mutable access to the convex hulls.
void setQuality(Size index, QualityType q)
Set the quality in dimension c.
void setOverallQuality(QualityType q)
Set the overall quality.
Linear Resampling of raw data with alignment.
Definition: LinearResamplerAlign.h:33
void raster(PeakContainerT &container)
Applies the resampling algorithm to a container (MSSpectrum or MSChromatogram).
Definition: LinearResamplerAlign.h:54
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:26
void addFeature(const Feature &feature, const String &key)
Adds an feature from a single chromatogram into the feature.
void addPrecursorFeature(const Feature &feature, const String &key)
Adds a precursor feature from a single chromatogram into the feature.
The MRMTransitionGroupPicker finds peaks in chromatograms that belong to the same precursors.
Definition: MRMTransitionGroupPicker.h:55
void findLargestPeak(const std::vector< MSChromatogram > &picked_chroms, int &chr_idx, int &peak_idx)
Find largest peak in a vector of chromatograms.
void prepareMasterContainer_(const PeakContainerT &ref_chromatogram, PeakContainerT &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:1056
String boundary_selection_method_
Which method to use for selecting peaks' boundaries.
Definition: MRMTransitionGroupPicker.h:1141
String peak_integration_
Definition: MRMTransitionGroupPicker.h:1120
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1127
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:764
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1121
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1125
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1123
PeakContainerT resampleChromatogram_(const PeakContainerT &chromatogram, const PeakContainerT &master_peak_container, double left_boundary, double right_boundary)
Resample a container at the positions indicated by the master peak container.
Definition: MRMTransitionGroupPicker.h:1094
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:325
double min_peak_width_
Definition: MRMTransitionGroupPicker.h:1132
void remove_overlapping_features(std::vector< SpectrumT > &picked_chroms, double best_left, double best_right)
Remove overlapping features.
Definition: MRMTransitionGroupPicker.h:701
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:948
MRMTransitionGroupPicker & operator=(const MRMTransitionGroupPicker &rhs)
Assignment operator is protected for algorithm.
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1130
PeakPickerChromatogram picker_
Definition: MRMTransitionGroupPicker.h:1143
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1133
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1144
void findWidestPeakIndices(const std::vector< MSChromatogram > &picked_chroms, Int &chrom_idx, Int &point_idx) const
Given a vector of chromatograms, find the indices of the chromatogram containing the widest peak and ...
double min_qual_
Definition: MRMTransitionGroupPicker.h:1128
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1134
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:200
~MRMTransitionGroupPicker() override
Destructor.
void updateMembers_() override
Synchronize members with param class.
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:797
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1131
void pickTransitionGroup(MRMTransitionGroup< SpectrumT, TransitionT > &transition_group)
Pick a group of chromatograms belonging to the same peptide.
Definition: MRMTransitionGroupPicker.h:87
bool use_consensus_
Definition: MRMTransitionGroupPicker.h:1124
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:366
bool recalculate_peaks_
Definition: MRMTransitionGroupPicker.h:1122
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:569
bool compute_peak_shape_metrics_
Definition: MRMTransitionGroupPicker.h:1126
std::vector< ChromatogramType > & getPrecursorChromatograms()
Definition: MRMTransitionGroup.h:211
ChromatogramType & getChromatogram(const String &key)
Definition: MRMTransitionGroup.h:193
bool hasPrecursorChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:242
const String & getTransitionGroupID() const
Definition: MRMTransitionGroup.h:104
bool isInternallyConsistent() const
Check whether internal state is consistent, e.g. same number of chromatograms and transitions are pre...
Definition: MRMTransitionGroup.h:291
bool chromatogramIdsMatch() const
Ensure that chromatogram native ids match their keys in the map.
Definition: MRMTransitionGroup.h:300
bool hasTransition(const String &key) const
Definition: MRMTransitionGroup.h:144
const TransitionType & getTransition(const String &key)
Definition: MRMTransitionGroup.h:149
bool hasChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:188
void addFeature(const MRMFeature &feature)
Definition: MRMTransitionGroup.h:275
const std::vector< TransitionType > & getTransitions() const
Definition: MRMTransitionGroup.h:116
ChromatogramType & getPrecursorChromatogram(const String &key)
Definition: MRMTransitionGroup.h:247
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:160
The representation of a chromatogram.
Definition: MSChromatogram.h:30
void sortByIntensity(bool reverse=false)
Lexicographically sorts the peaks by their intensity.
void setMetaValue(const String &name, const DataValue &value)
Sets the DataValue corresponding to a name.
const DataValue & getMetaValue(const String &name) const
Returns the value corresponding to a string, or DataValue::EMPTY if not found.
void setMZ(CoordinateType coordinate)
Mutable access to the m/z coordinate (index 1)
Definition: Peak2D.h:178
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:190
IntensityType getIntensity() const
Definition: Peak2D.h:142
void setIntensity(IntensityType intensity)
Sets data point intensity (height)
Definition: Peak2D.h:148
Compute the area, background and shape metrics of a peak.
Definition: PeakIntegrator.h:48
Int points_across_half_height
Definition: PeakIntegrator.h:185
double apex_pos
Definition: PeakIntegrator.h:72
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:76
double width_at_5
Definition: PeakIntegrator.h:107
Int points_across_baseline
Definition: PeakIntegrator.h:181
double width_at_50
Definition: PeakIntegrator.h:115
double end_position_at_50
Definition: PeakIntegrator.h:139
double baseline_delta_2_height
Definition: PeakIntegrator.h:177
double height
Definition: PeakIntegrator.h:68
double end_position_at_10
Definition: PeakIntegrator.h:135
double width_at_10
Definition: PeakIntegrator.h:111
double total_width
Definition: PeakIntegrator.h:143
double end_position_at_5
Definition: PeakIntegrator.h:131
double start_position_at_50
Definition: PeakIntegrator.h:127
double tailing_factor
Definition: PeakIntegrator.h:157
double slope_of_baseline
Definition: PeakIntegrator.h:172
double start_position_at_10
Definition: PeakIntegrator.h:123
double asymmetry_factor
Definition: PeakIntegrator.h:167
double start_position_at_5
Definition: PeakIntegrator.h:119
double area
Definition: PeakIntegrator.h:64
Definition: PeakIntegrator.h:60
Definition: PeakIntegrator.h:85
Definition: PeakIntegrator.h:103
The PeakPickerChromatogram finds peaks a single chromatogram.
Definition: PeakPickerChromatogram.h:44
@ IDX_LEFTBORDER
Definition: PeakPickerChromatogram.h:57
@ IDX_RIGHTBORDER
Definition: PeakPickerChromatogram.h:57
@ IDX_ABUNDANCE
Definition: PeakPickerChromatogram.h:57
A more convenient string class.
Definition: String.h:34
Size ensureUniqueId()
Assigns a valid unique id, but only if the present one is invalid. Returns 1 if the unique id was cha...
functor to compute the mean and stddev of sequence using the std::foreach algorithm
Definition: StatsHelpers.h:134
double mean() const
Definition: StatsHelpers.h:171
int Int
Signed integer type.
Definition: Types.h:72
unsigned int UInt
Unsigned integer type.
Definition: Types.h:64
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:97
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:94
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:94
const double k
Definition: Constants.h:132
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19
OPENSWATHALGO_DLLAPI double rankedMutualInformation(std::vector< unsigned int > &ranked_data1, std::vector< unsigned int > &ranked_data2, const unsigned int max_rank1, const unsigned int max_rank2)
OPENSWATHALGO_DLLAPI XCorrArrayType::const_iterator xcorrArrayGetMaxPeak(const XCorrArrayType &array)
Find best peak in an cross-correlation (highest apex)
OPENSWATHALGO_DLLAPI XCorrArrayType normalizedCrossCorrelation(std::vector< double > &data1, std::vector< double > &data2, const int maxdelay, const int lag)
OPENSWATHALGO_DLLAPI unsigned int computeAndAppendRank(const std::vector< double > &v, std::vector< unsigned int > &ranks)
Definition: Scoring.h:46