86 template <
typename SpectrumT,
typename TransitionT>
92 std::vector<MSChromatogram > picked_chroms;
93 std::vector<MSChromatogram > smoothed_chroms;
104 !transition_group.
getTransition(native_id).isDetectingTransition() )
111 picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
113 picked_chroms.push_back(std::move(picked_chrom));
114 smoothed_chroms.push_back(std::move(smoothed_chrom));
122 SpectrumT picked_chrom, smoothed_chrom;
125 picker_.pickChromatogram(chromatogram, picked_chrom, smoothed_chrom);
126 picked_chrom.sortByIntensity();
127 picked_chroms.push_back(picked_chrom);
128 smoothed_chroms.push_back(smoothed_chrom);
136 int chr_idx, peak_idx, cnt = 0;
137 std::vector<MRMFeature> features;
140 chr_idx = -1; peak_idx = -1;
142 if (boundary_selection_method_ ==
"largest")
144 findLargestPeak(picked_chroms, chr_idx, peak_idx);
146 else if (boundary_selection_method_ ==
"widest")
148 findWidestPeakIndices(picked_chroms, chr_idx, peak_idx);
151 if (chr_idx == -1 && peak_idx == -1)
153 OPENMS_LOG_DEBUG <<
"**** No more peaks left. Nr. chroms: " << picked_chroms.size() << std::endl;
158 MRMFeature mrm_feature = createMRMFeature(transition_group, picked_chroms, smoothed_chroms, chr_idx, peak_idx);
159 double total_xic = 0;
164 features.push_back(std::move(mrm_feature));
168 if (stop_after_feature_ > 0 && cnt > stop_after_feature_) {
172 if (intensity > 0 && intensity / total_xic < stop_after_intensity_ratio_)
174 OPENMS_LOG_DEBUG <<
"**** Minimum intensity ratio reached. Nr. chroms: " << picked_chroms.size() << std::endl;
180 for (
Size i = 0; i < features.size(); i++)
184 for (
Size j = 0; j < i; j++)
186 if ((
double)mrm_feature.
getMetaValue(
"leftWidth") >= (
double)features[j].getMetaValue(
"leftWidth") &&
187 (
double)mrm_feature.
getMetaValue(
"rightWidth") <= (
double)features[j].getMetaValue(
"rightWidth") )
199 template <
typename SpectrumT,
typename TransitionT>
201 std::vector<SpectrumT>& picked_chroms,
202 const std::vector<SpectrumT>& smoothed_chroms,
213 double peak_apex = picked_chroms[chr_idx][peak_idx].getRT();
214 OPENMS_LOG_DEBUG <<
"**** Creating MRMFeature for peak " << peak_idx <<
" in chrom. " << chr_idx <<
" with " <<
215 picked_chroms[chr_idx][peak_idx] <<
" and borders " << best_left <<
" " <<
216 best_right <<
" (" << best_right - best_left <<
")" << std::endl;
218 if (use_consensus_ && recalculate_peaks_)
221 recalculatePeakBorders_(picked_chroms, best_left, best_right, recalculate_peaks_max_z_);
222 if (peak_apex < best_left || peak_apex > best_right)
225 peak_apex = (best_left + best_right) / 2.0;
229 std::vector< double > left_edges;
230 std::vector< double > right_edges;
231 double min_left = best_left;
232 double max_right = best_right;
238 remove_overlapping_features(picked_chroms, best_left, best_right);
242 pickApex(picked_chroms, best_left, best_right, peak_apex,
243 min_left, max_right, left_edges, right_edges);
246 picked_chroms[chr_idx][peak_idx].setIntensity(0.0);
251 if (min_peak_width_ > 0.0 && std::fabs(best_right - best_left) < min_peak_width_)
256 if (compute_peak_quality_)
259 double qual = computeQuality_(transition_group, picked_chroms, chr_idx, best_left, best_right, outlier);
260 if (qual < min_qual_)
274 SpectrumT master_peak_container;
275 const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
276 prepareMasterContainer_(ref_chromatogram, master_peak_container, min_left, max_right);
281 double total_intensity = 0;
double total_peak_apices = 0;
double total_xic = 0;
double total_mi = 0;
282 pickFragmentChromatograms(transition_group, picked_chroms, mrmFeature, smoothed_chroms,
283 best_left, best_right, use_consensus_,
284 total_intensity, total_xic, total_mi, total_peak_apices,
285 master_peak_container, left_edges, right_edges,
290 pickPrecursorChromatograms(transition_group,
291 picked_chroms, mrmFeature, smoothed_chroms,
292 best_left, best_right, use_consensus_,
293 total_intensity, master_peak_container, left_edges, right_edges,
296 mrmFeature.
setRT(peak_apex);
302 if (compute_total_mi_)
306 mrmFeature.
setMetaValue(
"peak_apices_sum", total_peak_apices);
324 template <
typename SpectrumT>
325 void pickApex(std::vector<SpectrumT>& picked_chroms,
326 const double best_left,
const double best_right,
const double peak_apex,
327 double &min_left,
double &max_right,
328 std::vector< double > & left_edges, std::vector< double > & right_edges)
330 for (
Size k = 0;
k < picked_chroms.size();
k++)
332 double peak_apex_dist_min = std::numeric_limits<double>::max();
334 for (
Size i = 0; i < picked_chroms[
k].size(); i++)
340 if (pa_tmp.
apex_pos > 1e-11 && std::fabs(pa_tmp.
apex_pos - peak_apex) < peak_apex_dist_min)
342 peak_apex_dist_min = std::fabs(pa_tmp.
apex_pos - peak_apex);
348 double l = best_left;
349 double r = best_right;
354 picked_chroms[
k][min_dist].setIntensity(0.0);
357 left_edges.push_back(l);
358 right_edges.push_back(r);
360 if (l < min_left) {min_left = l;}
361 if (r > max_right) {max_right = r;}
365 template <
typename SpectrumT,
typename TransitionT>
367 const std::vector<SpectrumT>& picked_chroms,
369 const std::vector<SpectrumT>& smoothed_chroms,
370 const double best_left,
const double best_right,
371 const bool use_consensus_,
372 double & total_intensity,
375 double & total_peak_apices,
376 const SpectrumT & master_peak_container,
377 const std::vector< double > & left_edges,
378 const std::vector< double > & right_edges,
385 double local_left = best_left;
386 double local_right = best_right;
395 "When using non-consensus peak picker, all transitions need to be detecting transitions.");
397 local_left = left_edges[
k];
398 local_right = right_edges[
k];
401 const SpectrumT& chromatogram = selectChromHelper_(transition_group, transition_group.
getTransitions()[
k].getNativeID());
404 for (
typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
406 total_xic += it->getIntensity();
411 double transition_total_xic = 0;
413 for (
typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
415 transition_total_xic += it->getIntensity();
419 double transition_total_mi = 0;
420 if (compute_total_mi_)
422 std::vector<unsigned int> chrom_vect_id_ranked, chrom_vect_det_ranked;
423 std::vector<double> chrom_vect_id, chrom_vect_det;
424 for (
typename SpectrumT::const_iterator it = chromatogram.begin(); it != chromatogram.end(); it++)
426 chrom_vect_id.push_back(it->getIntensity());
430 int transition_total_mi_norm = 0;
435 const SpectrumT& chromatogram_det = selectChromHelper_(transition_group, transition_group.
getTransitions()[m].getNativeID());
436 chrom_vect_det.clear();
437 for (
typename SpectrumT::const_iterator it = chromatogram_det.begin(); it != chromatogram_det.end(); it++)
439 chrom_vect_det.push_back(it->getIntensity());
443 transition_total_mi_norm++;
446 if (transition_total_mi_norm > 0) { transition_total_mi /= transition_total_mi_norm; }
451 total_mi += transition_total_mi / transition_total_mi_norm;
455 SpectrumT used_chromatogram;
457 if (peak_integration_ ==
"original")
459 used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
461 else if (peak_integration_ ==
"smoothed")
463 if (smoothed_chroms.size() <=
k)
466 "Tried to calculate peak area and height without any smoothed chromatograms");
468 used_chromatogram = resampleChromatogram_(smoothed_chroms[
k], master_peak_container, local_left, local_right);
473 String(
"Peak integration chromatogram ") + peak_integration_ +
" is not a valid method for MRMTransitionGroupPicker");
482 double peak_integral = pa.
area;
483 double peak_apex_int = pa.
height;
485 if (background_subtraction_ !=
"none")
487 double background{0};
488 double avg_noise_level{0};
489 if (background_subtraction_ ==
"original")
491 const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
492 const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
493 const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
494 avg_noise_level = (intensity_right + intensity_left) / 2;
495 background = avg_noise_level * n_points;
497 else if (background_subtraction_ ==
"exact")
500 background = pb.
area;
501 avg_noise_level = pb.
height;
503 peak_integral -= background;
504 peak_apex_int -= avg_noise_level;
505 if (peak_integral < 0) {peak_integral = 0;}
506 if (peak_apex_int < 0) {peak_apex_int = 0;}
509 f.
setMetaValue(
"noise_background_level", avg_noise_level);
512 f.
setRT(picked_chroms[chr_idx][peak_idx].getMZ());
518 f.
setMZ(chromatogram.getProduct().getMZ());
519 mrmFeature.
setMZ(chromatogram.getPrecursor().getMZ());
521 if (chromatogram.metaValueExists(
"product_mz"))
523 f.
setMetaValue(
"MZ", chromatogram.getMetaValue(
"product_mz"));
524 f.
setMZ(chromatogram.getMetaValue(
"product_mz"));
527 f.
setMetaValue(
"native_id", chromatogram.getNativeID());
530 if (compute_total_mi_)
537 total_intensity += peak_integral;
538 total_peak_apices += peak_apex_int;
545 if (compute_peak_shape_metrics_)
564 mrmFeature.
addFeature(f, chromatogram.getNativeID());
568 template <
typename SpectrumT,
typename TransitionT>
570 const std::vector<SpectrumT>& picked_chroms,
572 const std::vector<SpectrumT>& smoothed_chroms,
573 const double best_left,
const double best_right,
574 const bool use_consensus_,
575 double & total_intensity,
576 const SpectrumT & master_peak_container,
577 const std::vector< double > & left_edges,
578 const std::vector< double > & right_edges,
590 double local_left = best_left;
591 double local_right = best_right;
592 if (!use_consensus_ && right_edges.size() > prec_idx && left_edges.size() > prec_idx)
594 local_left = left_edges[prec_idx];
595 local_right = right_edges[prec_idx];
598 SpectrumT used_chromatogram;
600 if (peak_integration_ ==
"original")
602 used_chromatogram = resampleChromatogram_(chromatogram, master_peak_container, local_left, local_right);
605 else if (peak_integration_ ==
"smoothed" && smoothed_chroms.size() <= prec_idx)
608 "Tried to calculate peak area and height without any smoothed chromatograms for precursors");
610 else if (peak_integration_ ==
"smoothed")
612 used_chromatogram = resampleChromatogram_(smoothed_chroms[prec_idx], master_peak_container, local_left, local_right);
617 String(
"Peak integration chromatogram ") + peak_integration_ +
" is not a valid method for MRMTransitionGroupPicker");
626 double peak_integral = pa.
area;
627 double peak_apex_int = pa.
height;
629 if (background_subtraction_ !=
"none")
631 double background{0};
632 double avg_noise_level{0};
633 if ((peak_integration_ ==
"smoothed") && smoothed_chroms.size() <= prec_idx)
636 "Tried to calculate background estimation without any smoothed chromatograms");
638 else if (background_subtraction_ ==
"original")
640 const double intensity_left = chromatogram.PosBegin(local_left)->getIntensity();
641 const double intensity_right = (chromatogram.PosEnd(local_right) - 1)->getIntensity();
642 const UInt n_points = std::distance(chromatogram.PosBegin(local_left), chromatogram.PosEnd(local_right));
643 avg_noise_level = (intensity_right + intensity_left) / 2;
644 background = avg_noise_level * n_points;
646 else if (background_subtraction_ ==
"exact")
649 background = pb.
area;
650 avg_noise_level = pb.
height;
652 peak_integral -= background;
653 peak_apex_int -= avg_noise_level;
654 if (peak_integral < 0) {peak_integral = 0;}
655 if (peak_apex_int < 0) {peak_apex_int = 0;}
658 f.
setMetaValue(
"noise_background_level", avg_noise_level);
661 f.
setMZ(chromatogram.getPrecursor().getMZ());
662 if (
k == 0) {mrmFeature.
setMZ(chromatogram.getPrecursor().getMZ());}
664 if (chromatogram.metaValueExists(
"precursor_mz"))
666 f.
setMZ(chromatogram.getMetaValue(
"precursor_mz"));
667 if (
k == 0) {mrmFeature.
setMZ(chromatogram.getMetaValue(
"precursor_mz"));}
670 f.
setRT(picked_chroms[chr_idx][peak_idx].getMZ());
675 f.
setMetaValue(
"native_id", chromatogram.getNativeID());
680 total_intensity += peak_integral;
700 template <
typename SpectrumT>
704 Size count_inside = 0;
705 for (
Size k = 0;
k < picked_chroms.size();
k++)
707 for (
Size i = 0; i < picked_chroms[
k].size(); i++)
709 if (picked_chroms[
k][i].getMZ() >= best_left && picked_chroms[
k][i].getMZ() <= best_right)
711 picked_chroms[
k][i].setIntensity(0.0);
717 Size count_overlap = 0;
719 for (
Size k = 0;
k < picked_chroms.size();
k++)
721 for (
Size i = 0; i < picked_chroms[
k].size(); i++)
723 if (picked_chroms[
k][i].getIntensity() <= 1e-11) {
continue; }
727 if ((left > best_left && left < best_right)
728 || (right > best_left && right < best_right))
730 picked_chroms[
k][i].setIntensity(0.0);
735 OPENMS_LOG_DEBUG <<
" ** Removed " << count_inside <<
" peaks enclosed in and " <<
736 count_overlap <<
" peaks overlapping with added feature." << std::endl;
740 void findLargestPeak(
const std::vector<MSChromatogram >& picked_chroms,
int& chr_idx,
int& peak_idx);
763 template <
typename SpectrumT,
typename TransitionT>
776 throw Exception::IllegalArgument(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
"Did not find chromatogram for id '" + native_id +
"'.");
796 template <
typename SpectrumT,
typename TransitionT>
798 const std::vector<SpectrumT>& picked_chroms,
800 const double best_left,
801 const double best_right,
807 double resample_boundary = resample_boundary_;
808 SpectrumT master_peak_container;
809 const SpectrumT& ref_chromatogram = selectChromHelper_(transition_group, picked_chroms[chr_idx].getNativeID());
810 prepareMasterContainer_(ref_chromatogram, master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
811 std::vector<std::vector<double> > all_ints;
812 for (
Size k = 0;
k < picked_chroms.size();
k++)
814 const SpectrumT& chromatogram = selectChromHelper_(transition_group, picked_chroms[
k].getNativeID());
815 const SpectrumT used_chromatogram = resampleChromatogram_(chromatogram,
816 master_peak_container, best_left - resample_boundary, best_right + resample_boundary);
818 std::vector<double> int_here;
819 for (
const auto& peak : used_chromatogram) int_here.push_back(peak.getIntensity());
821 double tic = std::accumulate(int_here.begin(), int_here.end(), 0.0);
822 if (tic > 1e-11) all_ints.push_back(int_here);
826 std::vector<double> mean_shapes;
827 std::vector<double> mean_coel;
828 for (
Size k = 0;
k < all_ints.size();
k++)
830 std::vector<double> shapes;
831 std::vector<double> coel;
832 for (
Size i = 0; i < all_ints.size(); i++)
834 if (i ==
k) {
continue;}
836 all_ints[
k], all_ints[i], boost::numeric_cast<int>(all_ints[i].size()), 1);
842 shapes.push_back(res_shape);
843 coel.push_back(res_coelution);
849 msc = std::for_each(shapes.begin(), shapes.end(), msc);
850 double shapes_mean = msc.
mean();
851 msc = std::for_each(coel.begin(), coel.end(), msc);
852 double coel_mean = msc.
mean();
856 mean_shapes.push_back(shapes_mean);
857 mean_coel.push_back(coel_mean);
863 int min_index_shape = std::distance(mean_shapes.begin(), std::min_element(mean_shapes.begin(), mean_shapes.end()));
864 int max_index_coel = std::distance(mean_coel.begin(), std::max_element(mean_coel.begin(), mean_coel.end()));
867 int missing_peaks = 0;
868 int multiple_peaks = 0;
871 std::vector<double> left_borders;
872 std::vector<double> right_borders;
873 for (
Size k = 0;
k < picked_chroms.size();
k++)
882 for (
Size i = 0; i < picked_chroms[
k].size(); i++)
884 if (picked_chroms[
k][i].getMZ() >= best_left && picked_chroms[
k][i].getMZ() <= best_right)
887 if (picked_chroms[
k][i].getIntensity() > max_int)
889 max_int = picked_chroms[
k][i].getIntensity() > max_int;
896 if (l_tmp > 1e-11) left_borders.push_back(l_tmp);
897 if (r_tmp > 1e-11) right_borders.push_back(r_tmp);
899 if (pfound == 0) missing_peaks++;
900 if (pfound > 1) multiple_peaks++;
905 OPENMS_LOG_DEBUG <<
" Overall found missing : " << missing_peaks <<
" and multiple : " << multiple_peaks << std::endl;
911 if (min_index_shape == max_index_coel)
913 OPENMS_LOG_DEBUG <<
" Element " << min_index_shape <<
" is a candidate for removal ... " << std::endl;
914 outlier =
String(picked_chroms[min_index_shape].getNativeID());
926 double shape_score = std::accumulate(mean_shapes.begin(), mean_shapes.end(), 0.0) / mean_shapes.size();
927 double coel_score = std::accumulate(mean_coel.begin(), mean_coel.end(), 0.0) / mean_coel.size();
928 coel_score = (coel_score - 1.0) / 2.0;
930 double score = shape_score - coel_score - 1.0 * missing_peaks / picked_chroms.size();
932 OPENMS_LOG_DEBUG <<
" Computed score " << score <<
" (from " << shape_score <<
933 " - " << coel_score <<
" - " << 1.0 * missing_peaks / picked_chroms.size() <<
")" << std::endl;
947 template <
typename SpectrumT>
948 void recalculatePeakBorders_(
const std::vector<SpectrumT>& picked_chroms,
double& best_left,
double& best_right,
double max_z)
956 std::vector<double> left_borders;
957 std::vector<double> right_borders;
958 for (
Size k = 0;
k < picked_chroms.size();
k++)
963 for (
Size i = 0; i < picked_chroms[
k].size(); i++)
965 if (picked_chroms[
k][i].getMZ() >= best_left && picked_chroms[
k][i].getMZ() <= best_right)
977 left_borders.push_back(left);
978 right_borders.push_back(right);
979 OPENMS_LOG_DEBUG <<
" * peak " <<
k <<
" left boundary " << left_borders.back() <<
" with inty " << max_int << std::endl;
980 OPENMS_LOG_DEBUG <<
" * peak " <<
k <<
" right boundary " << right_borders.back() <<
" with inty " << max_int << std::endl;
985 if (right_borders.empty())
1003 mean = std::accumulate(right_borders.begin(), right_borders.end(), 0.0) / (double) right_borders.size();
1004 stdev = std::sqrt(std::inner_product(right_borders.begin(), right_borders.end(), right_borders.begin(), 0.0)
1005 / right_borders.size() -
mean *
mean);
1006 std::sort(right_borders.begin(), right_borders.end());
1009 << best_right <<
" std " << stdev <<
" : " << std::fabs(best_right -
mean) / stdev
1010 <<
" coefficient of variation" << std::endl;
1013 if (std::fabs(best_right -
mean) / stdev > max_z)
1015 best_right = right_borders[right_borders.size() / 2];
1016 OPENMS_LOG_DEBUG <<
" - CV too large: correcting right boundary to " << best_right << std::endl;
1020 mean = std::accumulate(left_borders.begin(), left_borders.end(), 0.0) / (double) left_borders.size();
1021 stdev = std::sqrt(std::inner_product(left_borders.begin(), left_borders.end(), left_borders.begin(), 0.0)
1022 / left_borders.size() -
mean *
mean);
1023 std::sort(left_borders.begin(), left_borders.end());
1026 << best_left <<
" std " << stdev <<
" : " << std::fabs(best_left -
mean) / stdev
1027 <<
" coefficient of variation" << std::endl;
1030 if (std::fabs(best_left -
mean) / stdev > max_z)
1032 best_left = left_borders[left_borders.size() / 2];
1033 OPENMS_LOG_DEBUG <<
" - CV too large: correcting left boundary to " << best_left << std::endl;
1055 template <
typename SpectrumT>
1057 SpectrumT& master_peak_container,
double left_boundary,
double right_boundary)
1064 typename SpectrumT::const_iterator begin = ref_chromatogram.begin();
1065 while (begin != ref_chromatogram.end() && begin->getMZ() < left_boundary) {begin++; }
1066 if (begin != ref_chromatogram.begin()) {begin--; }
1068 typename SpectrumT::const_iterator end = begin;
1069 while (end != ref_chromatogram.end() && end->getMZ() < right_boundary) {end++; }
1070 if (end != ref_chromatogram.end()) {end++; }
1073 master_peak_container.resize(distance(begin, end));
1074 typename SpectrumT::iterator it = master_peak_container.begin();
1075 for (
typename SpectrumT::const_iterator chrom_it = begin; chrom_it != end; chrom_it++, it++)
1077 it->setMZ(chrom_it->getMZ());
1091 template <
typename SpectrumT>
1093 const SpectrumT& master_peak_container,
double left_boundary,
double right_boundary)
1098 typename SpectrumT::const_iterator begin = chromatogram.begin();
1099 while (begin != chromatogram.end() && begin->getMZ() < left_boundary) {begin++;}
1100 if (begin != chromatogram.begin()) {begin--;}
1102 typename SpectrumT::const_iterator end = begin;
1103 while (end != chromatogram.end() && end->getMZ() < right_boundary) {end++;}
1104 if (end != chromatogram.end()) {end++;}
1106 SpectrumT resampled_peak_container = master_peak_container;
1108 lresampler.
raster(begin, end, resampled_peak_container.begin(), resampled_peak_container.end());
1110 return resampled_peak_container;
#define OPENMS_LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:454
const String & getNativeID() const
returns the native identifier for the spectrum, used by the acquisition software.
void setNativeID(const String &native_id)
sets the native identifier for the spectrum, used by the acquisition software.
Definition: ConvexHull2D.h:47
void setHullPoints(const PointArrayType &points)
accessor for the outer(!) points (no checking is performed if this is actually a convex hull)
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:66
A method or algorithm argument contains illegal values.
Definition: Exception.h:616
An LC-MS feature.
Definition: Feature.h:46
const std::vector< ConvexHull2D > & getConvexHulls() const
Non-mutable access to the convex hulls.
void setQuality(Size index, QualityType q)
Set the quality in dimension c.
void setOverallQuality(QualityType q)
Set the overall quality.
Linear Resampling of raw data with alignment.
Definition: LinearResamplerAlign.h:33
void raster(SpecT &container)
Applies the resampling algorithm to a container (MSSpectrum or MSChromatogram).
Definition: LinearResamplerAlign.h:54
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:26
void addFeature(const Feature &feature, const String &key)
Adds an feature from a single chromatogram into the feature.
void addPrecursorFeature(const Feature &feature, const String &key)
Adds a precursor feature from a single chromatogram into the feature.
The MRMTransitionGroupPicker finds peaks in chromatograms that belong to the same precursors.
Definition: MRMTransitionGroupPicker.h:55
void findLargestPeak(const std::vector< MSChromatogram > &picked_chroms, int &chr_idx, int &peak_idx)
Find largest peak in a vector of chromatograms.
String boundary_selection_method_
Which method to use for selecting peaks' boundaries.
Definition: MRMTransitionGroupPicker.h:1137
String peak_integration_
Definition: MRMTransitionGroupPicker.h:1116
bool compute_total_mi_
Definition: MRMTransitionGroupPicker.h:1123
const SpectrumT & selectChromHelper_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const String &native_id)
Select matching precursor or fragment ion chromatogram.
Definition: MRMTransitionGroupPicker.h:764
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:1056
String background_subtraction_
Definition: MRMTransitionGroupPicker.h:1117
bool compute_peak_quality_
Definition: MRMTransitionGroupPicker.h:1121
bool use_precursors_
Definition: MRMTransitionGroupPicker.h:1119
void pickApex(std::vector< SpectrumT > &picked_chroms, const double best_left, const double best_right, const double peak_apex, double &min_left, double &max_right, std::vector< double > &left_edges, std::vector< double > &right_edges)
Apex-based peak picking.
Definition: MRMTransitionGroupPicker.h:325
double min_peak_width_
Definition: MRMTransitionGroupPicker.h:1128
void remove_overlapping_features(std::vector< SpectrumT > &picked_chroms, double best_left, double best_right)
Remove overlapping features.
Definition: MRMTransitionGroupPicker.h:701
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:948
MRMTransitionGroupPicker & operator=(const MRMTransitionGroupPicker &rhs)
Assignment operator is protected for algorithm.
int stop_after_feature_
Definition: MRMTransitionGroupPicker.h:1126
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:1092
PeakPickerChromatogram picker_
Definition: MRMTransitionGroupPicker.h:1139
double recalculate_peaks_max_z_
Definition: MRMTransitionGroupPicker.h:1129
PeakIntegrator pi_
Definition: MRMTransitionGroupPicker.h:1140
void findWidestPeakIndices(const std::vector< MSChromatogram > &picked_chroms, Int &chrom_idx, Int &point_idx) const
Given a vector of chromatograms, find the indices of the chromatogram containing the widest peak and ...
double min_qual_
Definition: MRMTransitionGroupPicker.h:1124
double resample_boundary_
Definition: MRMTransitionGroupPicker.h:1130
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:200
~MRMTransitionGroupPicker() override
Destructor.
void updateMembers_() override
Synchronize members with param class.
double computeQuality_(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, const int chr_idx, const double best_left, const double best_right, String &outlier)
Compute transition group quality (higher score is better)
Definition: MRMTransitionGroupPicker.h:797
double stop_after_intensity_ratio_
Definition: MRMTransitionGroupPicker.h:1127
void pickTransitionGroup(MRMTransitionGroup< SpectrumT, TransitionT > &transition_group)
Pick a group of chromatograms belonging to the same peptide.
Definition: MRMTransitionGroupPicker.h:87
bool use_consensus_
Definition: MRMTransitionGroupPicker.h:1120
void pickFragmentChromatograms(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, MRMFeature &mrmFeature, const std::vector< SpectrumT > &smoothed_chroms, const double best_left, const double best_right, const bool use_consensus_, double &total_intensity, double &total_xic, double &total_mi, double &total_peak_apices, const SpectrumT &master_peak_container, const std::vector< double > &left_edges, const std::vector< double > &right_edges, const int chr_idx, const int peak_idx)
Definition: MRMTransitionGroupPicker.h:366
bool recalculate_peaks_
Definition: MRMTransitionGroupPicker.h:1118
void pickPrecursorChromatograms(const MRMTransitionGroup< SpectrumT, TransitionT > &transition_group, const std::vector< SpectrumT > &picked_chroms, MRMFeature &mrmFeature, const std::vector< SpectrumT > &smoothed_chroms, const double best_left, const double best_right, const bool use_consensus_, double &total_intensity, const SpectrumT &master_peak_container, const std::vector< double > &left_edges, const std::vector< double > &right_edges, const int chr_idx, const int peak_idx)
Definition: MRMTransitionGroupPicker.h:569
bool compute_peak_shape_metrics_
Definition: MRMTransitionGroupPicker.h:1122
std::vector< ChromatogramType > & getPrecursorChromatograms()
Definition: MRMTransitionGroup.h:211
ChromatogramType & getChromatogram(const String &key)
Definition: MRMTransitionGroup.h:193
bool hasPrecursorChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:242
const String & getTransitionGroupID() const
Definition: MRMTransitionGroup.h:104
bool isInternallyConsistent() const
Check whether internal state is consistent, e.g. same number of chromatograms and transitions are pre...
Definition: MRMTransitionGroup.h:291
bool chromatogramIdsMatch() const
Ensure that chromatogram native ids match their keys in the map.
Definition: MRMTransitionGroup.h:300
bool hasTransition(const String &key) const
Definition: MRMTransitionGroup.h:144
const TransitionType & getTransition(const String &key)
Definition: MRMTransitionGroup.h:149
bool hasChromatogram(const String &key) const
Definition: MRMTransitionGroup.h:188
void addFeature(const MRMFeature &feature)
Definition: MRMTransitionGroup.h:275
const std::vector< TransitionType > & getTransitions() const
Definition: MRMTransitionGroup.h:116
ChromatogramType & getPrecursorChromatogram(const String &key)
Definition: MRMTransitionGroup.h:247
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:160
The representation of a chromatogram.
Definition: MSChromatogram.h:30
void sortByIntensity(bool reverse=false)
Lexicographically sorts the peaks by their intensity.
void setMZ(CoordinateType coordinate)
Mutable access to the m/z coordinate (index 1)
Definition: Peak2D.h:178
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:190
IntensityType getIntensity() const
Definition: Peak2D.h:142
void setIntensity(IntensityType intensity)
Sets data point intensity (height)
Definition: Peak2D.h:148
Compute the area, background and shape metrics of a peak.
Definition: PeakIntegrator.h:48
Int points_across_half_height
Definition: PeakIntegrator.h:185
double apex_pos
Definition: PeakIntegrator.h:72
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:76
double width_at_5
Definition: PeakIntegrator.h:107
Int points_across_baseline
Definition: PeakIntegrator.h:181
double width_at_50
Definition: PeakIntegrator.h:115
double end_position_at_50
Definition: PeakIntegrator.h:139
double baseline_delta_2_height
Definition: PeakIntegrator.h:177
double height
Definition: PeakIntegrator.h:68
double end_position_at_10
Definition: PeakIntegrator.h:135
double width_at_10
Definition: PeakIntegrator.h:111
double total_width
Definition: PeakIntegrator.h:143
double end_position_at_5
Definition: PeakIntegrator.h:131
double start_position_at_50
Definition: PeakIntegrator.h:127
double tailing_factor
Definition: PeakIntegrator.h:157
double slope_of_baseline
Definition: PeakIntegrator.h:172
double start_position_at_10
Definition: PeakIntegrator.h:123
double asymmetry_factor
Definition: PeakIntegrator.h:167
double start_position_at_5
Definition: PeakIntegrator.h:119
double area
Definition: PeakIntegrator.h:64
Definition: PeakIntegrator.h:60
Definition: PeakIntegrator.h:85
Definition: PeakIntegrator.h:103
The PeakPickerChromatogram finds peaks a single chromatogram.
Definition: PeakPickerChromatogram.h:44
@ IDX_LEFTBORDER
Definition: PeakPickerChromatogram.h:57
@ IDX_RIGHTBORDER
Definition: PeakPickerChromatogram.h:57
@ IDX_ABUNDANCE
Definition: PeakPickerChromatogram.h:57
A more convenient string class.
Definition: String.h:34
Size ensureUniqueId()
Assigns a valid unique id, but only if the present one is invalid. Returns 1 if the unique id was cha...
functor to compute the mean and stddev of sequence using the std::foreach algorithm
Definition: StatsHelpers.h:134
double mean() const
Definition: StatsHelpers.h:171
int Int
Signed integer type.
Definition: Types.h:72
unsigned int UInt
Unsigned integer type.
Definition: Types.h:64
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:97
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:94
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:94
const double k
Definition: Constants.h:132
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19
OPENSWATHALGO_DLLAPI double rankedMutualInformation(std::vector< unsigned int > &ranked_data1, std::vector< unsigned int > &ranked_data2, const unsigned int max_rank1, const unsigned int max_rank2)
OPENSWATHALGO_DLLAPI XCorrArrayType::const_iterator xcorrArrayGetMaxPeak(const XCorrArrayType &array)
Find best peak in an cross-correlation (highest apex)
OPENSWATHALGO_DLLAPI XCorrArrayType normalizedCrossCorrelation(std::vector< double > &data1, std::vector< double > &data2, const int maxdelay, const int lag)
OPENSWATHALGO_DLLAPI unsigned int computeAndAppendRank(const std::vector< double > &v, std::vector< unsigned int > &ranks)