00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef FIT_H
00030 #define FIT_H
00031
00032 #include <QObject>
00033
00034 #include "ApplicationWindow.h"
00035 #include "Filter.h"
00036
00037 #include <gsl/gsl_multifit_nlin.h>
00038 #include <gsl/gsl_multimin.h>
00039
00040 class Table;
00041 class Matrix;
00042
00044 class Fit : public Filter
00045 {
00046 Q_OBJECT
00047
00048 public:
00049
00050 typedef double (*fit_function_simplex)(const gsl_vector *, void *);
00051 typedef int (*fit_function)(const gsl_vector *, void *, gsl_vector *);
00052 typedef int (*fit_function_df)(const gsl_vector *, void *, gsl_matrix *);
00053 typedef int (*fit_function_fdf)(const gsl_vector *, void *, gsl_vector *, gsl_matrix *);
00054
00055 enum Algorithm{ScaledLevenbergMarquardt, UnscaledLevenbergMarquardt, NelderMeadSimplex};
00056 enum WeightingMethod{NoWeighting, Instrumental, Statistical, Dataset};
00057 enum FitType{BuiltIn = 0, Plugin = 1, User = 2};
00058
00059 Fit(ApplicationWindow *parent, Graph *g = 0, const QString& name = QString());
00060 Fit(ApplicationWindow *parent, Table *t, const QString& name = QString());
00061 ~Fit();
00062
00064 virtual void fit();
00065 virtual bool run(){fit(); return true;};
00066
00068 bool setWeightingData(WeightingMethod w, const QString& colName = QString::null);
00069
00070 void setDataCurve(int curve, double start, double end);
00071 bool setDataFromTable(Table *t, const QString& xColName, const QString& yColName, int from = 1, int to = -1);
00072
00073 QString resultFormula(){return d_result_formula;};
00074 QString formula(){return d_formula;};
00075 virtual void setFormula(const QString&){};
00076
00077 int numParameters(){return d_p;};
00078 QStringList parameterNames(){return d_param_names;};
00079 virtual void setParametersList(const QStringList&){};
00080 void setParameterExplanations(const QStringList& lst){d_param_explain = lst;};
00081
00082 double initialGuess(int parIndex){return gsl_vector_get(d_param_init, parIndex);};
00083 void setInitialGuess(int parIndex, double val){gsl_vector_set(d_param_init, parIndex, val);};
00084 void setInitialGuesses(double *x_init);
00085
00086 virtual void guessInitialValues(){};
00087
00088 void setParameterRange(int parIndex, double left, double right);
00089 void setAlgorithm(Algorithm s){d_solver = s;};
00090
00092 void generateFunction(bool yes, int points = 100);
00093
00095 virtual QString legendInfo();
00096
00098 double* results(){return d_results;};
00099
00101 double* errors();
00102
00104 double chiSquare() {return chi_2;};
00105
00107 double rSquare();
00108
00110 void scaleErrors(bool yes = true){d_scale_errors = yes;};
00111
00112 Table* parametersTable(const QString& tableName);
00113 void writeParametersToTable(Table *t, bool append = false);
00114
00115 Matrix* covarianceMatrix(const QString& matrixName);
00116
00117 bool save(const QString& fileName);
00118 bool load(const QString& fileName);
00119
00120 FitType type(){return d_fit_type;};
00121 void setType(FitType t){d_fit_type = t;};
00122
00123 QString fileName(){return d_file_name;};
00124 void setFileName(const QString& fn){d_file_name = fn;};
00125
00127 void freeMemory();
00128
00130 virtual double eval(double *, double){return 0.0;};
00131
00132 private:
00133 void init();
00134
00136 gsl_multimin_fminimizer * fitSimplex(gsl_multimin_function f, int &iterations, int &status);
00137
00139 gsl_multifit_fdfsolver * fitGSL(gsl_multifit_function_fdf f, int &iterations, int &status);
00140
00142 virtual void customizeFitResults(){};
00143
00144 protected:
00146 void initWorkspace(int par);
00148 void freeWorkspace();
00150 void insertFitFunctionCurve(const QString& name, double *x, double *y, int penWidth = 1);
00151
00153 virtual void generateFitCurve();
00154
00156 virtual void calculateFitCurveData(double *X, double *Y) {Q_UNUSED(X) Q_UNUSED(Y)};
00157
00159 virtual QString logFitInfo(int iterations, int status);
00160
00161 fit_function d_f;
00162 fit_function_df d_df;
00163 fit_function_fdf d_fdf;
00164 fit_function_simplex d_fsimplex;
00165
00167 int d_p;
00168
00170 gsl_vector *d_param_init;
00171
00175 bool is_non_linear;
00176
00178 double *d_w;
00179
00181 QStringList d_param_names;
00182
00184 QStringList d_param_explain;
00185
00187 bool d_gen_function;
00188
00190 Algorithm d_solver;
00191
00193 QString d_formula;
00194
00196 QString d_result_formula;
00197
00199 gsl_matrix *covar;
00200
00202 WeightingMethod d_weighting;
00203
00205 QString weighting_dataset;
00206
00208 double *d_results;
00209
00211 double *d_errors;
00212
00214 double chi_2;
00215
00217 bool d_scale_errors;
00218
00220 Table *d_param_table;
00221
00223 Matrix *d_cov_matrix;
00224
00225 FitType d_fit_type;
00226
00228 QString d_file_name;
00229
00231 double *d_param_range_left;
00232
00234 double *d_param_range_right;
00235 };
00236
00237 #endif