38 #include <boost/math/special_functions/fpclassify.hpp>
82 template<
typename SpectrumType,
class TransitionT>
86 bool smooth_data =
false;
101 avg_score /= transition_group.
size();
110 if (current_section.size() < 2)
117 typedef Peak1D LocalPeakType;
120 std::vector<LocalPeakType> data_to_fit;
121 prepareFit_(current_section, data_to_fit, smooth_data);
122 std::unique_ptr<InterpolationModel> model_rt;
123 double quality =
fitRT_(data_to_fit, model_rt);
130 template<
class LocalPeakType>
131 double fitRT_(std::vector<LocalPeakType>& rt_input_data, std::unique_ptr<InterpolationModel>& model)
const
137 return fitter_emg1D.
fit1d(rt_input_data, model);
142 template<
class LocalPeakType>
148 for (ConvexHull2D::PointArrayType::const_iterator it = current_section.begin(); it != current_section.end(); it++)
152 p.setIntensity(it->getY());
153 filter_spec.push_back(p);
158 std::vector<double> distances;
159 for (
Size j = 1; j < filter_spec.size(); ++j)
161 distances.push_back(filter_spec[j].getMZ() - filter_spec[j - 1].getMZ());
163 double dist_average = std::accumulate(distances.begin(), distances.end(), 0.0) / (double) distances.size();
168 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
169 filter_spec.push_back(new_peak);
170 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
171 filter_spec.push_back(new_peak);
172 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
173 filter_spec.push_back(new_peak);
176 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
177 filter_spec.insert(filter_spec.begin(), new_peak);
178 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
179 filter_spec.insert(filter_spec.begin(), new_peak);
180 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
181 filter_spec.insert(filter_spec.begin(), new_peak);
190 filter_param.
setValue(
"gaussian_width", 4 * dist_average);
192 filter.
filter(filter_spec);
196 for (
Size j = 0; j != filter_spec.size(); ++j)
199 p.setPosition(filter_spec[j].getMZ());
200 p.setIntensity(filter_spec[j].getIntensity());
201 data_to_fit.push_back(p);
std::vector< PointType > PointArrayType
Definition: ConvexHull2D.h:76
const Param & getParameters() const
Non-mutable access to the parameters.
void setParameters(const Param ¶m)
Sets the parameters.
const Param & getDefaults() const
Non-mutable access to the default parameters.
Exponentially modified gaussian distribution fitter (1-dim.) using Levenberg-Marquardt algorithm (Eig...
Definition: EmgFitter1D.h:49
QualityType fit1d(const RawDataArrayType &range, std::unique_ptr< InterpolationModel > &model) override
return interpolation model
Scoring of an elution peak using an exponentially modified gaussian distribution model.
Definition: EmgScoring.h:60
void setFitterParam(const Param ¶m)
Definition: EmgScoring.h:70
double elutionModelFit(const ConvexHull2D::PointArrayType ¤t_section, bool smooth_data) const
Definition: EmgScoring.h:107
Param fitter_emg1D_params_
Definition: EmgScoring.h:205
Param getDefaults()
Get default params for the Emg1D fitting.
Definition: EmgScoring.h:76
double calcElutionFitScore(MRMFeature &mrmfeature, MRMTransitionGroup< SpectrumType, TransitionT > &transition_group) const
calculate the elution profile fit score
Definition: EmgScoring.h:83
double fitRT_(std::vector< LocalPeakType > &rt_input_data, std::unique_ptr< InterpolationModel > &model) const
Definition: EmgScoring.h:131
void prepareFit_(const ConvexHull2D::PointArrayType ¤t_section, std::vector< LocalPeakType > &data_to_fit, bool smooth_data) const
Definition: EmgScoring.h:143
An LC-MS feature.
Definition: Feature.h:72
const std::vector< ConvexHull2D > & getConvexHulls() const
Non-mutable access to the convex hulls.
This class represents a Gaussian lowpass-filter which works on uniform as well as on non-uniform prof...
Definition: GaussFilter.h:77
void filter(MSSpectrum &spectrum)
Smoothes an MSSpectrum containing profile data.
Definition: GaussFilter.h:92
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:52
Feature & getFeature(const String &key)
get a specified feature
The representation of a group of transitions in a targeted proteomics experiment.
Definition: MRMTransitionGroup.h:68
Size size() const
Definition: MRMTransitionGroup.h:125
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:186
The representation of a 1D spectrum.
Definition: MSSpectrum.h:71
Management and storage of parameters / INI files.
Definition: Param.h:70
void setValue(const std::string &key, const ParamValue &value, const std::string &description="", const std::vector< std::string > &tags=std::vector< std::string >())
Sets a value.
A 1-dimensional raw data point or peak.
Definition: Peak1D.h:54
void setIntensity(IntensityType intensity)
Mutable access to the data point intensity (height)
Definition: Peak1D.h:104
void setMZ(CoordinateType mz)
Mutable access to m/z.
Definition: Peak1D.h:113
A more convenient string class.
Definition: String.h:61
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:120
const double k
Definition: Constants.h:153
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47