Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
InterpolationModel.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2017.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Timo Sachsenberg $
32 // $Authors: $
33 // --------------------------------------------------------------------------
34 
35 
36 #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_INTERPOLATIONMODEL_H
37 #define OPENMS_TRANSFORMATIONS_FEATUREFINDER_INTERPOLATIONMODEL_H
38 
41 
42 namespace OpenMS
43 {
55  class OPENMS_DLLAPI InterpolationModel :
56  public BaseModel<1>
57  {
58 
59 public:
60  typedef double IntensityType;
62  typedef double CoordinateType;
64 
67  BaseModel<1>(),
68  interpolation_()
69  {
70  this->defaults_.setValue("interpolation_step", 0.1, "Sampling rate for the interpolation of the model function ");
71  this->defaults_.setValue("intensity_scaling", 1.0, "Scaling factor used to adjust the model distribution to the intensities of the data");
72  }
73 
76  BaseModel<1>(source),
77  interpolation_(source.interpolation_),
78  interpolation_step_(source.interpolation_step_),
79  scaling_(source.scaling_)
80  {
81  }
82 
85  {
86  }
87 
90  {
91  if (&source == this) return *this;
92 
94  interpolation_step_ = source.interpolation_step_;
95  interpolation_ = source.interpolation_;
96  scaling_ = source.scaling_;
97 
98  return *this;
99  }
100 
102  IntensityType getIntensity(const PositionType & pos) const
103  {
104  return interpolation_.value(pos[0]);
105  }
106 
108  IntensityType getIntensity(CoordinateType coord) const
109  {
110  return interpolation_.value(coord);
111  }
112 
114  const LinearInterpolation & getInterpolation() const
115  {
116  return interpolation_;
117  }
118 
124  CoordinateType getScalingFactor() const
125  {
126  return scaling_;
127  }
128 
134  virtual void setOffset(CoordinateType offset)
135  {
136  interpolation_.setOffset(offset);
137  }
138 
140  void getSamples(SamplesType & cont) const
141  {
142  cont = SamplesType();
144  for (Size i = 0; i < interpolation_.getData().size(); ++i)
145  {
146 #pragma clang diagnostic push
147 #pragma clang diagnostic ignored "-Wconversion"
148  peak.setIntensity(interpolation_.getData()[i]);
149 #pragma clang diagnostic pop
150  peak.getPosition()[0] = interpolation_.index2key(i);
151  cont.push_back(peak);
152  }
153  }
154 
156  virtual CoordinateType getCenter() const
157  {
158  throw Exception::NotImplemented(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION);
159  }
160 
162  virtual void setSamples()
163  {
164  throw Exception::NotImplemented(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION);
165  }
166 
172  void setInterpolationStep(CoordinateType interpolation_step)
173  {
174  interpolation_step_ = interpolation_step;
175  this->param_.setValue("interpolation_step", interpolation_step_);
176  }
177 
178  void setScalingFactor(CoordinateType scaling)
179  {
180  scaling_ = scaling;
181  this->param_.setValue("intensity_scaling", scaling_);
182  }
183 
184 protected:
185  LinearInterpolation interpolation_;
186  CoordinateType interpolation_step_;
187  CoordinateType scaling_;
188 
190  {
192  interpolation_step_ = this->param_.getValue("interpolation_step");
193  scaling_ = this->param_.getValue("intensity_scaling");
194  }
195 
196  };
197 }
198 
199 #endif // OPENMS_TRANSFORMATIONS_FEATUREFINDER_INTERPOLATIONMODEL_H
const LinearInterpolation & getInterpolation() const
Returns the interpolation class.
Definition: InterpolationModel.h:114
Math::LinearInterpolation< double > LinearInterpolation
Definition: InterpolationModel.h:63
double CoordinateType
Definition: InterpolationModel.h:62
InterpolationModel(const InterpolationModel &source)
copy constructor
Definition: InterpolationModel.h:75
void getSamples(SamplesType &cont) const
get reasonable set of samples from the model (i.e. for printing)
Definition: InterpolationModel.h:140
DPosition< 1 > PositionType
Definition: InterpolationModel.h:61
CoordinateType getScalingFactor() const
get the scaling for the model
Definition: InterpolationModel.h:124
virtual CoordinateType getCenter() const
"center" of the model, particular definition (depends on the derived model)
Definition: InterpolationModel.h:156
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
IntensityType getIntensity(const PositionType &pos) const
access model predicted intensity at position pos
Definition: InterpolationModel.h:102
LinearInterpolation interpolation_
Definition: InterpolationModel.h:185
double IntensityType
Definition: InterpolationModel.h:60
Abstract class for 1D-models that are approximated using linear interpolation.
Definition: InterpolationModel.h:55
CoordinateType scaling_
Definition: InterpolationModel.h:187
virtual void updateMembers_()
This method is used to update extra member variables at the end of the setParameters() method...
Definition: BaseModel.h:158
void setInterpolationStep(CoordinateType interpolation_step)
Set the interpolation step for the linear interpolation of the model.
Definition: InterpolationModel.h:172
InterpolationModel()
Default constructor.
Definition: InterpolationModel.h:66
virtual void setSamples()
set sample/supporting points of interpolation wrt params.
Definition: InterpolationModel.h:162
void setScalingFactor(CoordinateType scaling)
Definition: InterpolationModel.h:178
virtual ~InterpolationModel()
destructor
Definition: InterpolationModel.h:84
virtual InterpolationModel & operator=(const InterpolationModel &source)
assignment operator
Definition: InterpolationModel.h:89
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:128
virtual void setOffset(CoordinateType offset)
set the offset of the model
Definition: InterpolationModel.h:134
CoordinateType interpolation_step_
Definition: InterpolationModel.h:186
Abstract base class for all D-dimensional models.
Definition: BaseModel.h:51
void updateMembers_()
This method is used to update extra member variables at the end of the setParameters() method...
Definition: InterpolationModel.h:189
std::vector< PeakType > SamplesType
Definition: BaseModel.h:61
virtual BaseModel & operator=(const BaseModel &source)
assignment operator
Definition: BaseModel.h:84
IntensityType getIntensity(CoordinateType coord) const
access model predicted intensity at position pos
Definition: InterpolationModel.h:108
Not implemented exception.
Definition: Exception.h:437

OpenMS / TOPP release 2.3.0 Documentation generated on Tue Jan 9 2018 18:22:01 using doxygen 1.8.13