Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
OptimizePick.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-2013.
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: Alexandra Zerck $
32 // $Authors: Eva Lange $
33 // --------------------------------------------------------------------------
34 
35 #ifndef OPENMS_TRANSFORMATIONS_RAW2PEAK_OPTIMIZEPICK_H
36 #define OPENMS_TRANSFORMATIONS_RAW2PEAK_OPTIMIZEPICK_H
37 
39 #include <OpenMS/KERNEL/Peak1D.h>
40 
41 #include <gsl/gsl_vector.h>
42 #include <gsl/gsl_multifit_nlin.h>
43 #include <gsl/gsl_blas.h>
44 
45 #include <iostream>
46 #include <fstream>
47 #include <vector>
48 
49 namespace OpenMS
50 {
56  namespace OptimizationFunctions
57  {
59  typedef std::vector<Peak1D> RawDataVector;
61  typedef RawDataVector::iterator PeakIterator;
62 
71  struct OPENMS_DLLAPI PenaltyFactors
72  {
74  pos(0), lWidth(0), rWidth(0) {}
76  pos(p.pos), lWidth(p.lWidth), rWidth(p.rWidth) {}
78  {
79  pos = p.pos;
80  lWidth = p.lWidth;
81  rWidth = p.rWidth;
82 
83  return *this;
84  }
85 
87 
89  double pos;
91  double lWidth;
93  double rWidth;
94  };
95 
97  int residual(const gsl_vector * x, void * params, gsl_vector * f);
98 
100  int jacobian(const gsl_vector * x, void * params, gsl_matrix * J);
101 
103  int evaluate(const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix * J);
104 
106  void printSignal(const gsl_vector * x, void * param, float resolution = 0.25);
107  }
108 
109 
116  class OPENMS_DLLAPI OptimizePick
117  {
118 public:
119 
120  struct Data
121  {
123  std::vector<double> positions;
124  std::vector<double> signal;
126  std::vector<PeakShape> peaks;
127 
129 
130  };
131 
132 
134  typedef std::vector<Peak1D> RawDataVector;
136  typedef RawDataVector::iterator PeakIterator;
137 
138 
141  max_iteration_(0),
142  eps_abs_(0),
143  eps_rel_(0) {}
144 
146  OptimizePick(const struct OptimizationFunctions::PenaltyFactors & penalties_,
147  const int max_iteration_,
148  const double eps_abs_,
149  const double eps_rel_);
150 
152  ~OptimizePick();
153 
155  inline const struct OptimizationFunctions::PenaltyFactors & getPenalties() const { return penalties_; }
157  inline struct OptimizationFunctions::PenaltyFactors & getPenalties() { return penalties_; }
159  inline void setPenalties(const struct OptimizationFunctions::PenaltyFactors & penalties) { penalties_ = penalties; }
160 
162  inline UInt getNumberIterations() const { return max_iteration_; }
164  inline unsigned int & getNumberIterations() { return max_iteration_; }
166  inline void setNumberIterations(const int max_iteration) { max_iteration_ = max_iteration; }
167 
169  inline DoubleReal getMaxAbsError() const { return eps_abs_; }
171  inline double & getMaxAbsError() { return eps_abs_; }
173  inline void setMaxAbsError(double eps_abs) { eps_abs_ = eps_abs; }
174 
176  inline DoubleReal getMaxRelError() const { return eps_rel_; }
178  inline double & getMaxRelError() { return eps_rel_; }
180  inline void setMaxRelError(double eps_rel) { eps_rel_ = eps_rel; }
181 
183  void optimize(std::vector<PeakShape> & peaks, Data & data);
184 
185 
186 protected:
189 
191  unsigned int max_iteration_;
192 
194  double eps_abs_;
195  double eps_rel_;
196 
197 // /** @brief Returns the squared pearson coefficient.
198 
199 // Computes the correlation of the peak and the original data given by the peak enpoints.
200 // If the value is near 1, the fitted peakshape and the raw data are expected to be very similar.
201 // */
202 // double correlate_(const PeakShape& peak,
203 // double left_endpoint,
204 // double right_endpoint,Data& data);
205 
206  };
207 }
208 
209 #endif
std::vector< double > positions
Positions and intensity values of the raw data.
Definition: OptimizePick.h:123
unsigned int max_iteration_
Maximum number of iterations during optimization.
Definition: OptimizePick.h:191
Definition: OptimizePick.h:120
void setMaxRelError(double eps_rel)
Mutable access to the maximum relative error.
Definition: OptimizePick.h:180
void setMaxAbsError(double eps_abs)
Mutable access to the maximum absolute error.
Definition: OptimizePick.h:173
void setNumberIterations(const int max_iteration)
Mutable access to the number of iterations.
Definition: OptimizePick.h:166
OptimizePick()
Constructor.
Definition: OptimizePick.h:140
~PenaltyFactors()
Definition: OptimizePick.h:86
DoubleReal getMaxAbsError() const
Non-mutable access to the maximum absolute error.
Definition: OptimizePick.h:169
PenaltyFactors & operator=(const PenaltyFactors &p)
Definition: OptimizePick.h:77
struct OptimizationFunctions::PenaltyFactors & getPenalties() const
Non-mutable access to the penalty factors.
Definition: OptimizePick.h:155
OptimizationFunctions::PenaltyFactors penalties
Definition: OptimizePick.h:128
PenaltyFactors(const PenaltyFactors &p)
Definition: OptimizePick.h:75
std::vector< Peak1D > RawDataVector
Raw data vector type.
Definition: OptimizePick.h:59
double lWidth
Penalty factor for the peak shape&#39;s left width parameter.
Definition: OptimizePick.h:91
DoubleReal getMaxRelError() const
Non-mutable access to the maximum relative error.
Definition: OptimizePick.h:176
RawDataVector::iterator PeakIterator
Raw data iterator type.
Definition: OptimizePick.h:61
std::vector< PeakShape > peaks
This container contains the peak shapes to be optimized.
Definition: OptimizePick.h:126
int evaluate(const gsl_vector *x, void *params, gsl_vector *f, gsl_matrix *J)
Driver function for the evaluation of function and jacobian.
double & getMaxRelError()
Mutable access to the maximum relative error.
Definition: OptimizePick.h:178
double eps_abs_
Maximum absolute and relative error used in the optimization.
Definition: OptimizePick.h:194
double & getMaxAbsError()
Mutable access to the maximum absolute error.
Definition: OptimizePick.h:171
int residual(const gsl_vector *x, void *params, gsl_vector *f)
Evaluation of the target function for nonlinear optimization.
double pos
Penalty factor for the peak shape&#39;s position.
Definition: OptimizePick.h:89
unsigned int & getNumberIterations()
Mutable access to the number of iterations.
Definition: OptimizePick.h:164
double rWidth
Penalty factor for the peak shape&#39;s right width parameter.
Definition: OptimizePick.h:93
int jacobian(const gsl_vector *x, void *params, gsl_matrix *J)
Compute the Jacobian of the residual, where each row of the matrix corresponds to a point in the data...
Class for the penalty factors used during the optimization.
Definition: OptimizePick.h:71
UInt getNumberIterations() const
Non-mutable access to the number of iterations.
Definition: OptimizePick.h:162
std::vector< double > signal
Definition: OptimizePick.h:124
RawDataVector::iterator PeakIterator
Raw data iterator type.
Definition: OptimizePick.h:136
double eps_rel_
Definition: OptimizePick.h:195
void printSignal(const gsl_vector *x, void *param, float resolution=0.25)
Print all peak shapes.
This class provides the non-linear optimization of the peak parameters.
Definition: OptimizePick.h:116
struct OptimizationFunctions::PenaltyFactors & getPenalties()
Mutable access to the penalty factors.
Definition: OptimizePick.h:157
void setPenalties(const struct OptimizationFunctions::PenaltyFactors &penalties)
Mutable access to the penalty factors.
Definition: OptimizePick.h:159
PenaltyFactors()
Definition: OptimizePick.h:73
std::vector< Peak1D > RawDataVector
Raw data vector type.
Definition: OptimizePick.h:134

OpenMS / TOPP release 1.11.1 Documentation generated on Thu Nov 14 2013 11:19:19 using doxygen 1.8.5