OpenMS
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-2023.
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 
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 << "**** No more peaks left. Nr. chroms: " << 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_) {
195  // If you set this, you only expect one feature anyway. No logging necessary why it stopped.
196  break;
197  }
198  if (intensity > 0 && intensity / total_xic < stop_after_intensity_ratio_)
199  {
200  OPENMS_LOG_DEBUG << "**** Minimum intensity ratio reached. Nr. chroms: " << picked_chroms.size() << std::endl;
201  break;
202  }
203  }
204 
205  // Check for completely overlapping features
206  for (Size i = 0; i < features.size(); i++)
207  {
208  MRMFeature& mrm_feature = features[i];
209  bool skip = false;
210  for (Size j = 0; j < i; j++)
211  {
212  if ((double)mrm_feature.getMetaValue("leftWidth") >= (double)features[j].getMetaValue("leftWidth") &&
213  (double)mrm_feature.getMetaValue("rightWidth") <= (double)features[j].getMetaValue("rightWidth") )
214  { skip = true; }
215  }
216  if (mrm_feature.getIntensity() > 0 && !skip)
217  {
218  transition_group.addFeature(mrm_feature);
219  }
220  }
221 
222  }
223 
225  template <typename SpectrumT, typename TransitionT>
227  std::vector<SpectrumT>& picked_chroms,
228  const std::vector<SpectrumT>& smoothed_chroms,
229  const int chr_idx,
230  const int peak_idx)
231  {
232  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
233  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
234 
235  MRMFeature mrmFeature;
236  mrmFeature.setIntensity(0.0);
237  double best_left = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][peak_idx];
238  double best_right = picked_chroms[chr_idx].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][peak_idx];
239  double peak_apex = picked_chroms[chr_idx][peak_idx].getRT();
240  OPENMS_LOG_DEBUG << "**** Creating MRMFeature for peak " << peak_idx << " in chrom. " << chr_idx << " with " <<
241  picked_chroms[chr_idx][peak_idx] << " and borders " << best_left << " " <<
242  best_right << " (" << best_right - best_left << ")" << std::endl;
243 
244  if (use_consensus_ && recalculate_peaks_)
245  {
246  // This may change best_left / best_right
247  recalculatePeakBorders_(picked_chroms, best_left, best_right, recalculate_peaks_max_z_);
248  if (peak_apex < best_left || peak_apex > best_right)
249  {
250  // apex fell out of range, lets correct it
251  peak_apex = (best_left + best_right) / 2.0;
252  }
253  }
254 
255  std::vector< double > left_edges;
256  std::vector< double > right_edges;
257  double min_left = best_left;
258  double max_right = best_right;
259  if (use_consensus_)
260  {
261  // Remove other, overlapping, picked peaks (in this and other
262  // chromatograms) and then ensure that at least one peak is set to zero
263  // (the currently best peak).
264  remove_overlapping_features(picked_chroms, best_left, best_right);
265  }
266  else
267  {
268  pickApex(picked_chroms, best_left, best_right, peak_apex,
269  min_left, max_right, left_edges, right_edges);
270 
271  } // end !use_consensus_
272  picked_chroms[chr_idx][peak_idx].setIntensity(0.0); // ensure that we set at least one peak to zero
273 
274  // Check for minimal peak width -> return empty feature (Intensity zero)
275  if (use_consensus_)
276  {
277  if (min_peak_width_ > 0.0 && std::fabs(best_right - best_left) < min_peak_width_)
278  {
279  return mrmFeature;
280  }
281 
282  if (compute_peak_quality_)
283  {
284  String outlier = "none";
285  double qual = computeQuality_(transition_group, picked_chroms, chr_idx, best_left, best_right, outlier);
286  if (qual < min_qual_)
287  {
288  return mrmFeature;
289  }
290  mrmFeature.setMetaValue("potentialOutlier", outlier);
291  mrmFeature.setMetaValue("initialPeakQuality", qual);
292  mrmFeature.setOverallQuality(qual);
293  }
294  }
295 
296  // Prepare linear resampling of all the chromatograms, here creating the
297  // empty master_peak_container with the same RT (m/z) values as the
298  // reference chromatogram. We use the overall minimal left boundary and
299  // maximal right boundary to prepare the container.
300  SpectrumT master_peak_container;
301  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
302  prepareMasterContainer_(ref_chromatogram, master_peak_container, min_left, max_right);
303 
304  // Iterate over initial transitions / chromatograms (note that we may
305  // have a different number of picked chromatograms than total transitions
306  // as not all are detecting transitions).
307  double total_intensity = 0; double total_peak_apices = 0; double total_xic = 0; double total_mi = 0;
308  pickFragmentChromatograms(transition_group, picked_chroms, mrmFeature, smoothed_chroms,
309  best_left, best_right, use_consensus_,
310  total_intensity, total_xic, total_mi, total_peak_apices,
311  master_peak_container, left_edges, right_edges,
312  chr_idx, peak_idx);
313 
314  // Also pick the precursor chromatogram(s); note total_xic is not
315  // extracted here, only for fragment traces
316  pickPrecursorChromatograms(transition_group,
317  picked_chroms, mrmFeature, smoothed_chroms,
318  best_left, best_right, use_consensus_,
319  total_intensity, master_peak_container, left_edges, right_edges,
320  chr_idx, peak_idx);
321 
322  mrmFeature.setRT(peak_apex);
323  mrmFeature.setIntensity(total_intensity);
324  mrmFeature.setMetaValue("PeptideRef", transition_group.getTransitionGroupID());
325  mrmFeature.setMetaValue("leftWidth", best_left);
326  mrmFeature.setMetaValue("rightWidth", best_right);
327  mrmFeature.setMetaValue("total_xic", total_xic);
328  if (compute_total_mi_)
329  {
330  mrmFeature.setMetaValue("total_mi", total_mi);
331  }
332  mrmFeature.setMetaValue("peak_apices_sum", total_peak_apices);
333 
334  mrmFeature.ensureUniqueId();
335  return mrmFeature;
336  }
337 
350  template <typename SpectrumT>
351  void pickApex(std::vector<SpectrumT>& picked_chroms,
352  const double best_left, const double best_right, const double peak_apex,
353  double &min_left, double &max_right,
354  std::vector< double > & left_edges, std::vector< double > & right_edges)
355  {
356  for (Size k = 0; k < picked_chroms.size(); k++)
357  {
358  double peak_apex_dist_min = std::numeric_limits<double>::max();
359  int min_dist = -1;
360  for (Size i = 0; i < picked_chroms[k].size(); i++)
361  {
362  PeakIntegrator::PeakArea pa_tmp = pi_.integratePeak( // get the peak apex
363  picked_chroms[k],
364  picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i],
365  picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i]);
366  if (pa_tmp.apex_pos > 1e-11 && std::fabs(pa_tmp.apex_pos - peak_apex) < peak_apex_dist_min)
367  { // update best candidate
368  peak_apex_dist_min = std::fabs(pa_tmp.apex_pos - peak_apex);
369  min_dist = (int)i;
370  }
371  }
372 
373  // Select master peak boundaries, or in the case we found at least one peak, the local peak boundaries
374  double l = best_left;
375  double r = best_right;
376  if (min_dist >= 0)
377  {
378  l = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][min_dist];
379  r = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][min_dist];
380  picked_chroms[k][min_dist].setIntensity(0.0); // only remove one peak per transition
381  }
382 
383  left_edges.push_back(l);
384  right_edges.push_back(r);
385  // ensure we remember the overall maxima / minima
386  if (l < min_left) {min_left = l;}
387  if (r > max_right) {max_right = r;}
388  }
389  }
390 
391  template <typename SpectrumT, typename TransitionT>
393  const std::vector<SpectrumT>& picked_chroms,
394  MRMFeature& mrmFeature,
395  const std::vector<SpectrumT>& smoothed_chroms,
396  const double best_left, const double best_right,
397  const bool use_consensus_,
398  double & total_intensity,
399  double & total_xic,
400  double & total_mi,
401  double & total_peak_apices,
402  const SpectrumT & master_peak_container,
403  const std::vector< double > & left_edges,
404  const std::vector< double > & right_edges,
405  const int chr_idx,
406  const int peak_idx)
407  {
408  for (Size k = 0; k < transition_group.getTransitions().size(); k++)
409  {
410 
411  double local_left = best_left;
412  double local_right = best_right;
413  if (!use_consensus_)
414  {
415  // We cannot have any non-detecting transitions (otherwise we have
416  // too few left / right edges) as we skipped those when doing peak
417  // picking and smoothing.
418  if (!transition_group.getTransitions()[k].isDetectingTransition())
419  {
420  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
421  "When using non-consensus peak picker, all transitions need to be detecting transitions.");
422  }
423  local_left = left_edges[k];
424  local_right = right_edges[k];
425  }
426 
427  const SpectrumT& chromatogram = selectChromHelper_(transition_group, transition_group.getTransitions()[k].getNativeID());
428  if (transition_group.getTransitions()[k].isDetectingTransition())
429  {
430  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
431  {
432  total_xic += it->getIntensity();
433  }
434  }
435 
436  // Compute total intensity on transition-level
437  double transition_total_xic = 0;
438 
439  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
440  {
441  transition_total_xic += it->getIntensity();
442  }
443 
444  // Compute total mutual information on transition-level.
445  double transition_total_mi = 0;
446  if (compute_total_mi_)
447  {
448  std::vector<unsigned int> chrom_vect_id_ranked, chrom_vect_det_ranked;
449  std::vector<double> chrom_vect_id, chrom_vect_det;
450  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
451  {
452  chrom_vect_id.push_back(it->getIntensity());
453  }
454  unsigned int max_rank_det = OpenSwath::Scoring::computeAndAppendRank(chrom_vect_id, chrom_vect_det_ranked);
455  // compute baseline mutual information
456  int transition_total_mi_norm = 0;
457  for (Size m = 0; m < transition_group.getTransitions().size(); m++)
458  {
459  if (transition_group.getTransitions()[m].isDetectingTransition())
460  {
461  const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.getTransitions()[m].getNativeID());
462  chrom_vect_det.clear();
463  for (typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
464  {
465  chrom_vect_det.push_back(it->getIntensity());
466  }
467  unsigned int max_rank_id = OpenSwath::Scoring::computeAndAppendRank(chrom_vect_det, chrom_vect_id_ranked);
468  transition_total_mi += OpenSwath::Scoring::rankedMutualInformation(chrom_vect_id_ranked, chrom_vect_det_ranked, max_rank_id, max_rank_det);
469  transition_total_mi_norm++;
470  }
471  }
472  if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
473 
474  if (transition_group.getTransitions()[k].isDetectingTransition())
475  {
476  // sum up all transition-level total MI and divide by the number of detection transitions to have peak group level total MI
477  total_mi += transition_total_mi / transition_total_mi_norm;
478  }
479  }
480 
481  SpectrumT used_chromatogram;
482  // resample the current chromatogram
483  if (peak_integration_ == "original")
484  {
485  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
486  }
487  else if (peak_integration_ == "smoothed")
488  {
489  if (smoothed_chroms.size() <= k)
490  {
491  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
492  "Tried to calculate peak area and height without any smoothed chromatograms");
493  }
494  used_chromatogram = resampleChromatogram_(smoothed_chroms[k], master_peak_container, local_left, local_right);
495  }
496  else
497  {
498  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
499  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
500  }
501 
502  Feature f;
503  double quality = 0;
504  f.setQuality(0, quality);
505  f.setOverallQuality(quality);
506 
507  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
508  double peak_integral = pa.area;
509  double peak_apex_int = pa.height;
510  f.setMetaValue("peak_apex_position", pa.apex_pos);
511  if (background_subtraction_ != "none")
512  {
513  double background{0};
514  double avg_noise_level{0};
515  if (background_subtraction_ == "original")
516  {
517  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
518  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
519  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
520  avg_noise_level = (intensity_right + intensity_left) / 2;
521  background = avg_noise_level * n_points;
522  }
523  else if (background_subtraction_ == "exact")
524  {
525  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
526  background = pb.area;
527  avg_noise_level = pb.height;
528  }
529  peak_integral -= background;
530  peak_apex_int -= avg_noise_level;
531  if (peak_integral < 0) {peak_integral = 0;}
532  if (peak_apex_int < 0) {peak_apex_int = 0;}
533 
534  f.setMetaValue("area_background_level", background);
535  f.setMetaValue("noise_background_level", avg_noise_level);
536  } // end background
537 
538  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
539  f.setIntensity(peak_integral);
540  ConvexHull2D hull;
541  hull.setHullPoints(pa.hull_points);
542  f.getConvexHulls().push_back(hull);
543 
544  f.setMZ(chromatogram.getProduct().getMZ());
545  mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());
546 
547  if (chromatogram.metaValueExists("product_mz")) // legacy code (ensures that old tests still work)
548  {
549  f.setMetaValue("MZ", chromatogram.getMetaValue("product_mz"));
550  f.setMZ(chromatogram.getMetaValue("product_mz"));
551  }
552 
553  f.setMetaValue("native_id", chromatogram.getNativeID());
554  f.setMetaValue("peak_apex_int", peak_apex_int);
555  f.setMetaValue("total_xic", transition_total_xic);
556  if (compute_total_mi_)
557  {
558  f.setMetaValue("total_mi", transition_total_mi);
559  }
560 
561  if (transition_group.getTransitions()[k].isQuantifyingTransition())
562  {
563  total_intensity += peak_integral;
564  total_peak_apices += peak_apex_int;
565  }
566 
567  // for backwards compatibility with TOPP tests
568  // Calculate peak shape metrics that will be used for later QC
569  PeakIntegrator::PeakShapeMetrics psm = pi_.calculatePeakShapeMetrics(used_chromatogram, local_left, local_right, peak_apex_int, pa.apex_pos);
570  f.setMetaValue("width_at_50", psm.width_at_50);
571  if (compute_peak_shape_metrics_)
572  {
573  f.setMetaValue("width_at_5", psm.width_at_5);
574  f.setMetaValue("width_at_10", psm.width_at_10);
575  f.setMetaValue("start_position_at_5", psm.start_position_at_5);
576  f.setMetaValue("start_position_at_10", psm.start_position_at_10);
577  f.setMetaValue("start_position_at_50", psm.start_position_at_50);
578  f.setMetaValue("end_position_at_5", psm.end_position_at_5);
579  f.setMetaValue("end_position_at_10", psm.end_position_at_10);
580  f.setMetaValue("end_position_at_50", psm.end_position_at_50);
581  f.setMetaValue("total_width", psm.total_width);
582  f.setMetaValue("tailing_factor", psm.tailing_factor);
583  f.setMetaValue("asymmetry_factor", psm.asymmetry_factor);
584  f.setMetaValue("slope_of_baseline", psm.slope_of_baseline);
585  f.setMetaValue("baseline_delta_2_height", psm.baseline_delta_2_height);
586  f.setMetaValue("points_across_baseline", psm.points_across_baseline);
587  f.setMetaValue("points_across_half_height", psm.points_across_half_height);
588  }
589 
590  mrmFeature.addFeature(f, chromatogram.getNativeID()); //map index and feature
591  }
592  }
593 
594  template <typename SpectrumT, typename TransitionT>
596  const std::vector<SpectrumT>& picked_chroms,
597  MRMFeature& mrmFeature,
598  const std::vector<SpectrumT>& smoothed_chroms,
599  const double best_left, const double best_right,
600  const bool use_consensus_,
601  double & total_intensity,
602  const SpectrumT & master_peak_container,
603  const std::vector< double > & left_edges,
604  const std::vector< double > & right_edges,
605  const int chr_idx,
606  const int peak_idx)
607  {
608  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
609  {
610  const SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
611 
612  // Identify precursor index
613  // note: this is only valid if all transitions are detecting transitions
614  Size prec_idx = transition_group.getChromatograms().size() + k;
615 
616  double local_left = best_left;
617  double local_right = best_right;
618  if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
619  {
620  local_left = left_edges[prec_idx];
621  local_right = right_edges[prec_idx];
622  }
623 
624  SpectrumT used_chromatogram;
625  // resample the current chromatogram
626  if (peak_integration_ == "original")
627  {
628  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
629  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
630  }
631  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= prec_idx)
632  {
633  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
634  "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
635  }
636  else if (peak_integration_ == "smoothed")
637  {
638  used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
639  }
640  else
641  {
642  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
643  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
644  }
645 
646  Feature f;
647  double quality = 0;
648  f.setQuality(0, quality);
649  f.setOverallQuality(quality);
650 
651  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
652  double peak_integral = pa.area;
653  double peak_apex_int = pa.height;
654 
655  if (background_subtraction_ != "none")
656  {
657  double background{0};
658  double avg_noise_level{0};
659  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= prec_idx)
660  {
661  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
662  "Tried to calculate background estimation without any smoothed chromatograms");
663  }
664  else if (background_subtraction_ == "original")
665  {
666  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
667  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
668  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
669  avg_noise_level = (intensity_right + intensity_left) / 2;
670  background = avg_noise_level * n_points;
671  }
672  else if (background_subtraction_ == "exact")
673  {
674  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
675  background = pb.area;
676  avg_noise_level = pb.height;
677  }
678  peak_integral -= background;
679  peak_apex_int -= avg_noise_level;
680  if (peak_integral < 0) {peak_integral = 0;}
681  if (peak_apex_int < 0) {peak_apex_int = 0;}
682 
683  f.setMetaValue("area_background_level", background);
684  f.setMetaValue("noise_background_level", avg_noise_level);
685  }
686 
687  f.setMZ(chromatogram.getPrecursor().getMZ());
688  if (k == 0) {mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());} // only use m/z if first (monoisotopic) isotope
689 
690  if (chromatogram.metaValueExists("precursor_mz")) // legacy code (ensures that old tests still work)
691  {
692  f.setMZ(chromatogram.getMetaValue("precursor_mz"));
693  if (k == 0) {mrmFeature.setMZ(chromatogram.getMetaValue("precursor_mz"));} // only use m/z if first (monoisotopic) isotope
694  }
695 
696  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
697  f.setIntensity(peak_integral);
698  ConvexHull2D hull;
699  hull.setHullPoints(pa.hull_points);
700  f.getConvexHulls().push_back(hull);
701  f.setMetaValue("native_id", chromatogram.getNativeID());
702  f.setMetaValue("peak_apex_int", peak_apex_int);
703 
704  if (use_precursors_ && transition_group.getTransitions().empty())
705  {
706  total_intensity += peak_integral;
707  }
708 
709  mrmFeature.addPrecursorFeature(f, chromatogram.getNativeID());
710  }
711  }
712 
713  // maybe private, but we have tests
714 
726  template <typename SpectrumT>
727  void remove_overlapping_features(std::vector<SpectrumT>& picked_chroms, double best_left, double best_right)
728  {
729  // delete all seeds that lie within the current seed
730  Size count_inside = 0;
731  for (Size k = 0; k < picked_chroms.size(); k++)
732  {
733  for (Size i = 0; i < picked_chroms[k].size(); i++)
734  {
735  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
736  {
737  picked_chroms[k][i].setIntensity(0.0);
738  count_inside++;
739  }
740  }
741  }
742 
743  Size count_overlap = 0;
744  // delete all seeds that overlap within the current seed
745  for (Size k = 0; k < picked_chroms.size(); k++)
746  {
747  for (Size i = 0; i < picked_chroms[k].size(); i++)
748  {
749  if (picked_chroms[k][i].getIntensity() <= 1e-11) {continue; }
750 
751  double left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
752  double right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
753  if ((left > best_left && left < best_right)
754  || (right > best_left && right < best_right))
755  {
756  picked_chroms[k][i].setIntensity(0.0);
757  count_overlap++;
758  }
759  }
760  }
761  OPENMS_LOG_DEBUG << " ** Removed " << count_inside << " peaks enclosed in and " <<
762  count_overlap << " peaks overlapping with added feature." << std::endl;
763  }
764 
766  void findLargestPeak(const std::vector<MSChromatogram >& picked_chroms, int& chr_idx, int& peak_idx);
767 
776  void findWidestPeakIndices(const std::vector<MSChromatogram>& picked_chroms, Int& chrom_idx, Int& point_idx) const;
777 
778 protected:
779 
781  void updateMembers_() override;
782 
785 
789  template <typename SpectrumT, typename TransitionT>
790  const SpectrumT& selectChromHelper_(const MRMTransitionGroup<SpectrumT, TransitionT>& transition_group, const String& native_id)
791  {
792  if (transition_group.hasChromatogram(native_id))
793  {
794  return transition_group.getChromatogram(native_id);
795  }
796  else if (transition_group.hasPrecursorChromatogram(native_id))
797  {
798  return transition_group.getPrecursorChromatogram(native_id);
799  }
800  else
801  {
802  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION, "Did not find chromatogram for id '" + native_id + "'.");
803  }
804  }
805 
822  template <typename SpectrumT, typename TransitionT>
824  const std::vector<SpectrumT>& picked_chroms,
825  const int chr_idx,
826  const double best_left,
827  const double best_right,
828  String& outlier)
829  {
830  // Resample all chromatograms around the current estimated peak and
831  // collect the raw intensities. For resampling, use a bit more on either
832  // side to correctly identify shoulders etc.
833  double resample_boundary = resample_boundary_; // sample 15 seconds more on each side
834  SpectrumT master_peak_container;
835  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
836  prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
837  std::vector<std::vector<double> > all_ints;
838  for (Size k = 0; k < picked_chroms.size(); k++)
839  {
840  const SpectrumT& chromatogram = selectChromHelper_(transition_group, picked_chroms[k].getNativeID());
841  const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
842  master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
843 
844  std::vector<double> int_here;
845  for (const auto& peak : used_chromatogram) int_here.push_back(peak.getIntensity());
846  // Remove chromatograms without a single peak
847  double tic = std::accumulate(int_here.begin(), int_here.end(), 0.0);
848  if (tic > 1e-11) all_ints.push_back(int_here);
849  }
850 
851  // Compute the cross-correlation for the collected intensities
852  std::vector<double> mean_shapes;
853  std::vector<double> mean_coel;
854  for (Size k = 0; k < all_ints.size(); k++)
855  {
856  std::vector<double> shapes;
857  std::vector<double> coel;
858  for (Size i = 0; i < all_ints.size(); i++)
859  {
860  if (i == k) {continue;}
862  all_ints[k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
863 
864  // the first value is the x-axis (retention time) and should be an int -> it show the lag between the two
865  double res_coelution = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->first);
866  double res_shape = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->second);
867 
868  shapes.push_back(res_shape);
869  coel.push_back(res_coelution);
870  }
871 
872  // We have computed the cross-correlation of chromatogram k against
873  // all others. Use the mean of these computations as the value for k.
875  msc = std::for_each(shapes.begin(), shapes.end(), msc);
876  double shapes_mean = msc.mean();
877  msc = std::for_each(coel.begin(), coel.end(), msc);
878  double coel_mean = msc.mean();
879 
880  // mean shape scores below 0.5-0.6 should be a real sign of trouble ... !
881  // mean coel scores above 3.5 should be a real sign of trouble ... !
882  mean_shapes.push_back(shapes_mean);
883  mean_coel.push_back(coel_mean);
884  }
885 
886  // find the chromatogram with the minimal shape score and the maximal
887  // coelution score -> if it is the same chromatogram, the chance is
888  // pretty good that it is different from the others...
889  int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
890  int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
891 
892  // Look at the picked peaks that are within the current left/right borders
893  int missing_peaks = 0;
894  int multiple_peaks = 0;
895 
896  // collect all seeds that lie within the current seed
897  std::vector<double> left_borders;
898  std::vector<double> right_borders;
899  for (Size k = 0; k < picked_chroms.size(); k++)
900  {
901  double l_tmp;
902  double r_tmp;
903  double max_int = -1;
904 
905  int pfound = 0;
906  l_tmp = -1;
907  r_tmp = -1;
908  for (Size i = 0; i < picked_chroms[k].size(); i++)
909  {
910  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
911  {
912  pfound++;
913  if (picked_chroms[k][i].getIntensity() > max_int)
914  {
915  max_int = picked_chroms[k][i].getIntensity() > max_int;
916  l_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
917  r_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
918  }
919  }
920  }
921 
922  if (l_tmp > 1e-11) left_borders.push_back(l_tmp);
923  if (r_tmp > 1e-11) right_borders.push_back(r_tmp);
924 
925  if (pfound == 0) missing_peaks++;
926  if (pfound > 1) multiple_peaks++;
927  }
928 
929  // Check how many chromatograms had exactly one peak picked between our
930  // current left/right borders -> this would be a sign of consistency.
931  OPENMS_LOG_DEBUG << " Overall found missing : " << missing_peaks << " and multiple : " << multiple_peaks << std::endl;
932 
934 
935  // Is there one transitions that is very different from the rest (e.g.
936  // the same element has a bad shape and a bad coelution score) -> potential outlier
937  if (min_index_shape == max_index_coel)
938  {
939  OPENMS_LOG_DEBUG << " Element " << min_index_shape << " is a candidate for removal ... " << std::endl;
940  outlier = String(picked_chroms[min_index_shape].getNativeID());
941  }
942  else
943  {
944  outlier = "none";
945  }
946 
947  // For the final score (larger is better), consider these scores:
948  // - missing_peaks (the more peaks are missing, the worse)
949  // - multiple_peaks
950  // - mean of the shapes (1 is very good, 0 is bad)
951  // - mean of the co-elution scores (0 is good, 1 is ok, above 1 is pretty bad)
952  double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
953  double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
954  coel_score = (coel_score - 1.0) / 2.0;
955 
956  double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
957 
958  OPENMS_LOG_DEBUG << " Computed score " << score << " (from " << shape_score <<
959  " - " << coel_score << " - " << 1.0 * missing_peaks / picked_chroms.size() << ")" << std::endl;
960 
961  return score;
962  }
963 
973  template <typename SpectrumT>
974  void recalculatePeakBorders_(const std::vector<SpectrumT>& picked_chroms, double& best_left, double& best_right, double max_z)
975  {
976  // 1. Collect all seeds that lie within the current seed
977  // - Per chromatogram only the most intense one counts, otherwise very
978  // - low intense peaks can contribute disproportionally to the voting
979  // - procedure.
980  // TODO Especially with DDA MS1 "transitions" from FFID, you often have exactly the same
981  // borders and the logs get confusing and you get -Nan CVs. You also might save time if you use a set here.
982  std::vector<double> left_borders;
983  std::vector<double> right_borders;
984  for (Size k = 0; k < picked_chroms.size(); k++)
985  {
986  double max_int = -1;
987  double left = -1;
988  double right = -1;
989  for (Size i = 0; i < picked_chroms[k].size(); i++)
990  {
991  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
992  {
993  if (picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i] > max_int)
994  {
995  max_int = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i];
996  left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
997  right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
998  }
999  }
1000  }
1001  if (max_int > -1 )
1002  {
1003  left_borders.push_back(left);
1004  right_borders.push_back(right);
1005  OPENMS_LOG_DEBUG << " * peak " << k << " left boundary " << left_borders.back() << " with inty " << max_int << std::endl;
1006  OPENMS_LOG_DEBUG << " * peak " << k << " right boundary " << right_borders.back() << " with inty " << max_int << std::endl;
1007  }
1008  }
1009 
1010  // Return for empty peak list
1011  if (right_borders.empty())
1012  {
1013  return;
1014  }
1015 
1016  // FEATURE IDEA: instead of Z-score use modified Z-score for small data sets
1017  // http://d-scholarship.pitt.edu/7948/1/Seo.pdf
1018  // http://www.itl.nist.gov/div898/handbook/eda/section3/eda35h.htm
1019  // 1. calculate median
1020  // 2. MAD = calculate difference to median for each value -> take median of that
1021  // 3. Mi = 0.6745*(xi - median) / MAD
1022 
1023  // 2. Calculate mean and standard deviation
1024  // If the coefficient of variation is too large for one border, we use a
1025  // "pseudo-median" instead of the border of the most intense peak.
1026  double mean, stdev;
1027 
1028  // Right borders
1029  mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
1030  stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
1031  / right_borders.size() - mean * mean);
1032  std::sort(right_borders.begin(), right_borders.end());
1033 
1034  OPENMS_LOG_DEBUG << " - Recalculating right peak boundaries " << mean << " mean / best "
1035  << best_right << " std " << stdev << " : " << std::fabs(best_right - mean) / stdev
1036  << " coefficient of variation" << std::endl;
1037 
1038  // Compare right borders of best transition with the mean
1039  if (std::fabs(best_right - mean) / stdev > max_z)
1040  {
1041  best_right = right_borders[right_borders.size() / 2]; // pseudo median
1042  OPENMS_LOG_DEBUG << " - CV too large: correcting right boundary to " << best_right << std::endl;
1043  }
1044 
1045  // Left borders
1046  mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
1047  stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
1048  / left_borders.size() - mean * mean);
1049  std::sort(left_borders.begin(), left_borders.end());
1050 
1051  OPENMS_LOG_DEBUG << " - Recalculating left peak boundaries " << mean << " mean / best "
1052  << best_left << " std " << stdev << " : " << std::fabs(best_left - mean) / stdev
1053  << " coefficient of variation" << std::endl;
1054 
1055  // Compare left borders of best transition with the mean
1056  if (std::fabs(best_left - mean) / stdev > max_z)
1057  {
1058  best_left = left_borders[left_borders.size() / 2]; // pseudo median
1059  OPENMS_LOG_DEBUG << " - CV too large: correcting left boundary to " << best_left << std::endl;
1060  }
1061 
1062  }
1063 
1065 
1066 
1081  template <typename SpectrumT>
1082  void prepareMasterContainer_(const SpectrumT& ref_chromatogram,
1083  SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1084  {
1085  OPENMS_PRECONDITION(master_peak_container.empty(), "Master peak container must be empty")
1086 
1087  // get the start / end point of this chromatogram => then add one more
1088  // point beyond the two boundaries to make the resampling accurate also
1089  // at the edge.
1090  typename SpectrumT::const_iterator begin = ref_chromatogram.begin();
1091  while (begin != ref_chromatogram.end() && begin->getMZ() < left_boundary) {begin++; }
1092  if (begin != ref_chromatogram.begin()) {begin--; }
1093 
1094  typename SpectrumT::const_iterator end = begin;
1095  while (end != ref_chromatogram.end() && end->getMZ() < right_boundary) {end++; }
1096  if (end != ref_chromatogram.end()) {end++; }
1097 
1098  // resize the master container and set the m/z values to the ones of the master container
1099  master_peak_container.resize(distance(begin, end)); // initialize to zero
1100  typename SpectrumT::iterator it = master_peak_container.begin();
1101  for (typename SpectrumT::const_iterator chrom_it = begin; chrom_it != end; chrom_it++, it++)
1102  {
1103  it->setMZ(chrom_it->getMZ());
1104  }
1105  }
1106 
1117  template <typename SpectrumT>
1118  SpectrumT resampleChromatogram_(const SpectrumT& chromatogram,
1119  const SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1120  {
1121  // get the start / end point of this chromatogram => then add one more
1122  // point beyond the two boundaries to make the resampling accurate also
1123  // at the edge.
1124  typename SpectrumT::const_iterator begin = chromatogram.begin();
1125  while (begin != chromatogram.end() && begin->getMZ() < left_boundary) {begin++;}
1126  if (begin != chromatogram.begin()) {begin--;}
1127 
1128  typename SpectrumT::const_iterator end = begin;
1129  while (end != chromatogram.end() && end->getMZ() < right_boundary) {end++;}
1130  if (end != chromatogram.end()) {end++;}
1131 
1132  SpectrumT resampled_peak_container = master_peak_container; // copy the master container, which contains the RT values
1133  LinearResamplerAlign lresampler;
1134  lresampler.raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1135 
1136  return resampled_peak_container;
1137  }
1138 
1140 
1141  // Members
1150  double min_qual_;
1151 
1157 
1164 
1167  };
1168 }
1169 
1170 
#define OPENMS_LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:480
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.
A 2-dimensional hull representation in [counter]clockwise direction - depending on axis labelling.
Definition: ConvexHull2D.h:73
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:92
A method or algorithm argument contains illegal values.
Definition: Exception.h:650
An LC-MS feature.
Definition: Feature.h:72
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:59
void raster(SpecT &container)
Applies the resampling algorithm to a container (MSSpectrum or MSChromatogram).
Definition: LinearResamplerAlign.h:80
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:52
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:81
void findLargestPeak(const std::vector< MSChromatogram > &picked_chroms, int &chr_idx, int &peak_idx)
Find largest peak in a vector of chromatograms.
String boundary_selection_method_
Which method to use for selecting peaks' boundaries.
Definition: MRMTransitionGroupPicker.h:1163
String peak_integration_
Definition: MRMTransitionGroupPicker.h:1142
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1149
PeakPickerMRM picker_
Definition: MRMTransitionGroupPicker.h:1165
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:790
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:1082
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1143
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1147
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1145
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:351
double min_peak_width_
Definition: MRMTransitionGroupPicker.h:1154
void remove_overlapping_features(std::vector< SpectrumT > &picked_chroms, double best_left, double best_right)
Remove overlapping features.
Definition: MRMTransitionGroupPicker.h:727
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:974
MRMTransitionGroupPicker & operator=(const MRMTransitionGroupPicker &rhs)
Assignment operator is protected for algorithm.
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1152
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:1118
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1155
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1166
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:1150
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1156
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:226
~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:823
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1153
void pickTransitionGroup(MRMTransitionGroup< SpectrumT, TransitionT > &transition_group)
Pick a group of chromatograms belonging to the same peptide.
Definition: MRMTransitionGroupPicker.h:113
bool use_consensus_
Definition: MRMTransitionGroupPicker.h:1146
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:392
bool recalculate_peaks_
Definition: MRMTransitionGroupPicker.h:1144
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:595
bool compute_peak_shape_metrics_
Definition: MRMTransitionGroupPicker.h:1148
std::vector< ChromatogramType > & getPrecursorChromatograms()
Definition: MRMTransitionGroup.h:237
ChromatogramType & getChromatogram(const String &key)
Definition: MRMTransitionGroup.h:219
bool hasPrecursorChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:267
const String & getTransitionGroupID() const
Definition: MRMTransitionGroup.h:130
bool isInternallyConsistent() const
Check whether internal state is consistent, e.g. same number of chromatograms and transitions are pre...
Definition: MRMTransitionGroup.h:316
bool chromatogramIdsMatch() const
Ensure that chromatogram native ids match their keys in the map.
Definition: MRMTransitionGroup.h:325
bool hasTransition(const String &key) const
Definition: MRMTransitionGroup.h:170
const TransitionType & getTransition(const String &key)
Definition: MRMTransitionGroup.h:175
bool hasChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:214
void addFeature(const MRMFeature &feature)
Definition: MRMTransitionGroup.h:300
const std::vector< TransitionType > & getTransitions() const
Definition: MRMTransitionGroup.h:142
ChromatogramType & getPrecursorChromatogram(const String &key)
Definition: MRMTransitionGroup.h:272
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:186
The representation of a chromatogram.
Definition: MSChromatogram.h:57
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:204
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:216
IntensityType getIntensity() const
Definition: Peak2D.h:168
void setIntensity(IntensityType intensity)
Sets data point intensity (height)
Definition: Peak2D.h:174
Compute the area, background and shape metrics of a peak.
Definition: PeakIntegrator.h:74
Int points_across_half_height
Definition: PeakIntegrator.h:211
double apex_pos
Definition: PeakIntegrator.h:98
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:102
double width_at_5
Definition: PeakIntegrator.h:133
Int points_across_baseline
Definition: PeakIntegrator.h:207
double width_at_50
Definition: PeakIntegrator.h:141
double end_position_at_50
Definition: PeakIntegrator.h:165
double baseline_delta_2_height
Definition: PeakIntegrator.h:203
double height
Definition: PeakIntegrator.h:94
double end_position_at_10
Definition: PeakIntegrator.h:161
double width_at_10
Definition: PeakIntegrator.h:137
double total_width
Definition: PeakIntegrator.h:169
double end_position_at_5
Definition: PeakIntegrator.h:157
double start_position_at_50
Definition: PeakIntegrator.h:153
double tailing_factor
Definition: PeakIntegrator.h:183
double slope_of_baseline
Definition: PeakIntegrator.h:198
double start_position_at_10
Definition: PeakIntegrator.h:149
double asymmetry_factor
Definition: PeakIntegrator.h:193
double start_position_at_5
Definition: PeakIntegrator.h:145
double area
Definition: PeakIntegrator.h:90
Definition: PeakIntegrator.h:86
Definition: PeakIntegrator.h:111
Definition: PeakIntegrator.h:129
The PeakPickerMRM finds peaks a single chromatogram.
Definition: PeakPickerMRM.h:70
@ IDX_LEFTBORDER
Definition: PeakPickerMRM.h:83
@ IDX_RIGHTBORDER
Definition: PeakPickerMRM.h:83
@ IDX_ABUNDANCE
Definition: PeakPickerMRM.h:83
A more convenient string class.
Definition: String.h:60
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:160
double mean() const
Definition: StatsHelpers.h:197
int Int
Signed integer type.
Definition: Types.h:102
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:120
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:120
const double k
Definition: Constants.h:158
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:48
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:72