|  | OpenMS
    2.6.0
    | 
 
 
  
  
 
Go to the documentation of this file.
   38 #include <boost/math/special_functions/fpclassify.hpp>  
   79     template<
typename SpectrumType, 
class TransitionT>
 
   82       std::vector<double> fit_scores;
 
   84       bool smooth_data = 
false;
 
   94         fit_scores.push_back(fscore);
 
   98       avg_score /= transition_group.
size();
 
  107       if (current_section.size() < 2)
 
  114       typedef Peak1D LocalPeakType;
 
  117       std::vector<LocalPeakType> data_to_fit;
 
  118       prepareFit_(current_section, data_to_fit, smooth_data);
 
  120       double quality = 
fitRT_(data_to_fit, model_rt);
 
  129     template<
class LocalPeakType>
 
  151       if (boost::math::isnan(quality)) quality = -1.0;
 
  157     template<
class LocalPeakType>
 
  163       for (ConvexHull2D::PointArrayType::const_iterator it = current_section.begin(); it != current_section.end(); it++)
 
  167         p.setIntensity(it->getY());
 
  168         filter_spec.push_back(p);
 
  173       std::vector<double> distances;
 
  174       for (
Size j = 1; j < filter_spec.size(); ++j)
 
  176         distances.push_back(filter_spec[j].getMZ() - filter_spec[j - 1].getMZ());
 
  178       double dist_average = std::accumulate(distances.begin(), distances.end(), 0.0) / (
double) distances.size();
 
  183       new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
 
  184       filter_spec.push_back(new_peak);
 
  185       new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
 
  186       filter_spec.push_back(new_peak);
 
  187       new_peak.
setMZ(filter_spec.back().getMZ() + dist_average);
 
  188       filter_spec.push_back(new_peak);
 
  191       new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
 
  192       filter_spec.insert(filter_spec.begin(), new_peak);
 
  193       new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
 
  194       filter_spec.insert(filter_spec.begin(), new_peak);
 
  195       new_peak.
setMZ(filter_spec.front().getMZ() - dist_average);
 
  196       filter_spec.insert(filter_spec.begin(), new_peak);
 
  205         filter_param.setValue(
"gaussian_width", 4 * dist_average);
 
  207         filter.
filter(filter_spec);
 
  211       for (
Size j = 0; j != filter_spec.size(); ++j)
 
  214         p.setPosition(filter_spec[j].getMZ());
 
  215         p.setIntensity(filter_spec[j].getIntensity());
 
  216         data_to_fit.push_back(p);
 
  
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:121
A more convenient string class.
Definition: String.h:59
void filter(MSSpectrum &spectrum)
Smoothes an MSSpectrum containing profile data.
Definition: GaussFilter.h:92
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
~EmgScoring()
Definition: EmgScoring.h:66
void setIntensity(IntensityType intensity)
Mutable access to the data point intensity (height)
Definition: Peak1D.h:112
#define OPENMS_PRECONDITION(condition, message)
Precondition macro.
Definition: openms/include/OpenMS/CONCEPT/Macros.h:136
A multi-chromatogram MRM feature.
Definition: MRMFeature.h:50
std::vector< PointType > PointArrayType
Definition: ConvexHull2D.h:76
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
Feature & getFeature(const String &key)
get a specified feature
void setFitterParam(Param param)
Definition: EmgScoring.h:68
Abstract class for 1D-models that are approximated using linear interpolation.
Definition: InterpolationModel.h:54
Exponentially modified gaussian distribution fitter (1-dim.) using Levenberg-Marquardt algorithm (Eig...
Definition: EmgFitter1D.h:46
void prepareFit_(const ConvexHull2D::PointArrayType ¤t_section, std::vector< LocalPeakType > &data_to_fit, bool smooth_data)
Definition: EmgScoring.h:158
void setParameters(const Param ¶m)
Sets the parameters.
const Param & getDefaults() const
Non-mutable access to the default parameters.
const Param & getParameters() const
Non-mutable access to the parameters.
A 1-dimensional raw data point or peak.
Definition: Peak1D.h:54
double fitRT_(std::vector< LocalPeakType > &rt_input_data, InterpolationModel *&model)
Definition: EmgScoring.h:130
double elutionModelFit(ConvexHull2D::PointArrayType current_section, bool smooth_data)
Definition: EmgScoring.h:104
Scoring of an elution peak using an exponentially modified gaussian distribution model.
Definition: EmgScoring.h:59
An LC-MS feature.
Definition: Feature.h:70
EmgScoring()
Definition: EmgScoring.h:64
The representation of a group of transitions in a targeted proteomics experiment.
Definition: MRMTransitionGroup.h:67
QualityType fit1d(const RawDataArrayType &range, InterpolationModel *&model) override
return interpolation model
Management and storage of parameters / INI files.
Definition: Param.h:73
Param getDefaults()
Definition: EmgScoring.h:73
std::vector< ChromatogramType > & getChromatograms()
Definition: MRMTransitionGroup.h:174
The representation of a 1D spectrum.
Definition: MSSpectrum.h:67
double calcElutionFitScore(MRMFeature &mrmfeature, MRMTransitionGroup< SpectrumType, TransitionT > &transition_group)
calculate the elution profile fit score
Definition: EmgScoring.h:80
Size size() const
Definition: MRMTransitionGroup.h:125
This class represents a Gaussian lowpass-filter which works on uniform as well as on non-uniform prof...
Definition: GaussFilter.h:74
EmgFitter1D fitter_emg1D_
Definition: EmgScoring.h:220