35 #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_EMGSCORING_H 36 #define OPENMS_TRANSFORMATIONS_FEATUREFINDER_EMGSCORING_H 39 #include <boost/math/special_functions/fpclassify.hpp> 80 template<
typename SpectrumType,
class TransitionT>
83 std::vector<double> fit_scores;
85 bool smooth_data =
false;
95 fit_scores.push_back(fscore);
99 avg_score /= transition_group.
size();
108 if (current_section.size() < 2)
115 typedef Peak1D LocalPeakType;
118 std::vector<LocalPeakType> data_to_fit;
119 prepareFit_(current_section, data_to_fit, smooth_data);
121 double quality =
fitRT_(data_to_fit, model_rt);
130 template<
class LocalPeakType>
152 if (boost::math::isnan(quality)) quality = -1.0;
158 template<
class LocalPeakType>
164 for (ConvexHull2D::PointArrayType::const_iterator it = current_section.begin(); it != current_section.end(); it++)
168 p.setIntensity(it->getY());
169 filter_spec.push_back(p);
174 std::vector<double> distances;
175 for (
Size j = 1; j < filter_spec.size(); ++j)
177 distances.push_back(filter_spec[j].getMZ() - filter_spec[j - 1].getMZ());
179 double dist_average = std::accumulate(distances.begin(), distances.end(), 0.0) / (
double) distances.size();
184 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
185 filter_spec.push_back(new_peak);
186 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
187 filter_spec.push_back(new_peak);
188 new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
189 filter_spec.push_back(new_peak);
192 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
193 filter_spec.insert(filter_spec.begin(), new_peak);
194 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
195 filter_spec.insert(filter_spec.begin(), new_peak);
196 new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
197 filter_spec.insert(filter_spec.begin(), new_peak);
206 filter_param.setValue(
"gaussian_width", 4 * dist_average);
208 filter.
filter(filter_spec);
212 for (
Size j = 0; j != filter_spec.size(); ++j)
215 p.setPosition(filter_spec[j].getMZ());
216 p.setIntensity(filter_spec[j].getIntensity());
217 data_to_fit.push_back(p);
A more convenient string class.
Definition: String.h:57
Scoring of an elution peak using an exponentially modified gaussian distribution model.
Definition: EmgScoring.h:60
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:107
Size size() const
Definition: MRMTransitionGroup.h:126
std::vector< PointType > PointArrayType
Definition: ConvexHull2D.h:77
EmgFitter1D fitter_emg1D_
Definition: EmgScoring.h:221
Exponentially modified gaussian distribution fitter (1-dim.) using Levenberg-Marquardt algorithm (Eig...
Definition: EmgFitter1D.h:47
QualityType fit1d(const RawDataArrayType &range, InterpolationModel *&model)
return interpolation model
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
const std::vector< ConvexHull2D > & getConvexHulls() const
Non-mutable access to the convex hulls.
void setMZ(CoordinateType mz)
Mutable access to m/z.
Definition: Peak1D.h:120
void setIntensity(IntensityType intensity)
Mutable access to the data point intensity (height)
Definition: Peak1D.h:111
void setParameters(const Param ¶m)
Sets the parameters.
The representation of a group of transitions in a targeted proteomics experiment. ...
Definition: MRMTransitionGroup.h:68
Abstract class for 1D-models that are approximated using linear interpolation.
Definition: InterpolationModel.h:55
const Param & getDefaults() const
Non-mutable access to the default parameters.
The representation of a 1D spectrum.
Definition: MSSpectrum.h:67
const std::vector< ChromatogramType > & getChromatograms() const
Definition: MRMTransitionGroup.h:175
A 1-dimensional raw data point or peak.
Definition: Peak1D.h:55
void setFitterParam(Param param)
Definition: EmgScoring.h:69
An LC-MS feature.
Definition: Feature.h:70
const Param & getParameters() const
Non-mutable access to the parameters.
~EmgScoring()
Definition: EmgScoring.h:67
double elutionModelFit(ConvexHull2D::PointArrayType current_section, bool smooth_data)
Definition: EmgScoring.h:105
Management and storage of parameters / INI files.
Definition: Param.h:75
This class represents a Gaussian lowpass-filter which works on uniform as well as on non-uniform prof...
Definition: GaussFilter.h:75
Feature & getFeature(String key)
get a specified feature
double calcElutionFitScore(MRMFeature &mrmfeature, MRMTransitionGroup< SpectrumType, TransitionT > &transition_group)
calculate the elution profile fit score
Definition: EmgScoring.h:81
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:128
void prepareFit_(const ConvexHull2D::PointArrayType ¤t_section, std::vector< LocalPeakType > &data_to_fit, bool smooth_data)
Definition: EmgScoring.h:159
EmgScoring()
Definition: EmgScoring.h:65
void filter(MSSpectrum &spectrum)
Smoothes an MSSpectrum containing profile data.
Definition: GaussFilter.h:93
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:50
Param getDefaults()
Definition: EmgScoring.h:74
double fitRT_(std::vector< LocalPeakType > &rt_input_data, InterpolationModel *&model)
Definition: EmgScoring.h:131