OpenMS  2.4.0
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-2018.
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 #pragma once
37 
40 
41 namespace OpenMS
42 {
54  class OPENMS_DLLAPI InterpolationModel :
55  public BaseModel<1>
56  {
57 
58 public:
59  typedef double IntensityType;
61  typedef double CoordinateType;
63 
66  BaseModel<1>(),
67  interpolation_()
68  {
69  this->defaults_.setValue("interpolation_step", 0.1, "Sampling rate for the interpolation of the model function ");
70  this->defaults_.setValue("intensity_scaling", 1.0, "Scaling factor used to adjust the model distribution to the intensities of the data");
71  }
72 
75  BaseModel<1>(source),
76  interpolation_(source.interpolation_),
77  interpolation_step_(source.interpolation_step_),
78  scaling_(source.scaling_)
79  {
80  }
81 
84  {
85  }
86 
89  {
90  if (&source == this) return *this;
91 
93  interpolation_step_ = source.interpolation_step_;
94  interpolation_ = source.interpolation_;
95  scaling_ = source.scaling_;
96 
97  return *this;
98  }
99 
101  IntensityType getIntensity(const PositionType & pos) const override
102  {
103  return interpolation_.value(pos[0]);
104  }
105 
108  {
109  return interpolation_.value(coord);
110  }
111 
114  {
115  return interpolation_;
116  }
117 
124  {
125  return scaling_;
126  }
127 
133  virtual void setOffset(CoordinateType offset)
134  {
135  interpolation_.setOffset(offset);
136  }
137 
139  void getSamples(SamplesType & cont) const override
140  {
141  cont = SamplesType();
143  for (Size i = 0; i < interpolation_.getData().size(); ++i)
144  {
145 #pragma clang diagnostic push
146 #pragma clang diagnostic ignored "-Wconversion"
147  peak.setIntensity(interpolation_.getData()[i]);
148 #pragma clang diagnostic pop
149  peak.getPosition()[0] = interpolation_.index2key(i);
150  cont.push_back(peak);
151  }
152  }
153 
155  virtual CoordinateType getCenter() const
156  {
157  throw Exception::NotImplemented(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION);
158  }
159 
161  virtual void setSamples()
162  {
163  throw Exception::NotImplemented(__FILE__, __LINE__, OPENMS_PRETTY_FUNCTION);
164  }
165 
171  void setInterpolationStep(CoordinateType interpolation_step)
172  {
173  interpolation_step_ = interpolation_step;
174  this->param_.setValue("interpolation_step", interpolation_step_);
175  }
176 
178  {
179  scaling_ = scaling;
180  this->param_.setValue("intensity_scaling", scaling_);
181  }
182 
183 protected:
187 
188  void updateMembers_() override
189  {
191  interpolation_step_ = this->param_.getValue("interpolation_step");
192  scaling_ = this->param_.getValue("intensity_scaling");
193  }
194 
195  };
196 }
197 
const LinearInterpolation & getInterpolation() const
Returns the interpolation class.
Definition: InterpolationModel.h:113
Math::LinearInterpolation< double > LinearInterpolation
Definition: InterpolationModel.h:62
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method...
Definition: InterpolationModel.h:188
void getSamples(SamplesType &cont) const override
get reasonable set of samples from the model (i.e. for printing)
Definition: InterpolationModel.h:139
~InterpolationModel() override
destructor
Definition: InterpolationModel.h:83
double CoordinateType
Definition: InterpolationModel.h:61
InterpolationModel(const InterpolationModel &source)
copy constructor
Definition: InterpolationModel.h:74
DPosition< 1 > PositionType
Definition: InterpolationModel.h:60
CoordinateType getScalingFactor() const
get the scaling for the model
Definition: InterpolationModel.h:123
virtual CoordinateType getCenter() const
"center" of the model, particular definition (depends on the derived model)
Definition: InterpolationModel.h:155
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
LinearInterpolation interpolation_
Definition: InterpolationModel.h:184
double IntensityType
Definition: InterpolationModel.h:59
Abstract class for 1D-models that are approximated using linear interpolation.
Definition: InterpolationModel.h:54
CoordinateType scaling_
Definition: InterpolationModel.h:186
void setInterpolationStep(CoordinateType interpolation_step)
Set the interpolation step for the linear interpolation of the model.
Definition: InterpolationModel.h:171
InterpolationModel()
Default constructor.
Definition: InterpolationModel.h:65
virtual void setSamples()
set sample/supporting points of interpolation wrt params.
Definition: InterpolationModel.h:161
void setScalingFactor(CoordinateType scaling)
Definition: InterpolationModel.h:177
IntensityType getIntensity(const PositionType &pos) const override
access model predicted intensity at position pos
Definition: InterpolationModel.h:101
virtual InterpolationModel & operator=(const InterpolationModel &source)
assignment operator
Definition: InterpolationModel.h:88
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
virtual void setOffset(CoordinateType offset)
set the offset of the model
Definition: InterpolationModel.h:133
CoordinateType interpolation_step_
Definition: InterpolationModel.h:185
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method...
Definition: BaseModel.h:157
Abstract base class for all D-dimensional models.
Definition: BaseModel.h:50
std::vector< PeakType > SamplesType
Definition: BaseModel.h:60
virtual BaseModel & operator=(const BaseModel &source)
assignment operator
Definition: BaseModel.h:83
IntensityType getIntensity(CoordinateType coord) const
access model predicted intensity at position pos
Definition: InterpolationModel.h:107
Not implemented exception.
Definition: Exception.h:436