OpenMS
TransformationModelInterpolated.h
Go to the documentation of this file.
1 // Copyright (c) 2002-2023, The OpenMS Team -- EKU Tuebingen, ETH Zurich, and FU Berlin
2 // SPDX-License-Identifier: BSD-3-Clause
3 //
4 // --------------------------------------------------------------------------
5 // $Maintainer: Timo Sachsenberg $
6 // $Authors: Stephan Aiche, Hendrik Weisser $
7 // --------------------------------------------------------------------------
8 
9 #pragma once
10 
11 #include <OpenMS/config.h>
12 
15 
16 namespace OpenMS
17 {
43  class OPENMS_DLLAPI TransformationModelInterpolated :
44  public TransformationModel
45  {
46 public:
47 
56  TransformationModelInterpolated(const DataPoints& data, const Param& params);
57  TransformationModelInterpolated(const std::vector<std::pair<double,double>>& data, const Param& params, bool preprocess);
58 
61 
69  double evaluate(double value) const override;
70 
72  static void getDefaultParameters(Param& params);
73 
82  {
83 public:
84 
91  virtual void init(std::vector<double>& x, std::vector<double>& y) = 0;
92 
100  virtual double eval(const double& x) const = 0;
101 
105  virtual ~Interpolator() {}
106  };
107 
108 private:
110  std::vector<double> x_;
111 
113  std::vector<double> y_;
114 
117 
120 
123 
126 
128  void preprocessDataPoints_(const std::vector<std::pair<double,double>>& data);
129  };
130 
131 } // namespace
132 
Management and storage of parameters / INI files.
Definition: Param.h:44
The class defines a generic interpolation technique used in the TransformationModelInterpolated.
Definition: TransformationModelInterpolated.h:82
virtual ~Interpolator()
d'tor.
Definition: TransformationModelInterpolated.h:105
virtual double eval(const double &x) const =0
Evaluate the underlying interpolation at a specific position x.
virtual void init(std::vector< double > &x, std::vector< double > &y)=0
Initialize the Interpolator.
Interpolation model for transformations.
Definition: TransformationModelInterpolated.h:45
~TransformationModelInterpolated() override
Destructor.
static void getDefaultParameters(Param &params)
Gets the default parameters.
void preprocessDataPoints_(const DataPoints &data)
Preprocesses the incoming data and fills the (private) vectors x_ and y_.
double evaluate(double value) const override
Evaluate the interpolation model at the given value.
TransformationModelLinear * lm_back_
Linear model for extrapolation (back)
Definition: TransformationModelInterpolated.h:122
std::vector< double > y_
Data coordinates y.
Definition: TransformationModelInterpolated.h:113
std::vector< double > x_
Data coordinates x.
Definition: TransformationModelInterpolated.h:110
Interpolator * interp_
Interpolation function.
Definition: TransformationModelInterpolated.h:116
TransformationModelLinear * lm_front_
Linear model for extrapolation (front)
Definition: TransformationModelInterpolated.h:119
void preprocessDataPoints_(const std::vector< std::pair< double, double >> &data)
Preprocesses the incoming data and fills the (private) vectors x_ and y_.
TransformationModelInterpolated(const std::vector< std::pair< double, double >> &data, const Param &params, bool preprocess)
TransformationModelInterpolated(const DataPoints &data, const Param &params)
Constructor.
Linear model for transformations.
Definition: TransformationModelLinear.h:34
Base class for transformation models.
Definition: TransformationModel.h:29
std::vector< DataPoint > DataPoints
Vector of coordinate pairs.
Definition: TransformationModel.h:65
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:22