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