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