OpenMS  2.7.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-2021.
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 > 0.0 && 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<double> chrom_vect_id, chrom_vect_det;
449  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
450  {
451  chrom_vect_id.push_back(it->getIntensity());
452  }
453 
454  // compute baseline mutual information
455  int transition_total_mi_norm = 0;
456  for (Size m = 0; m < transition_group.getTransitions().size(); m++)
457  {
458  if (transition_group.getTransitions()[m].isDetectingTransition())
459  {
460  const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.getTransitions()[m].getNativeID());
461  chrom_vect_det.clear();
462  for (typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
463  {
464  chrom_vect_det.push_back(it->getIntensity());
465  }
466  transition_total_mi += OpenSwath::Scoring::rankedMutualInformation(chrom_vect_det, chrom_vect_id);
467  transition_total_mi_norm++;
468  }
469  }
470  if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
471 
472  if (transition_group.getTransitions()[k].isDetectingTransition())
473  {
474  // sum up all transition-level total MI and divide by the number of detection transitions to have peak group level total MI
475  total_mi += transition_total_mi / transition_total_mi_norm;
476  }
477  }
478 
479  SpectrumT used_chromatogram;
480  // resample the current chromatogram
481  if (peak_integration_ == "original")
482  {
483  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
484  }
485  else if (peak_integration_ == "smoothed")
486  {
487  if (smoothed_chroms.size() <= k)
488  {
489  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
490  "Tried to calculate peak area and height without any smoothed chromatograms");
491  }
492  used_chromatogram = resampleChromatogram_(smoothed_chroms[k], master_peak_container, local_left, local_right);
493  }
494  else
495  {
496  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
497  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
498  }
499 
500  Feature f;
501  double quality = 0;
502  f.setQuality(0, quality);
503  f.setOverallQuality(quality);
504 
505  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
506  double peak_integral = pa.area;
507  double peak_apex_int = pa.height;
508  f.setMetaValue("peak_apex_position", pa.apex_pos);
509  if (background_subtraction_ != "none")
510  {
511  double background{0};
512  double avg_noise_level{0};
513  if (background_subtraction_ == "original")
514  {
515  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
516  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
517  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
518  avg_noise_level = (intensity_right + intensity_left) / 2;
519  background = avg_noise_level * n_points;
520  }
521  else if (background_subtraction_ == "exact")
522  {
523  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
524  background = pb.area;
525  avg_noise_level = pb.height;
526  }
527  peak_integral -= background;
528  peak_apex_int -= avg_noise_level;
529  if (peak_integral < 0) {peak_integral = 0;}
530  if (peak_apex_int < 0) {peak_apex_int = 0;}
531 
532  f.setMetaValue("area_background_level", background);
533  f.setMetaValue("noise_background_level", avg_noise_level);
534  } // end background
535 
536  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
537  f.setIntensity(peak_integral);
538  ConvexHull2D hull;
539  hull.setHullPoints(pa.hull_points);
540  f.getConvexHulls().push_back(hull);
541 
542  f.setMZ(chromatogram.getProduct().getMZ());
543  mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());
544 
545  if (chromatogram.metaValueExists("product_mz")) // legacy code (ensures that old tests still work)
546  {
547  f.setMetaValue("MZ", chromatogram.getMetaValue("product_mz"));
548  f.setMZ(chromatogram.getMetaValue("product_mz"));
549  }
550 
551  f.setMetaValue("native_id", chromatogram.getNativeID());
552  f.setMetaValue("peak_apex_int", peak_apex_int);
553  f.setMetaValue("total_xic", transition_total_xic);
554  if (compute_total_mi_)
555  {
556  f.setMetaValue("total_mi", transition_total_mi);
557  }
558 
559  if (transition_group.getTransitions()[k].isQuantifyingTransition())
560  {
561  total_intensity += peak_integral;
562  total_peak_apices += peak_apex_int;
563  }
564 
565  // for backwards compatibility with TOPP tests
566  // Calculate peak shape metrics that will be used for later QC
567  PeakIntegrator::PeakShapeMetrics psm = pi_.calculatePeakShapeMetrics(used_chromatogram, local_left, local_right, peak_apex_int, pa.apex_pos);
568  f.setMetaValue("width_at_50", psm.width_at_50);
569  if (compute_peak_shape_metrics_)
570  {
571  f.setMetaValue("width_at_5", psm.width_at_5);
572  f.setMetaValue("width_at_10", psm.width_at_10);
573  f.setMetaValue("start_position_at_5", psm.start_position_at_5);
574  f.setMetaValue("start_position_at_10", psm.start_position_at_10);
575  f.setMetaValue("start_position_at_50", psm.start_position_at_50);
576  f.setMetaValue("end_position_at_5", psm.end_position_at_5);
577  f.setMetaValue("end_position_at_10", psm.end_position_at_10);
578  f.setMetaValue("end_position_at_50", psm.end_position_at_50);
579  f.setMetaValue("total_width", psm.total_width);
580  f.setMetaValue("tailing_factor", psm.tailing_factor);
581  f.setMetaValue("asymmetry_factor", psm.asymmetry_factor);
582  f.setMetaValue("slope_of_baseline", psm.slope_of_baseline);
583  f.setMetaValue("baseline_delta_2_height", psm.baseline_delta_2_height);
584  f.setMetaValue("points_across_baseline", psm.points_across_baseline);
585  f.setMetaValue("points_across_half_height", psm.points_across_half_height);
586  }
587 
588  mrmFeature.addFeature(f, chromatogram.getNativeID()); //map index and feature
589  }
590  }
591 
592  template <typename SpectrumT, typename TransitionT>
594  const std::vector<SpectrumT>& picked_chroms,
595  MRMFeature& mrmFeature,
596  const std::vector<SpectrumT>& smoothed_chroms,
597  const double best_left, const double best_right,
598  const bool use_consensus_,
599  double & total_intensity,
600  const SpectrumT & master_peak_container,
601  const std::vector< double > & left_edges,
602  const std::vector< double > & right_edges,
603  const int chr_idx,
604  const int peak_idx)
605  {
606  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
607  {
608  const SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
609 
610  // Identify precursor index
611  // note: this is only valid if all transitions are detecting transitions
612  Size prec_idx = transition_group.getChromatograms().size() + k;
613 
614  double local_left = best_left;
615  double local_right = best_right;
616  if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
617  {
618  local_left = left_edges[prec_idx];
619  local_right = right_edges[prec_idx];
620  }
621 
622  SpectrumT used_chromatogram;
623  // resample the current chromatogram
624  if (peak_integration_ == "original")
625  {
626  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
627  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
628  }
629  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= prec_idx)
630  {
631  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
632  "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
633  }
634  else if (peak_integration_ == "smoothed")
635  {
636  used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
637  }
638  else
639  {
640  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
641  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
642  }
643 
644  Feature f;
645  double quality = 0;
646  f.setQuality(0, quality);
647  f.setOverallQuality(quality);
648 
649  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
650  double peak_integral = pa.area;
651  double peak_apex_int = pa.height;
652 
653  if (background_subtraction_ != "none")
654  {
655  double background{0};
656  double avg_noise_level{0};
657  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= prec_idx)
658  {
659  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
660  "Tried to calculate background estimation without any smoothed chromatograms");
661  }
662  else if (background_subtraction_ == "original")
663  {
664  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
665  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
666  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
667  avg_noise_level = (intensity_right + intensity_left) / 2;
668  background = avg_noise_level * n_points;
669  }
670  else if (background_subtraction_ == "exact")
671  {
672  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
673  background = pb.area;
674  avg_noise_level = pb.height;
675  }
676  peak_integral -= background;
677  peak_apex_int -= avg_noise_level;
678  if (peak_integral < 0) {peak_integral = 0;}
679  if (peak_apex_int < 0) {peak_apex_int = 0;}
680 
681  f.setMetaValue("area_background_level", background);
682  f.setMetaValue("noise_background_level", avg_noise_level);
683  }
684 
685  f.setMZ(chromatogram.getPrecursor().getMZ());
686  if (k == 0) {mrmFeature.setMZ(chromatogram.getPrecursor().getMZ());} // only use m/z if first (monoisotopic) isotope
687 
688  if (chromatogram.metaValueExists("precursor_mz")) // legacy code (ensures that old tests still work)
689  {
690  f.setMZ(chromatogram.getMetaValue("precursor_mz"));
691  if (k == 0) {mrmFeature.setMZ(chromatogram.getMetaValue("precursor_mz"));} // only use m/z if first (monoisotopic) isotope
692  }
693 
694  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
695  f.setIntensity(peak_integral);
696  ConvexHull2D hull;
697  hull.setHullPoints(pa.hull_points);
698  f.getConvexHulls().push_back(hull);
699  f.setMetaValue("native_id", chromatogram.getNativeID());
700  f.setMetaValue("peak_apex_int", peak_apex_int);
701 
702  if (use_precursors_ && transition_group.getTransitions().empty())
703  {
704  total_intensity += peak_integral;
705  }
706 
707  mrmFeature.addPrecursorFeature(f, chromatogram.getNativeID());
708  }
709  }
710 
711  // maybe private, but we have tests
712 
724  template <typename SpectrumT>
725  void remove_overlapping_features(std::vector<SpectrumT>& picked_chroms, double best_left, double best_right)
726  {
727  // delete all seeds that lie within the current seed
728  Size count_inside = 0;
729  for (Size k = 0; k < picked_chroms.size(); k++)
730  {
731  for (Size i = 0; i < picked_chroms[k].size(); i++)
732  {
733  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
734  {
735  picked_chroms[k][i].setIntensity(0.0);
736  count_inside++;
737  }
738  }
739  }
740 
741  Size count_overlap = 0;
742  // delete all seeds that overlap within the current seed
743  for (Size k = 0; k < picked_chroms.size(); k++)
744  {
745  for (Size i = 0; i < picked_chroms[k].size(); i++)
746  {
747  if (picked_chroms[k][i].getIntensity() <= 0.0) {continue; }
748 
749  double left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
750  double right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
751  if ((left > best_left && left < best_right)
752  || (right > best_left && right < best_right))
753  {
754  picked_chroms[k][i].setIntensity(0.0);
755  count_overlap++;
756  }
757  }
758  }
759  OPENMS_LOG_DEBUG << " ** Removed " << count_inside << " peaks enclosed in and " <<
760  count_overlap << " peaks overlapping with added feature." << std::endl;
761  }
762 
764  void findLargestPeak(const std::vector<MSChromatogram >& picked_chroms, int& chr_idx, int& peak_idx);
765 
774  void findWidestPeakIndices(const std::vector<MSChromatogram>& picked_chroms, Int& chrom_idx, Int& point_idx) const;
775 
776 protected:
777 
779  void updateMembers_() override;
780 
783 
787  template <typename SpectrumT, typename TransitionT>
788  const SpectrumT& selectChromHelper_(const MRMTransitionGroup<SpectrumT, TransitionT>& transition_group, const String& native_id)
789  {
790  if (transition_group.hasChromatogram(native_id))
791  {
792  return transition_group.getChromatogram(native_id);
793  }
794  else if (transition_group.hasPrecursorChromatogram(native_id))
795  {
796  return transition_group.getPrecursorChromatogram(native_id);
797  }
798  else
799  {
800  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION, "Did not find chromatogram for id '" + native_id + "'.");
801  }
802  }
803 
820  template <typename SpectrumT, typename TransitionT>
822  const std::vector<SpectrumT>& picked_chroms,
823  const int chr_idx,
824  const double best_left,
825  const double best_right,
826  String& outlier)
827  {
828  // Resample all chromatograms around the current estimated peak and
829  // collect the raw intensities. For resampling, use a bit more on either
830  // side to correctly identify shoulders etc.
831  double resample_boundary = resample_boundary_; // sample 15 seconds more on each side
832  SpectrumT master_peak_container;
833  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
834  prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
835  std::vector<std::vector<double> > all_ints;
836  for (Size k = 0; k < picked_chroms.size(); k++)
837  {
838  const SpectrumT& chromatogram = selectChromHelper_(transition_group, picked_chroms[k].getNativeID());
839  const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
840  master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
841 
842  std::vector<double> int_here;
843  for (const auto& peak : used_chromatogram) int_here.push_back(peak.getIntensity());
844  // Remove chromatograms without a single peak
845  double tic = std::accumulate(int_here.begin(), int_here.end(), 0.0);
846  if (tic > 0.0) all_ints.push_back(int_here);
847  }
848 
849  // Compute the cross-correlation for the collected intensities
850  std::vector<double> mean_shapes;
851  std::vector<double> mean_coel;
852  for (Size k = 0; k < all_ints.size(); k++)
853  {
854  std::vector<double> shapes;
855  std::vector<double> coel;
856  for (Size i = 0; i < all_ints.size(); i++)
857  {
858  if (i == k) {continue;}
860  all_ints[k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
861 
862  // the first value is the x-axis (retention time) and should be an int -> it show the lag between the two
863  double res_coelution = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->first);
864  double res_shape = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->second);
865 
866  shapes.push_back(res_shape);
867  coel.push_back(res_coelution);
868  }
869 
870  // We have computed the cross-correlation of chromatogram k against
871  // all others. Use the mean of these computations as the value for k.
873  msc = std::for_each(shapes.begin(), shapes.end(), msc);
874  double shapes_mean = msc.mean();
875  msc = std::for_each(coel.begin(), coel.end(), msc);
876  double coel_mean = msc.mean();
877 
878  // mean shape scores below 0.5-0.6 should be a real sign of trouble ... !
879  // mean coel scores above 3.5 should be a real sign of trouble ... !
880  mean_shapes.push_back(shapes_mean);
881  mean_coel.push_back(coel_mean);
882  }
883 
884  // find the chromatogram with the minimal shape score and the maximal
885  // coelution score -> if it is the same chromatogram, the chance is
886  // pretty good that it is different from the others...
887  int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
888  int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
889 
890  // Look at the picked peaks that are within the current left/right borders
891  int missing_peaks = 0;
892  int multiple_peaks = 0;
893 
894  // collect all seeds that lie within the current seed
895  std::vector<double> left_borders;
896  std::vector<double> right_borders;
897  for (Size k = 0; k < picked_chroms.size(); k++)
898  {
899  double l_tmp;
900  double r_tmp;
901  double max_int = -1;
902 
903  int pfound = 0;
904  l_tmp = -1;
905  r_tmp = -1;
906  for (Size i = 0; i < picked_chroms[k].size(); i++)
907  {
908  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
909  {
910  pfound++;
911  if (picked_chroms[k][i].getIntensity() > max_int)
912  {
913  max_int = picked_chroms[k][i].getIntensity() > max_int;
914  l_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
915  r_tmp = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
916  }
917  }
918  }
919 
920  if (l_tmp > 0.0) left_borders.push_back(l_tmp);
921  if (r_tmp > 0.0) right_borders.push_back(r_tmp);
922 
923  if (pfound == 0) missing_peaks++;
924  if (pfound > 1) multiple_peaks++;
925  }
926 
927  // Check how many chromatograms had exactly one peak picked between our
928  // current left/right borders -> this would be a sign of consistency.
929  OPENMS_LOG_DEBUG << " Overall found missing : " << missing_peaks << " and multiple : " << multiple_peaks << std::endl;
930 
932 
933  // Is there one transitions that is very different from the rest (e.g.
934  // the same element has a bad shape and a bad coelution score) -> potential outlier
935  if (min_index_shape == max_index_coel)
936  {
937  OPENMS_LOG_DEBUG << " Element " << min_index_shape << " is a candidate for removal ... " << std::endl;
938  outlier = String(picked_chroms[min_index_shape].getNativeID());
939  }
940  else
941  {
942  outlier = "none";
943  }
944 
945  // For the final score (larger is better), consider these scores:
946  // - missing_peaks (the more peaks are missing, the worse)
947  // - multiple_peaks
948  // - mean of the shapes (1 is very good, 0 is bad)
949  // - mean of the co-elution scores (0 is good, 1 is ok, above 1 is pretty bad)
950  double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
951  double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
952  coel_score = (coel_score - 1.0) / 2.0;
953 
954  double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
955 
956  OPENMS_LOG_DEBUG << " Computed score " << score << " (from " << shape_score <<
957  " - " << coel_score << " - " << 1.0 * missing_peaks / picked_chroms.size() << ")" << std::endl;
958 
959  return score;
960  }
961 
971  template <typename SpectrumT>
972  void recalculatePeakBorders_(const std::vector<SpectrumT>& picked_chroms, double& best_left, double& best_right, double max_z)
973  {
974  // 1. Collect all seeds that lie within the current seed
975  // - Per chromatogram only the most intense one counts, otherwise very
976  // - low intense peaks can contribute disproportionally to the voting
977  // - procedure.
978  // TODO Especially with DDA MS1 "transitions" from FFID, you often have exactly the same
979  // borders and the logs get confusing and you get -Nan CVs. You also might save time if you use a set here.
980  std::vector<double> left_borders;
981  std::vector<double> right_borders;
982  for (Size k = 0; k < picked_chroms.size(); k++)
983  {
984  double max_int = -1;
985  double left = -1;
986  double right = -1;
987  for (Size i = 0; i < picked_chroms[k].size(); i++)
988  {
989  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
990  {
991  if (picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i] > max_int)
992  {
993  max_int = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_ABUNDANCE][i];
994  left = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_LEFTBORDER][i];
995  right = picked_chroms[k].getFloatDataArrays()[PeakPickerMRM::IDX_RIGHTBORDER][i];
996  }
997  }
998  }
999  if (max_int > -1 )
1000  {
1001  left_borders.push_back(left);
1002  right_borders.push_back(right);
1003  OPENMS_LOG_DEBUG << " * peak " << k << " left boundary " << left_borders.back() << " with inty " << max_int << std::endl;
1004  OPENMS_LOG_DEBUG << " * peak " << k << " right boundary " << right_borders.back() << " with inty " << max_int << std::endl;
1005  }
1006  }
1007 
1008  // Return for empty peak list
1009  if (right_borders.empty())
1010  {
1011  return;
1012  }
1013 
1014  // FEATURE IDEA: instead of Z-score use modified Z-score for small data sets
1015  // http://d-scholarship.pitt.edu/7948/1/Seo.pdf
1016  // http://www.itl.nist.gov/div898/handbook/eda/section3/eda35h.htm
1017  // 1. calculate median
1018  // 2. MAD = calculate difference to median for each value -> take median of that
1019  // 3. Mi = 0.6745*(xi - median) / MAD
1020 
1021  // 2. Calculate mean and standard deviation
1022  // If the coefficient of variation is too large for one border, we use a
1023  // "pseudo-median" instead of the border of the most intense peak.
1024  double mean, stdev;
1025 
1026  // Right borders
1027  mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
1028  stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
1029  / right_borders.size() - mean * mean);
1030  std::sort(right_borders.begin(), right_borders.end());
1031 
1032  OPENMS_LOG_DEBUG << " - Recalculating right peak boundaries " << mean << " mean / best "
1033  << best_right << " std " << stdev << " : " << std::fabs(best_right - mean) / stdev
1034  << " coefficient of variation" << std::endl;
1035 
1036  // Compare right borders of best transition with the mean
1037  if (std::fabs(best_right - mean) / stdev > max_z)
1038  {
1039  best_right = right_borders[right_borders.size() / 2]; // pseudo median
1040  OPENMS_LOG_DEBUG << " - CV too large: correcting right boundary to " << best_right << std::endl;
1041  }
1042 
1043  // Left borders
1044  mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
1045  stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
1046  / left_borders.size() - mean * mean);
1047  std::sort(left_borders.begin(), left_borders.end());
1048 
1049  OPENMS_LOG_DEBUG << " - Recalculating left peak boundaries " << mean << " mean / best "
1050  << best_left << " std " << stdev << " : " << std::fabs(best_left - mean) / stdev
1051  << " coefficient of variation" << std::endl;
1052 
1053  // Compare left borders of best transition with the mean
1054  if (std::fabs(best_left - mean) / stdev > max_z)
1055  {
1056  best_left = left_borders[left_borders.size() / 2]; // pseudo median
1057  OPENMS_LOG_DEBUG << " - CV too large: correcting left boundary to " << best_left << std::endl;
1058  }
1059 
1060  }
1061 
1063 
1064 
1079  template <typename SpectrumT>
1080  void prepareMasterContainer_(const SpectrumT& ref_chromatogram,
1081  SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1082  {
1083  OPENMS_PRECONDITION(master_peak_container.empty(), "Master peak container must be empty")
1084 
1085  // get the start / end point of this chromatogram => then add one more
1086  // point beyond the two boundaries to make the resampling accurate also
1087  // at the edge.
1088  typename SpectrumT::const_iterator begin = ref_chromatogram.begin();
1089  while (begin != ref_chromatogram.end() && begin->getMZ() < left_boundary) {begin++; }
1090  if (begin != ref_chromatogram.begin()) {begin--; }
1091 
1092  typename SpectrumT::const_iterator end = begin;
1093  while (end != ref_chromatogram.end() && end->getMZ() < right_boundary) {end++; }
1094  if (end != ref_chromatogram.end()) {end++; }
1095 
1096  // resize the master container and set the m/z values to the ones of the master container
1097  master_peak_container.resize(distance(begin, end)); // initialize to zero
1098  typename SpectrumT::iterator it = master_peak_container.begin();
1099  for (typename SpectrumT::const_iterator chrom_it = begin; chrom_it != end; chrom_it++, it++)
1100  {
1101  it->setMZ(chrom_it->getMZ());
1102  }
1103  }
1104 
1115  template <typename SpectrumT>
1116  SpectrumT resampleChromatogram_(const SpectrumT& chromatogram,
1117  const SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1118  {
1119  // get the start / end point of this chromatogram => then add one more
1120  // point beyond the two boundaries to make the resampling accurate also
1121  // at the edge.
1122  typename SpectrumT::const_iterator begin = chromatogram.begin();
1123  while (begin != chromatogram.end() && begin->getMZ() < left_boundary) {begin++;}
1124  if (begin != chromatogram.begin()) {begin--;}
1125 
1126  typename SpectrumT::const_iterator end = begin;
1127  while (end != chromatogram.end() && end->getMZ() < right_boundary) {end++;}
1128  if (end != chromatogram.end()) {end++;}
1129 
1130  SpectrumT resampled_peak_container = master_peak_container; // copy the master container, which contains the RT values
1131  LinearResamplerAlign lresampler;
1132  lresampler.raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1133 
1134  return resampled_peak_container;
1135  }
1136 
1138 
1139  // Members
1148  double min_qual_;
1149 
1155 
1162 
1165  };
1166 }
1167 
1168 
#define OPENMS_LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:470
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:93
A method or algorithm argument contains illegal values.
Definition: Exception.h:656
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:1161
String peak_integration_
Definition: MRMTransitionGroupPicker.h:1140
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1147
PeakPickerMRM picker_
Definition: MRMTransitionGroupPicker.h:1163
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:788
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:1080
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1141
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1145
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1143
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:1152
void remove_overlapping_features(std::vector< SpectrumT > &picked_chroms, double best_left, double best_right)
Remove overlapping features.
Definition: MRMTransitionGroupPicker.h:725
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:972
MRMTransitionGroupPicker & operator=(const MRMTransitionGroupPicker &rhs)
Assignment operator is protected for algorithm.
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1150
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:1116
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1153
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1164
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:1148
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1154
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:821
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1151
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:1144
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:1142
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:593
bool compute_peak_shape_metrics_
Definition: MRMTransitionGroupPicker.h:1146
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:58
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 DataValue &default_value=DataValue::EMPTY) const
Returns the value corresponding to a string, or a default value (default: DataValue::EMPTY) if not fo...
void setMZ(CoordinateType coordinate)
Mutable access to the m/z coordinate (index 1)
Definition: Peak2D.h:202
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:214
IntensityType getIntensity() const
Definition: Peak2D.h:166
void setIntensity(IntensityType intensity)
Non-mutable access to the data point intensity (height)
Definition: Peak2D.h:172
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:61
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
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:133
const double k
Definition: Constants.h:153
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
OPENSWATHALGO_DLLAPI double rankedMutualInformation(std::vector< double > &data1, std::vector< double > &data2)
OPENSWATHALGO_DLLAPI XCorrArrayType normalizedCrossCorrelation(std::vector< double > &data1, std::vector< double > &data2, const int &maxdelay, const int &lag)
OPENSWATHALGO_DLLAPI XCorrArrayType::const_iterator xcorrArrayGetMaxPeak(const XCorrArrayType &array)
Find best peak in an cross-correlation (highest apex)
Definition: Scoring.h:62