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  String native_id = chromatogram.getNativeID();
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) break;
178 
179  // Compute a feature from the individual chromatograms and add non-zero features
180  MRMFeature mrm_feature = createMRMFeature(transition_group, picked_chroms_, smoothed_chroms_, chr_idx, peak_idx);
181  if (mrm_feature.getIntensity() > 0)
182  {
183  features.push_back(mrm_feature);
184  }
185 
186  cnt++;
187  if (stop_after_feature_ > 0 && cnt > stop_after_feature_) {break;}
188  if (mrm_feature.getIntensity() > 0 &&
189  mrm_feature.getIntensity() / (double)mrm_feature.getMetaValue("total_xic") < stop_after_intensity_ratio_)
190  {
191  break;
192  }
193  }
194 
195  // Check for completely overlapping features
196  for (Size i = 0; i < features.size(); i++)
197  {
198  MRMFeature& mrm_feature = features[i];
199  bool skip = false;
200  for (Size j = 0; j < i; j++)
201  {
202  if ((double)mrm_feature.getMetaValue("leftWidth") >= (double)features[j].getMetaValue("leftWidth") &&
203  (double)mrm_feature.getMetaValue("rightWidth") <= (double)features[j].getMetaValue("rightWidth") )
204  { skip = true; }
205  }
206  if (mrm_feature.getIntensity() > 0 && !skip)
207  {
208  transition_group.addFeature(mrm_feature);
209  }
210  }
211 
212  }
213 
215  template <typename SpectrumT, typename TransitionT>
217  std::vector<SpectrumT>& picked_chroms,
218  const std::vector<SpectrumT>& smoothed_chroms,
219  const int chr_idx,
220  const int peak_idx)
221  {
222  OPENMS_PRECONDITION(transition_group.isInternallyConsistent(), "Consistent state required")
223  OPENMS_PRECONDITION(transition_group.chromatogramIdsMatch(), "Chromatogram native IDs need to match keys in transition group")
224 
225  MRMFeature mrmFeature;
226  mrmFeature.setIntensity(0.0);
227  double best_left = picked_chroms[chr_idx].getFloatDataArrays()[1][peak_idx];
228  double best_right = picked_chroms[chr_idx].getFloatDataArrays()[2][peak_idx];
229  double peak_apex = picked_chroms[chr_idx][peak_idx].getRT();
230  LOG_DEBUG << "**** Creating MRMFeature for peak " << chr_idx << " " << peak_idx << " " <<
231  picked_chroms[chr_idx][peak_idx] << " with borders " << best_left << " " <<
232  best_right << " (" << best_right - best_left << ")" << std::endl;
233 
234  if (use_consensus_ && recalculate_peaks_)
235  {
236  // This may change best_left / best_right
237  recalculatePeakBorders_(picked_chroms, best_left, best_right, recalculate_peaks_max_z_);
238  if (peak_apex < best_left || peak_apex > best_right)
239  {
240  // apex fell out of range, lets correct it
241  peak_apex = (best_left + best_right) / 2.0;
242  }
243  }
244 
245  std::vector< double > left_edges;
246  std::vector< double > right_edges;
247  double min_left = best_left;
248  double max_right = best_right;
249  if (use_consensus_)
250  {
251  // Remove other, overlapping, picked peaks (in this and other
252  // chromatograms) and then ensure that at least one peak is set to zero
253  // (the currently best peak).
254  remove_overlapping_features(picked_chroms, best_left, best_right);
255  }
256  else
257  {
258  // Pick the peak with the closest apex to the consensus apex for each chromatogram.
259  // Use the closest peak for the current peak. Note that we will only set the closest peak
260  // per chromatogram to zero, so if there are two peaks for some transitions, we will get
261  // to them later. If there is no peak, then we transfer transition boundaries from "master" peak.
262  {
263  for (Size k = 0; k < picked_chroms.size(); k++)
264  {
265  double peak_apex_dist_min = 1e6;
266  int min_dist = -1;
267  for (Size i = 0; i < picked_chroms[k].size(); i++)
268  {
269  PeakIntegrator::PeakArea pa_tmp = pi_.integratePeak( // get the peak apex
270  picked_chroms[k], picked_chroms[k].getFloatDataArrays()[1][i], picked_chroms[k].getFloatDataArrays()[2][i]);
271  if (pa_tmp.apex_pos > 0.0 && std::fabs(pa_tmp.apex_pos - peak_apex) < peak_apex_dist_min)
272  {
273  min_dist = (int)i;
274  }
275  }
276 
277  // Select master peak boundaries, or in the case we found at least one peak, the local peak boundaries
278  double l = best_left;
279  double r = best_right;
280  if (min_dist >= 0)
281  {
282  l = picked_chroms[k].getFloatDataArrays()[1][min_dist];
283  r = picked_chroms[k].getFloatDataArrays()[2][min_dist];
284  picked_chroms[k][min_dist].setIntensity(0.0); // only remove one peak per transition
285  }
286 
287  left_edges.push_back(l);
288  right_edges.push_back(r);
289  // ensure we remember the overall maxima / minima
290  if (l < min_left) {min_left = l;}
291  if (r > max_right) {max_right = r;}
292  }
293  }
294  }
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_ && min_peak_width_ > 0.0 && std::fabs(best_right - best_left) < min_peak_width_)
299  {
300  return mrmFeature;
301  }
302 
303  if (use_consensus_ && compute_peak_quality_)
304  {
305  String outlier = "none";
306  double qual = computeQuality_(transition_group, picked_chroms, chr_idx, best_left, best_right, outlier);
307  if (qual < min_qual_)
308  {
309  return mrmFeature;
310  }
311  mrmFeature.setMetaValue("potentialOutlier", outlier);
312  mrmFeature.setMetaValue("initialPeakQuality", qual);
313  mrmFeature.setOverallQuality(qual);
314  }
315 
316  // Prepare linear resampling of all the chromatograms, here creating the
317  // empty master_peak_container with the same RT (m/z) values as the
318  // reference chromatogram. We use the overall minimal left boundary and
319  // maximal right boundary to prepare the container.
320  SpectrumT master_peak_container;
321  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
322  prepareMasterContainer_(ref_chromatogram, master_peak_container, min_left, max_right);
323 
324  // Iterate over initial transitions / chromatograms (note that we may
325  // have a different number of picked chromatograms than total transitions
326  // as not all are detecting transitions).
327  double total_intensity = 0; double total_peak_apices = 0; double total_xic = 0; double total_mi = 0;
328  for (Size k = 0; k < transition_group.getTransitions().size(); k++)
329  {
330 
331  double local_left = best_left;
332  double local_right = best_right;
333  if (!use_consensus_)
334  {
335  local_left = left_edges[k];
336  local_right = right_edges[k];
337  }
338 
339  const SpectrumT& chromatogram = selectChromHelper_(transition_group, transition_group.getTransitions()[k].getNativeID());
340  if (transition_group.getTransitions()[k].isDetectingTransition())
341  {
342  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
343  {
344  total_xic += it->getIntensity();
345  }
346  }
347 
348  // Compute total intensity on transition-level
349  double transition_total_xic = 0;
350 
351  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
352  {
353  transition_total_xic += it->getIntensity();
354  }
355 
356  // Compute total mutual information on transition-level.
357  double transition_total_mi = 0;
358  if (compute_total_mi_)
359  {
360  std::vector<double> chrom_vect_id, chrom_vect_det;
361  for (typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
362  {
363  chrom_vect_id.push_back(it->getIntensity());
364  }
365 
366  // compute baseline mutual information
367  int transition_total_mi_norm = 0;
368  for (Size m = 0; m < transition_group.getTransitions().size(); m++)
369  {
370  if (transition_group.getTransitions()[m].isDetectingTransition())
371  {
372  const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.getTransitions()[m].getNativeID());
373  chrom_vect_det.clear();
374  for (typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
375  {
376  chrom_vect_det.push_back(it->getIntensity());
377  }
378  transition_total_mi += OpenSwath::Scoring::rankedMutualInformation(chrom_vect_det, chrom_vect_id);
379  transition_total_mi_norm++;
380  }
381  }
382  if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
383 
384  if (transition_group.getTransitions()[k].isDetectingTransition())
385  {
386  // sum up all transition-level total MI and divide by the number of detection transitions to have peak group level total MI
387  total_mi += transition_total_mi / transition_total_mi_norm;
388  }
389  }
390 
391  SpectrumT used_chromatogram;
392  // resample the current chromatogram
393  if (peak_integration_ == "original")
394  {
395  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
396  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
397  }
398  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= k)
399  {
400  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
401  "Tried to calculate peak area and height without any smoothed chromatograms");
402  }
403  else if (peak_integration_ == "smoothed")
404  {
405  used_chromatogram = resampleChromatogram_(smoothed_chroms[k], master_peak_container, local_left, local_right);
406  }
407  else
408  {
409  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
410  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
411  }
412 
413  Feature f;
414  double quality = 0;
415  f.setQuality(0, quality);
416  f.setOverallQuality(quality);
417 
418  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
419  double peak_integral = pa.area;
420  double peak_apex_int = pa.height;
421  f.setMetaValue("peak_apex_position", pa.apex_pos);
422  if (background_subtraction_ != "none")
423  {
424  double background{0};
425  double avg_noise_level{0};
426  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= k)
427  {
428  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
429  "Tried to calculate background estimation without any smoothed chromatograms");
430  }
431  else if (background_subtraction_ == "original")
432  {
433  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
434  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
435  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
436  avg_noise_level = (intensity_right + intensity_left) / 2;
437  background = avg_noise_level * n_points;
438  }
439  else if (background_subtraction_ == "exact")
440  {
441  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
442  background = pb.area;
443  avg_noise_level = pb.height;
444  }
445  peak_integral -= background;
446  peak_apex_int -= avg_noise_level;
447  if (peak_integral < 0) {peak_integral = 0;}
448  if (peak_apex_int < 0) {peak_apex_int = 0;}
449 
450  f.setMetaValue("area_background_level", background);
451  f.setMetaValue("noise_background_level", avg_noise_level);
452  }
453 
454  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
455  f.setIntensity(peak_integral);
456  ConvexHull2D hull;
457  hull.setHullPoints(pa.hull_points);
458  f.getConvexHulls().push_back(hull);
459  if (chromatogram.metaValueExists("product_mz"))
460  {
461  f.setMetaValue("MZ", chromatogram.getMetaValue("product_mz"));
462  f.setMZ(chromatogram.getMetaValue("product_mz"));
463  }
464  else
465  {
466  LOG_WARN << "Please set meta value 'product_mz' on chromatogram to populate feature m/z value" << std::endl;
467  }
468  f.setMetaValue("native_id", chromatogram.getNativeID());
469  f.setMetaValue("peak_apex_int", peak_apex_int);
470  f.setMetaValue("total_xic", transition_total_xic);
471  if (compute_total_mi_)
472  {
473  f.setMetaValue("total_mi", transition_total_mi);
474  }
475 
476  if (transition_group.getTransitions()[k].isDetectingTransition())
477  {
478  total_intensity += peak_integral;
479  total_peak_apices += peak_apex_int;
480  }
481 
482  // for backwards compatibility with TOPP tests
483  // Calculate peak shape metrics that will be used for later QC
484  if (compute_peak_shape_metrics_)
485  {
486  PeakIntegrator::PeakShapeMetrics psm = pi_.calculatePeakShapeMetrics(used_chromatogram, local_left, local_right, peak_apex_int, pa.apex_pos);
487  f.setMetaValue("width_at_5", psm.width_at_5);
488  f.setMetaValue("width_at_10", psm.width_at_10);
489  f.setMetaValue("width_at_50", psm.width_at_50);
490  f.setMetaValue("start_position_at_5", psm.start_position_at_5);
491  f.setMetaValue("start_position_at_10", psm.start_position_at_10);
492  f.setMetaValue("start_position_at_50", psm.start_position_at_50);
493  f.setMetaValue("end_position_at_5", psm.end_position_at_5);
494  f.setMetaValue("end_position_at_10", psm.end_position_at_10);
495  f.setMetaValue("end_position_at_50", psm.end_position_at_50);
496  f.setMetaValue("total_width", psm.total_width);
497  f.setMetaValue("tailing_factor", psm.tailing_factor);
498  f.setMetaValue("asymmetry_factor", psm.asymmetry_factor);
499  f.setMetaValue("slope_of_baseline", psm.slope_of_baseline);
500  f.setMetaValue("baseline_delta_2_height", psm.baseline_delta_2_height);
501  f.setMetaValue("points_across_baseline", psm.points_across_baseline);
502  f.setMetaValue("points_across_half_height", psm.points_across_half_height);
503  }
504 
505  mrmFeature.addFeature(f, chromatogram.getNativeID()); //map index and feature
506  }
507 
508  // Also pick the precursor chromatogram(s); note total_xic is not
509  // extracted here, only for fragment traces
510  for (Size k = 0; k < transition_group.getPrecursorChromatograms().size(); k++)
511  {
512  const SpectrumT& chromatogram = transition_group.getPrecursorChromatograms()[k];
513  Size prec_idx = transition_group.getChromatograms().size() + k;
514 
515  double local_left = best_left;
516  double local_right = best_right;
517  if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
518  {
519  local_left = left_edges[prec_idx];
520  local_right = right_edges[prec_idx];
521  }
522 
523  SpectrumT used_chromatogram;
524  // resample the current chromatogram
525  if (peak_integration_ == "original")
526  {
527  used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
528  // const SpectrumT& used_chromatogram = chromatogram; // instead of resampling
529  }
530  else if (peak_integration_ == "smoothed" && smoothed_chroms.size() <= prec_idx)
531  {
532  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
533  "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
534  }
535  else if (peak_integration_ == "smoothed")
536  {
537  used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
538  }
539  else
540  {
541  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
542  String("Peak integration chromatogram ") + peak_integration_ + " is not a valid method for MRMTransitionGroupPicker");
543  }
544 
545  Feature f;
546  double quality = 0;
547  f.setQuality(0, quality);
548  f.setOverallQuality(quality);
549 
550  PeakIntegrator::PeakArea pa = pi_.integratePeak(used_chromatogram, local_left, local_right);
551  double peak_integral = pa.area;
552  double peak_apex_int = pa.height;
553 
554  if (background_subtraction_ != "none")
555  {
556  double background{0};
557  double avg_noise_level{0};
558  if ((peak_integration_ == "smoothed") && smoothed_chroms.size() <= prec_idx)
559  {
560  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
561  "Tried to calculate background estimation without any smoothed chromatograms");
562  }
563  else if (background_subtraction_ == "original")
564  {
565  const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
566  const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
567  const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
568  avg_noise_level = (intensity_right + intensity_left) / 2;
569  background = avg_noise_level * n_points;
570  }
571  else if (background_subtraction_ == "exact")
572  {
573  PeakIntegrator::PeakBackground pb = pi_.estimateBackground(used_chromatogram, local_left, local_right, pa.apex_pos);
574  background = pb.area;
575  avg_noise_level = pb.height;
576  }
577  peak_integral -= background;
578  peak_apex_int -= avg_noise_level;
579  if (peak_integral < 0) {peak_integral = 0;}
580  if (peak_apex_int < 0) {peak_apex_int = 0;}
581 
582  f.setMetaValue("area_background_level", background);
583  f.setMetaValue("noise_background_level", avg_noise_level);
584  }
585 
586  if (chromatogram.metaValueExists("precursor_mz"))
587  {
588  f.setMZ(chromatogram.getMetaValue("precursor_mz"));
589  mrmFeature.setMZ(chromatogram.getMetaValue("precursor_mz"));
590  }
591 
592  f.setRT(picked_chroms[chr_idx][peak_idx].getMZ());
593  f.setIntensity(peak_integral);
594  ConvexHull2D hull;
595  hull.setHullPoints(pa.hull_points);
596  f.getConvexHulls().push_back(hull);
597  f.setMetaValue("native_id", chromatogram.getNativeID());
598  f.setMetaValue("peak_apex_int", peak_apex_int);
599 
600  if (use_precursors_ && transition_group.getTransitions().empty())
601  {
602  total_intensity += peak_integral;
603  }
604 
605  mrmFeature.addPrecursorFeature(f, chromatogram.getNativeID());
606  }
607 
608  mrmFeature.setRT(peak_apex);
609  mrmFeature.setIntensity(total_intensity);
610  mrmFeature.setMetaValue("PeptideRef", transition_group.getTransitionGroupID());
611  mrmFeature.setMetaValue("leftWidth", best_left);
612  mrmFeature.setMetaValue("rightWidth", best_right);
613  mrmFeature.setMetaValue("total_xic", total_xic);
614  if (compute_total_mi_)
615  {
616  mrmFeature.setMetaValue("total_mi", total_mi);
617  }
618  mrmFeature.setMetaValue("peak_apices_sum", total_peak_apices);
619 
620  mrmFeature.ensureUniqueId();
621  return mrmFeature;
622  }
623 
624  // maybe private, but we have tests
625 
637  template <typename SpectrumT>
638  void remove_overlapping_features(std::vector<SpectrumT>& picked_chroms, double best_left, double best_right)
639  {
640  // delete all seeds that lie within the current seed
641  //std::cout << "Removing features for peak between " << best_left << " " << best_right << std::endl;
642  for (Size k = 0; k < picked_chroms.size(); k++)
643  {
644  for (Size i = 0; i < picked_chroms[k].size(); i++)
645  {
646  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
647  {
648  //std::cout << "For Chrom " << k << " removing peak " << picked_chroms[k][i].getMZ() << " l/r : " << picked_chroms[k].getFloatDataArrays()[1][i] << " " <<
649  // picked_chroms[k].getFloatDataArrays()[2][i] << " with int " << picked_chroms[k][i].getIntensity() <<std::endl;
650  picked_chroms[k][i].setIntensity(0.0);
651  }
652  }
653  }
654 
655  // delete all seeds that overlap within the current seed
656  for (Size k = 0; k < picked_chroms.size(); k++)
657  {
658  for (Size i = 0; i < picked_chroms[k].size(); i++)
659  {
660  if (picked_chroms[k][i].getIntensity() <= 0.0) {continue; }
661 
662  double left = picked_chroms[k].getFloatDataArrays()[1][i];
663  double right = picked_chroms[k].getFloatDataArrays()[2][i];
664  if ((left > best_left && left < best_right)
665  || (right > best_left && right < best_right))
666  {
667  //std::cout << "= For Chrom " << k << " removing contained peak " << picked_chroms[k][i].getMZ() << " l/r : " << picked_chroms[k].getFloatDataArrays()[1][i] << " " <<
668  // picked_chroms[k].getFloatDataArrays()[2][i] << " with int " << picked_chroms[k][i].getIntensity() <<std::endl;
669  picked_chroms[k][i].setIntensity(0.0);
670  }
671  }
672  }
673  }
674 
676  void findLargestPeak(const std::vector<MSChromatogram >& picked_chroms, int& chr_idx, int& peak_idx);
677 
686  void findWidestPeakIndices(const std::vector<MSChromatogram>& picked_chroms, Int& chrom_idx, Int& point_idx) const;
687 
688 protected:
689 
691  void updateMembers_() override;
692 
694  MRMTransitionGroupPicker& operator=(const MRMTransitionGroupPicker& rhs);
695 
699  template <typename SpectrumT, typename TransitionT>
700  const SpectrumT& selectChromHelper_(const MRMTransitionGroup<SpectrumT, TransitionT>& transition_group, const String& native_id)
701  {
702  if (transition_group.hasChromatogram(native_id))
703  {
704  return transition_group.getChromatogram(native_id);
705  }
706  else if (transition_group.hasPrecursorChromatogram(native_id))
707  {
708  return transition_group.getPrecursorChromatogram(native_id);
709  }
710  else
711  {
712  throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION, "Did not find chromatogram for id '" + native_id + "'.");
713  }
714  }
715 
732  template <typename SpectrumT, typename TransitionT>
734  const std::vector<SpectrumT>& picked_chroms,
735  const int chr_idx,
736  const double best_left,
737  const double best_right,
738  String& outlier)
739  {
740  // Resample all chromatograms around the current estimated peak and
741  // collect the raw intensities. For resampling, use a bit more on either
742  // side to correctly identify shoulders etc.
743  double resample_boundary = resample_boundary_; // sample 15 seconds more on each side
744  SpectrumT master_peak_container;
745  const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
746  prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
747  std::vector<std::vector<double> > all_ints;
748  for (Size k = 0; k < picked_chroms.size(); k++)
749  {
750  const SpectrumT chromatogram = selectChromHelper_(transition_group, picked_chroms[k].getNativeID());
751  const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
752  master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
753 
754  std::vector<double> int_here;
755  for (Size i = 0; i < used_chromatogram.size(); i++)
756  {
757  int_here.push_back(used_chromatogram[i].getIntensity());
758  }
759  all_ints.push_back(int_here);
760  }
761 
762  // Compute the cross-correlation for the collected intensities
763  std::vector<double> mean_shapes;
764  std::vector<double> mean_coel;
765  for (Size k = 0; k < all_ints.size(); k++)
766  {
767  std::vector<double> shapes;
768  std::vector<double> coel;
769  for (Size i = 0; i < all_ints.size(); i++)
770  {
771  if (i == k) {continue;}
773  all_ints[k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
774 
775  // the first value is the x-axis (retention time) and should be an int -> it show the lag between the two
776  double res_coelution = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->first);
777  double res_shape = std::abs(OpenSwath::Scoring::xcorrArrayGetMaxPeak(res)->second);
778 
779  shapes.push_back(res_shape);
780  coel.push_back(res_coelution);
781  }
782 
783  // We have computed the cross-correlation of chromatogram k against
784  // all others. Use the mean of these computations as the value for k.
786  msc = std::for_each(shapes.begin(), shapes.end(), msc);
787  double shapes_mean = msc.mean();
788  msc = std::for_each(coel.begin(), coel.end(), msc);
789  double coel_mean = msc.mean();
790 
791  // mean shape scores below 0.5-0.6 should be a real sign of trouble ... !
792  // mean coel scores above 3.5 should be a real sign of trouble ... !
793  mean_shapes.push_back(shapes_mean);
794  mean_coel.push_back(coel_mean);
795  }
796 
797  // find the chromatogram with the minimal shape score and the maximal
798  // coelution score -> if it is the same chromatogram, the chance is
799  // pretty good that it is different from the others...
800  int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
801  int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
802 
803  // Look at the picked peaks that are within the current left/right borders
804  int missing_peaks = 0;
805  int multiple_peaks = 0;
806 
807  // collect all seeds that lie within the current seed
808  std::vector<double> left_borders;
809  std::vector<double> right_borders;
810  for (Size k = 0; k < picked_chroms.size(); k++)
811  {
812  double l_tmp;
813  double r_tmp;
814  double max_int = -1;
815 
816  int pfound = 0;
817  l_tmp = -1;
818  r_tmp = -1;
819  for (Size i = 0; i < picked_chroms[k].size(); i++)
820  {
821  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
822  {
823  pfound++;
824  if (picked_chroms[k][i].getIntensity() > max_int)
825  {
826  max_int = picked_chroms[k][i].getIntensity() > max_int;
827  l_tmp = picked_chroms[k].getFloatDataArrays()[1][i];
828  r_tmp = picked_chroms[k].getFloatDataArrays()[2][i];
829  }
830  }
831  }
832 
833  if (l_tmp > 0.0) left_borders.push_back(l_tmp);
834  if (r_tmp > 0.0) right_borders.push_back(r_tmp);
835 
836  if (pfound == 0) missing_peaks++;
837  if (pfound > 1) multiple_peaks++;
838  }
839 
840  // Check how many chromatograms had exactly one peak picked between our
841  // current left/right borders -> this would be a sign of consistency.
842  LOG_DEBUG << " Overall found missing : " << missing_peaks << " and multiple : " << multiple_peaks << std::endl;
843 
845 
846  // Is there one transitions that is very different from the rest (e.g.
847  // the same element has a bad shape and a bad coelution score) -> potential outlier
848  if (min_index_shape == max_index_coel)
849  {
850  LOG_DEBUG << " element " << min_index_shape << " is a candidate for removal ... " << std::endl;
851  outlier = String(picked_chroms[min_index_shape].getNativeID());
852  }
853  else
854  {
855  outlier = "none";
856  }
857 
858  // For the final score (larger is better), consider these scores:
859  // - missing_peaks (the more peaks are missing, the worse)
860  // - multiple_peaks
861  // - mean of the shapes (1 is very good, 0 is bad)
862  // - mean of the co-elution scores (0 is good, 1 is ok, above 1 is pretty bad)
863  double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
864  double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
865  coel_score = (coel_score - 1.0) / 2.0;
866 
867  double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
868 
869  LOG_DEBUG << " computed score " << score << " (from " << shape_score <<
870  " - " << coel_score << " - " << 1.0 * missing_peaks / picked_chroms.size() << ")" << std::endl;
871 
872  return score;
873  }
874 
884  template <typename SpectrumT>
885  void recalculatePeakBorders_(const std::vector<SpectrumT>& picked_chroms, double& best_left, double& best_right, double max_z)
886  {
887  // 1. Collect all seeds that lie within the current seed
888  // - Per chromatogram only the most intense one counts, otherwise very
889  // - low intense peaks can contribute disproportionally to the voting
890  // - procedure.
891  std::vector<double> left_borders;
892  std::vector<double> right_borders;
893  for (Size k = 0; k < picked_chroms.size(); k++)
894  {
895  double max_int = -1;
896  double left = -1;
897  double right = -1;
898  for (Size i = 0; i < picked_chroms[k].size(); i++)
899  {
900  if (picked_chroms[k][i].getMZ() >= best_left && picked_chroms[k][i].getMZ() <= best_right)
901  {
902  if (picked_chroms[k].getFloatDataArrays()[0][i] > max_int)
903  {
904  max_int = picked_chroms[k].getFloatDataArrays()[0][i];
905  left = picked_chroms[k].getFloatDataArrays()[1][i];
906  right = picked_chroms[k].getFloatDataArrays()[2][i];
907  }
908  }
909  }
910  if (max_int > -1 )
911  {
912  left_borders.push_back(left);
913  right_borders.push_back(right);
914  LOG_DEBUG << " * " << k << " left boundary " << left_borders.back() << " with int " << max_int << std::endl;
915  LOG_DEBUG << " * " << k << " right boundary " << right_borders.back() << " with int " << max_int << std::endl;
916  }
917  }
918 
919  // Return for empty peak list
920  if (right_borders.empty())
921  {
922  return;
923  }
924 
925  // FEATURE IDEA: instead of Z-score use modified Z-score for small data sets
926  // http://d-scholarship.pitt.edu/7948/1/Seo.pdf
927  // http://www.itl.nist.gov/div898/handbook/eda/section3/eda35h.htm
928  // 1. calculate median
929  // 2. MAD = calculate difference to median for each value -> take median of that
930  // 3. Mi = 0.6745*(xi - median) / MAD
931 
932  // 2. Calculate mean and standard deviation
933  // If the coefficient of variation is too large for one border, we use a
934  // "pseudo-median" instead of the border of the most intense peak.
935  double mean, stdev;
936 
937  // Right borders
938  mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
939  stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
940  / right_borders.size() - mean * mean);
941  std::sort(right_borders.begin(), right_borders.end());
942 
943  LOG_DEBUG << " - Recalculating right peak boundaries " << mean << " mean / best "
944  << best_right << " std " << stdev << " : " << std::fabs(best_right - mean) / stdev
945  << " coefficient of variation" << std::endl;
946 
947  // Compare right borders of best transition with the mean
948  if (std::fabs(best_right - mean) / stdev > max_z)
949  {
950  best_right = right_borders[right_borders.size() / 2]; // pseudo median
951  LOG_DEBUG << " - Setting right boundary to " << best_right << std::endl;
952  }
953 
954  // Left borders
955  mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
956  stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
957  / left_borders.size() - mean * mean);
958  std::sort(left_borders.begin(), left_borders.end());
959 
960  LOG_DEBUG << " - Recalculating left peak boundaries " << mean << " mean / best "
961  << best_left << " std " << stdev << " : " << std::fabs(best_left - mean) / stdev
962  << " coefficient of variation" << std::endl;
963 
964  // Compare left borders of best transition with the mean
965  if (std::fabs(best_left - mean) / stdev > max_z)
966  {
967  best_left = left_borders[left_borders.size() / 2]; // pseudo median
968  LOG_DEBUG << " - Setting left boundary to " << best_left << std::endl;
969  }
970 
971  }
972 
974 
975 
990  template <typename SpectrumT>
991  void prepareMasterContainer_(const SpectrumT& ref_chromatogram,
992  SpectrumT& master_peak_container, double left_boundary, double right_boundary)
993  {
994  OPENMS_PRECONDITION(master_peak_container.empty(), "Master peak container must be empty")
995 
996  // get the start / end point of this chromatogram => then add one more
997  // point beyond the two boundaries to make the resampling accurate also
998  // at the edge.
999  typename SpectrumT::const_iterator begin = ref_chromatogram.begin();
1000  while (begin != ref_chromatogram.end() && begin->getMZ() < left_boundary) {begin++; }
1001  if (begin != ref_chromatogram.begin()) {begin--; }
1002 
1003  typename SpectrumT::const_iterator end = begin;
1004  while (end != ref_chromatogram.end() && end->getMZ() < right_boundary) {end++; }
1005  if (end != ref_chromatogram.end()) {end++; }
1006 
1007  // resize the master container and set the m/z values to the ones of the master container
1008  master_peak_container.resize(distance(begin, end)); // initialize to zero
1009  typename SpectrumT::iterator it = master_peak_container.begin();
1010  for (typename SpectrumT::const_iterator chrom_it = begin; chrom_it != end; chrom_it++, it++)
1011  {
1012  it->setMZ(chrom_it->getMZ());
1013  }
1014  }
1015 
1026  template <typename SpectrumT>
1027  SpectrumT resampleChromatogram_(const SpectrumT& chromatogram,
1028  const SpectrumT& master_peak_container, double left_boundary, double right_boundary)
1029  {
1030  // get the start / end point of this chromatogram => then add one more
1031  // point beyond the two boundaries to make the resampling accurate also
1032  // at the edge.
1033  typename SpectrumT::const_iterator begin = chromatogram.begin();
1034  while (begin != chromatogram.end() && begin->getMZ() < left_boundary) {begin++;}
1035  if (begin != chromatogram.begin()) {begin--;}
1036 
1037  typename SpectrumT::const_iterator end = begin;
1038  while (end != chromatogram.end() && end->getMZ() < right_boundary) {end++;}
1039  if (end != chromatogram.end()) {end++;}
1040 
1041  SpectrumT resampled_peak_container = master_peak_container; // copy the master container, which contains the RT values
1042  LinearResamplerAlign lresampler;
1043  lresampler.raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1044 
1045  return resampled_peak_container;
1046  }
1047 
1049 
1050  // Members
1059  double min_qual_;
1060 
1066 
1073 
1076  };
1077 }
1078 
1079 
const double k
double total_width
Definition: PeakIntegrator.h:165
void setMetaValue(const String &name, const DataValue &value)
Sets the DataValue corresponding to a name.
bool hasTransition(String key) const
Definition: MRMTransitionGroup.h:158
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:700
double min_qual_
Definition: MRMTransitionGroupPicker.h:1059
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1052
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1064
bool use_consensus_
Definition: MRMTransitionGroupPicker.h:1055
A more convenient string class.
Definition: String.h:57
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1056
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:106
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:207
double start_position_at_50
Definition: PeakIntegrator.h:149
The PeakPickerMRM finds peaks a single chromatogram.
Definition: PeakPickerMRM.h:68
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1075
double width_at_5
Definition: PeakIntegrator.h:129
double slope_of_baseline
Definition: PeakIntegrator.h:194
double start_position_at_5
Definition: PeakIntegrator.h:141
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1061
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:1051
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:458
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:134
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
#define LOG_WARN
Macro if a warning, a piece of information which should be read by the user, should be logged...
Definition: LogStream.h:450
Definition: Scoring.h:61
double start_position_at_10
Definition: PeakIntegrator.h:145
void sortByIntensity(bool reverse=false)
Lexicographically sorts the peaks by their intensity.
double min_peak_width_
Definition: MRMTransitionGroupPicker.h:1063
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1065
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:157
bool recalculate_peaks_
Definition: MRMTransitionGroupPicker.h:1053
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:98
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:137
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:1057
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:991
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:203
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:153
double area
Definition: PeakIntegrator.h:86
bool hasChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:192
Definition: PeakIntegrator.h:124
double asymmetry_factor
Definition: PeakIntegrator.h:189
void setOverallQuality(QualityType q)
Set the overall quality.
double end_position_at_50
Definition: PeakIntegrator.h:161
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:1027
double height
Definition: PeakIntegrator.h:115
double width_at_10
Definition: PeakIntegrator.h:133
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:68
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:216
PeakPickerMRM picker_
Definition: MRMTransitionGroupPicker.h:1074
Definition: PeakIntegrator.h:81
double area
Definition: PeakIntegrator.h:111
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1054
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1058
double apex_pos
Definition: PeakIntegrator.h:94
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:885
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:638
String boundary_selection_method_
Which method to use for selecting peaks&#39; boundaries.
Definition: MRMTransitionGroupPicker.h:1072
std::vector< ChromatogramType > & getPrecursorChromatograms()
Definition: MRMTransitionGroup.h:215
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1062
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:733
double height
Definition: PeakIntegrator.h:90
int Int
Signed integer type.
Definition: Types.h:102
double tailing_factor
Definition: PeakIntegrator.h:179
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:199