35 #ifndef OPENMS_TRANSFORMATIONS_RAW2PEAK_OPTIMIZEPICK_H 
   36 #define OPENMS_TRANSFORMATIONS_RAW2PEAK_OPTIMIZEPICK_H 
   41 #include <gsl/gsl_vector.h> 
   42 #include <gsl/gsl_multifit_nlin.h> 
   43 #include <gsl/gsl_blas.h> 
   56   namespace OptimizationFunctions
 
   74         pos(0), lWidth(0), rWidth(0) {}
 
   76         pos(p.pos), lWidth(p.lWidth), rWidth(p.rWidth) {}
 
   97     int residual(
const gsl_vector * x, 
void * params, gsl_vector * f);
 
  100     int jacobian(
const gsl_vector * x, 
void * params, gsl_matrix * J);
 
  103     int evaluate(
const gsl_vector * x, 
void * params, gsl_vector * f, gsl_matrix * J);
 
  106     void printSignal(
const gsl_vector * x, 
void * param, 
float resolution = 0.25);
 
  147                  const int max_iteration_,
 
  148                  const double eps_abs_,
 
  149                  const double eps_rel_);
 
  183     void optimize(std::vector<PeakShape> & peaks, Data & data);
 
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'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'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'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