11 #include <OpenMS/config.h>
72 double apex_pos = 0.0;
107 double width_at_5 = 0.0;
111 double width_at_10 = 0.0;
115 double width_at_50 = 0.0;
119 double start_position_at_5 = 0.0;
123 double start_position_at_10 = 0.0;
127 double start_position_at_50 = 0.0;
131 double end_position_at_5 = 0.0;
135 double end_position_at_10 = 0.0;
139 double end_position_at_50 = 0.0;
143 double total_width = 0.0;
157 double tailing_factor = 0.0;
167 double asymmetry_factor = 0.0;
172 double slope_of_baseline = 0.0;
177 double baseline_delta_2_height = 0.0;
181 Int points_across_baseline = 0;
185 Int points_across_half_height = 0;
196 static constexpr
const char* INTEGRATION_TYPE_INTENSITYSUM =
"intensity_sum";
198 static constexpr
const char* INTEGRATION_TYPE_TRAPEZOID =
"trapezoid";
200 static constexpr
const char* INTEGRATION_TYPE_SIMPSON =
"simpson";
202 static constexpr
const char* BASELINE_TYPE_BASETOBASE =
"base_to_base";
204 static constexpr
const char* BASELINE_TYPE_VERTICALDIVISION =
"vertical_division";
206 static constexpr
const char* BASELINE_TYPE_VERTICALDIVISION_MIN =
"vertical_division_min";
208 static constexpr
const char* BASELINE_TYPE_VERTICALDIVISION_MAX =
"vertical_division_max";
230 const MSChromatogram& chromatogram,
const double left,
const double right
274 const MSSpectrum& spectrum,
const double left,
const double right
324 const MSChromatogram& chromatogram,
const double left,
const double right,
325 const double peak_apex_pos
354 const double peak_apex_pos
382 const MSSpectrum& spectrum,
const double left,
const double right,
383 const double peak_apex_pos
412 const double peak_apex_pos
435 const MSChromatogram& chromatogram,
const double left,
const double right,
436 const double peak_height,
const double peak_apex_pos
460 const double peak_height,
const double peak_apex_pos
483 const MSSpectrum& spectrum,
const double left,
const double right,
484 const double peak_height,
const double peak_apex_pos
508 const double peak_height,
const double peak_apex_pos
516 template <
typename PeakContainerT>
519 OPENMS_PRECONDITION(left <= right,
"Left peak boundary must be smaller than right boundary!")
520 PeakContainerT emg_pc;
521 const PeakContainerT& p = EMGPreProcess_(pc, emg_pc, left, right);
523 std::function<double(
const double,
const double)>
524 compute_peak_area_trapezoid = [&p](
const double left,
const double right)
526 double peak_area { 0.0 };
527 for (
typename PeakContainerT::ConstIterator it = p.PosBegin(left); it != p.PosEnd(right) - 1; ++it)
529 peak_area += ((it + 1)->getPos() - it->getPos()) * ((it->getIntensity() + (it + 1)->getIntensity()) / 2.0);
534 std::function<double(
const double,
const double)>
535 compute_peak_area_intensity_sum = [&p](
const double left,
const double right)
538 double peak_area { 0.0 };
539 for (
typename PeakContainerT::ConstIterator it = p.PosBegin(left); it != p.PosEnd(right); ++it)
541 peak_area += it->getIntensity();
548 UInt n_points = std::distance(p.PosBegin(left), p.PosEnd(right));
549 for (
auto it = p.PosBegin(left); it != p.PosEnd(right); ++it)
552 if (pa.
height < it->getIntensity())
554 pa.
height = it->getIntensity();
559 if (integration_type_ == INTEGRATION_TYPE_TRAPEZOID)
563 pa.
area = compute_peak_area_trapezoid(left, right);
566 else if (integration_type_ == INTEGRATION_TYPE_SIMPSON)
571 "number of points is 2, falling back to `trapezoid`." << std::endl;
572 pa.
area = compute_peak_area_trapezoid(left, right);
574 else if (n_points > 2)
578 pa.
area = simpson_(p.PosBegin(left), p.PosEnd(right));
582 double areas[4] = {-1.0, -1.0, -1.0, -1.0};
583 areas[0] = simpson_(p.PosBegin(left), p.PosEnd(right) - 1);
584 areas[1] = simpson_(p.PosBegin(left) + 1, p.PosEnd(right));
585 if (p.begin() <= p.PosBegin(left) - 1)
587 areas[2] = simpson_(p.PosBegin(left) - 1, p.PosEnd(right));
589 if (p.PosEnd(right) < p.end())
591 areas[3] = simpson_(p.PosBegin(left), p.PosEnd(right) + 1);
594 for (
const auto& area : areas)
606 else if (integration_type_ == INTEGRATION_TYPE_INTENSITYSUM)
608 pa.
area = compute_peak_area_intensity_sum(left, right);
612 throw Exception::InvalidParameter(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
"Please set a valid value for the parameter \"integration_type\".");
620 template <
typename PeakContainerT>
622 const PeakContainerT& pc,
double left,
double right,
623 const double peak_apex_pos
626 PeakContainerT emg_pc;
627 const PeakContainerT& p = EMGPreProcess_(pc, emg_pc, left, right);
629 const double int_l = p.PosBegin(left)->getIntensity();
630 const double int_r = (p.PosEnd(right) - 1)->getIntensity();
631 const double delta_int = int_r - int_l;
632 const double delta_pos = (p.PosEnd(right) - 1)->getPos() - p.PosBegin(left)->getPos();
633 const double min_int_pos = int_r <= int_l ? (p.PosEnd(right) - 1)->getPos() : p.PosBegin(left)->getPos();
634 const double delta_int_apex = std::fabs(delta_int) * std::fabs(min_int_pos - peak_apex_pos) / delta_pos;
637 if (baseline_type_ == BASELINE_TYPE_BASETOBASE)
639 height = std::min(int_r, int_l) + delta_int_apex;
640 if (integration_type_ == INTEGRATION_TYPE_TRAPEZOID || integration_type_ == INTEGRATION_TYPE_SIMPSON)
644 area = delta_pos * (std::min(int_r, int_l) + 0.5 * std::fabs(delta_int));
646 else if (integration_type_ == INTEGRATION_TYPE_INTENSITYSUM)
653 double pos_sum = 0.0;
654 for (
auto it = p.PosBegin(left); it != p.PosEnd(right); ++it)
656 pos_sum += it->getPos();
658 UInt n_points = std::distance(p.PosBegin(left), p.PosEnd(right));
663 const double rectangle_area = n_points * int_l;
664 const double slope = delta_int / delta_pos;
665 const double triangle_area = (pos_sum - n_points * p.PosBegin(left)->getPos()) * slope;
666 area = triangle_area + rectangle_area;
669 else if (baseline_type_ == BASELINE_TYPE_VERTICALDIVISION || baseline_type_ == BASELINE_TYPE_VERTICALDIVISION_MIN)
671 height = std::min(int_r, int_l);
672 if (integration_type_ == INTEGRATION_TYPE_TRAPEZOID || integration_type_ == INTEGRATION_TYPE_SIMPSON)
674 area = delta_pos * std::min(int_r, int_l);
676 else if (integration_type_ == INTEGRATION_TYPE_INTENSITYSUM)
678 area = std::min(int_r, int_l) * std::distance(p.PosBegin(left), p.PosEnd(right));;
681 else if (baseline_type_ == BASELINE_TYPE_VERTICALDIVISION_MAX)
683 height = std::max(int_r, int_l);
684 if (integration_type_ == INTEGRATION_TYPE_TRAPEZOID || integration_type_ == INTEGRATION_TYPE_SIMPSON)
686 area = delta_pos * std::max(int_r, int_l);
688 else if (integration_type_ == INTEGRATION_TYPE_INTENSITYSUM)
690 area = std::max(int_r, int_l) * std::distance(p.PosBegin(left), p.PosEnd(right));
695 throw Exception::InvalidParameter(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION,
"Please set a valid value for the parameter \"baseline_type\".");
717 template <
typename PeakContainerConstIteratorT>
718 double simpson_(PeakContainerConstIteratorT it_begin, PeakContainerConstIteratorT it_end)
const
720 double integral = 0.0;
721 for (
auto it = it_begin + 1; it < it_end - 1; it = it + 2)
723 const double h = it->getPos() - (it - 1)->getPos();
724 const double k = (it + 1)->getPos() - it->getPos();
725 const double y_h = (it - 1)->getIntensity();
726 const double y_0 = it->getIntensity();
727 const double y_k = (it + 1)->getIntensity();
728 integral += (1.0 / 6.0) * (
h +
k) * ((2.0 -
k /
h) * y_h + (pow(
h +
k, 2) / (
h *
k)) * y_0 + (2.0 -
h /
k) * y_k);
734 template <
typename PeakContainerT>
736 const PeakContainerT& pc,
double left,
double right,
737 const double peak_height,
const double peak_apex_pos
742 if (pc.empty())
return psm;
745 if (!(left <= peak_apex_pos && peak_apex_pos <= right))
throw Exception::InvalidRange(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION);
747 PeakContainerT emg_pc;
748 const PeakContainerT& p = EMGPreProcess_(pc, emg_pc, left, right);
750 typename PeakContainerT::ConstIterator it_PosBegin_l = p.PosBegin(left);
751 typename PeakContainerT::ConstIterator it_PosEnd_apex = p.PosBegin(peak_apex_pos);
752 typename PeakContainerT::ConstIterator it_PosEnd_r = p.PosEnd(right);
753 for (
auto it = it_PosBegin_l; it != it_PosEnd_r; ++it)
757 if (it->getIntensity() >= 0.5 * peak_height)
763 psm.
start_position_at_5 = findPosAtPeakHeightPercent_(it_PosBegin_l, it_PosEnd_apex, p.end(), peak_height, 0.05,
true);
764 psm.
start_position_at_10 = findPosAtPeakHeightPercent_(it_PosBegin_l, it_PosEnd_apex, p.end(), peak_height, 0.1,
true);
765 psm.
start_position_at_50 = findPosAtPeakHeightPercent_(it_PosBegin_l, it_PosEnd_apex, p.end(), peak_height, 0.5,
true);
766 psm.
end_position_at_5 = findPosAtPeakHeightPercent_(it_PosEnd_apex, it_PosEnd_r, p.end(), peak_height, 0.05,
false);
767 psm.
end_position_at_10 = findPosAtPeakHeightPercent_(it_PosEnd_apex, it_PosEnd_r, p.end(), peak_height, 0.1,
false);
768 psm.
end_position_at_50 = findPosAtPeakHeightPercent_(it_PosEnd_apex, it_PosEnd_r, p.end(), peak_height, 0.5,
false);
773 psm.
total_width = (p.PosEnd(right) - 1)->getPos() - p.PosBegin(left)->getPos();
774 psm.
slope_of_baseline = (p.PosEnd(right) - 1)->getIntensity() - p.PosBegin(left)->getIntensity();
775 if (peak_height != 0.0)
816 template <
typename PeakContainerConstIteratorT>
818 PeakContainerConstIteratorT it_left,
819 PeakContainerConstIteratorT it_right,
820 PeakContainerConstIteratorT it_end,
821 const double peak_height,
822 const double percent,
823 const bool is_left_half)
const
829 if (it_left == it_right)
return it_left->getPos();
831 const double percent_intensity = peak_height * percent;
832 PeakContainerConstIteratorT closest;
837 PeakContainerConstIteratorT it = it_left;
838 it < it_right && it->getIntensity() <= percent_intensity;
844 closest = it_right - 1;
846 PeakContainerConstIteratorT it = it_right - 1;
847 it >= it_left && it->getIntensity() <= percent_intensity;
851 return closest->getPos();
864 String integration_type_ = INTEGRATION_TYPE_INTENSITYSUM;
869 String baseline_type_ = BASELINE_TYPE_BASETOBASE;
890 template <
typename PeakContainerT>
892 const PeakContainerT& pc,
893 PeakContainerT& emg_pc,
901 left = emg_pc.front().getPos();
902 right = emg_pc.back().getPos();
#define OPENMS_LOG_WARN
Macro if a warning, a piece of information which should be read by the user, should be logged.
Definition: LogStream.h:444
std::vector< PointType > PointArrayType
Definition: ConvexHull2D.h:50
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:66
Compute the area, background and shape metrics of a peak.
Definition: EmgGradientDescent.h:40
void fitEMGPeakModel(const PeakContainerT &input_peak, PeakContainerT &output_peak, const double left_pos=0.0, const double right_pos=0.0) const
Fit the given peak (either MSChromatogram or MSSpectrum) to the EMG peak model.
Exception indicating that an invalid parameter was handed over to an algorithm.
Definition: Exception.h:316
Invalid range exception.
Definition: Exception.h:257
The representation of a chromatogram.
Definition: MSChromatogram.h:30
ContainerType::const_iterator ConstIterator
Non-mutable iterator.
Definition: MSChromatogram.h:66
The representation of a 1D spectrum.
Definition: MSSpectrum.h:44
ContainerType::const_iterator ConstIterator
Non-mutable iterator.
Definition: MSSpectrum.h:110
Management and storage of parameters / INI files.
Definition: Param.h:44
Compute the area, background and shape metrics of a peak.
Definition: PeakIntegrator.h:48
Int points_across_half_height
Definition: PeakIntegrator.h:185
bool fit_EMG_
Enable/disable EMG peak model fitting.
Definition: PeakIntegrator.h:873
PeakArea integratePeak(const MSChromatogram &chromatogram, const double left, const double right) const
Compute the area of a peak contained in a MSChromatogram.
double apex_pos
Definition: PeakIntegrator.h:72
ConvexHull2D::PointArrayType hull_points
Definition: PeakIntegrator.h:76
PeakBackground estimateBackground_(const PeakContainerT &pc, double left, double right, const double peak_apex_pos) const
Definition: PeakIntegrator.h:621
PeakArea integratePeak(const MSChromatogram &chromatogram, MSChromatogram::ConstIterator &left, MSChromatogram::ConstIterator &right) const
Compute the area of a peak contained in a MSChromatogram.
double width_at_5
Definition: PeakIntegrator.h:107
Int points_across_baseline
Definition: PeakIntegrator.h:181
PeakShapeMetrics calculatePeakShapeMetrics(const MSSpectrum &spectrum, MSSpectrum::ConstIterator &left, MSSpectrum::ConstIterator &right, const double peak_height, const double peak_apex_pos) const
Calculate peak's shape metrics.
double width_at_50
Definition: PeakIntegrator.h:115
double end_position_at_50
Definition: PeakIntegrator.h:139
double findPosAtPeakHeightPercent_(PeakContainerConstIteratorT it_left, PeakContainerConstIteratorT it_right, PeakContainerConstIteratorT it_end, const double peak_height, const double percent, const bool is_left_half) const
Find the position (RT/MZ) at a given percentage of peak's height.
Definition: PeakIntegrator.h:817
double baseline_delta_2_height
Definition: PeakIntegrator.h:177
EmgGradientDescent emg_
Definition: PeakIntegrator.h:874
double height
Definition: PeakIntegrator.h:68
double end_position_at_10
Definition: PeakIntegrator.h:135
double simpson_(PeakContainerConstIteratorT it_begin, PeakContainerConstIteratorT it_end) const
Simpson's rule algorithm.
Definition: PeakIntegrator.h:718
~PeakIntegrator() override
Destructor.
PeakShapeMetrics calculatePeakShapeMetrics(const MSChromatogram &chromatogram, MSChromatogram::ConstIterator &left, MSChromatogram::ConstIterator &right, const double peak_height, const double peak_apex_pos) const
Calculate peak's shape metrics.
PeakShapeMetrics calculatePeakShapeMetrics_(const PeakContainerT &pc, double left, double right, const double peak_height, const double peak_apex_pos) const
Definition: PeakIntegrator.h:735
PeakBackground estimateBackground(const MSSpectrum &spectrum, MSSpectrum::ConstIterator &left, MSSpectrum::ConstIterator &right, const double peak_apex_pos) const
Estimate the background of a peak contained in a MSSpectrum.
double width_at_10
Definition: PeakIntegrator.h:111
PeakArea integratePeak(const MSSpectrum &spectrum, MSSpectrum::ConstIterator &left, MSSpectrum::ConstIterator &right) const
Compute the area of a peak contained in a MSSpectrum.
double total_width
Definition: PeakIntegrator.h:143
PeakShapeMetrics calculatePeakShapeMetrics(const MSChromatogram &chromatogram, const double left, const double right, const double peak_height, const double peak_apex_pos) const
Calculate peak's shape metrics.
double end_position_at_5
Definition: PeakIntegrator.h:131
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
double start_position_at_50
Definition: PeakIntegrator.h:127
void getDefaultParameters(Param ¶ms)
double tailing_factor
Definition: PeakIntegrator.h:157
PeakShapeMetrics calculatePeakShapeMetrics(const MSSpectrum &spectrum, const double left, const double right, const double peak_height, const double peak_apex_pos) const
Calculate peak's shape metrics.
double slope_of_baseline
Definition: PeakIntegrator.h:172
PeakBackground estimateBackground(const MSChromatogram &chromatogram, const double left, const double right, const double peak_apex_pos) const
Estimate the background of a peak contained in a MSChromatogram.
const PeakContainerT & EMGPreProcess_(const PeakContainerT &pc, PeakContainerT &emg_pc, double &left, double &right) const
Fit the peak to the EMG model.
Definition: PeakIntegrator.h:891
double start_position_at_10
Definition: PeakIntegrator.h:123
PeakArea integratePeak(const MSSpectrum &spectrum, const double left, const double right) const
Compute the area of a peak contained in a MSSpectrum.
double asymmetry_factor
Definition: PeakIntegrator.h:167
double start_position_at_5
Definition: PeakIntegrator.h:119
double area
Definition: PeakIntegrator.h:64
PeakBackground estimateBackground(const MSSpectrum &spectrum, const double left, const double right, const double peak_apex_pos) const
Estimate the background of a peak contained in a MSSpectrum.
PeakBackground estimateBackground(const MSChromatogram &chromatogram, MSChromatogram::ConstIterator &left, MSChromatogram::ConstIterator &right, const double peak_apex_pos) const
Estimate the background of a peak contained in a MSChromatogram.
PeakArea integratePeak_(const PeakContainerT &pc, double left, double right) const
Definition: PeakIntegrator.h:517
PeakIntegrator()
Constructor.
Definition: PeakIntegrator.h:60
Definition: PeakIntegrator.h:85
Definition: PeakIntegrator.h:103
A more convenient string class.
Definition: String.h:34
int Int
Signed integer type.
Definition: Types.h:72
unsigned int UInt
Unsigned integer type.
Definition: Types.h:64
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:94
const double k
Definition: Constants.h:132
const double h
Definition: Constants.h:141
Main OpenMS namespace.
Definition: openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19