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 #ifndef ROOT_Math_GSLMultiFit
00028 #define ROOT_Math_GSLMultiFit
00029
00030 #include "gsl/gsl_vector.h"
00031 #include "gsl/gsl_matrix.h"
00032 #include "gsl/gsl_multifit_nlin.h"
00033 #include "gsl/gsl_blas.h"
00034 #include "GSLMultiFitFunctionWrapper.h"
00035
00036 #include "Math/IFunction.h"
00037 #include <string>
00038 #include <cassert>
00039
00040
00041 namespace ROOT {
00042
00043 namespace Math {
00044
00045
00046
00047
00048
00049
00050
00051 class GSLMultiFit {
00052
00053 public:
00054
00055
00056
00057
00058
00059 GSLMultiFit () :
00060 fSolver(0),
00061 fVec(0),
00062 fCov(0),
00063 fType(gsl_multifit_fdfsolver_lmsder)
00064 {}
00065
00066
00067
00068
00069 ~GSLMultiFit () {
00070 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
00071 if (fVec != 0) gsl_vector_free(fVec);
00072 if (fCov != 0) gsl_matrix_free(fCov);
00073 }
00074
00075 private:
00076
00077
00078
00079
00080
00081 GSLMultiFit(const GSLMultiFit &) {}
00082
00083
00084
00085
00086 GSLMultiFit & operator = (const GSLMultiFit & rhs) {
00087 if (this == &rhs) return *this;
00088 return *this;
00089 }
00090
00091
00092 public:
00093
00094
00095 void CreateSolver(unsigned int npoints, unsigned int npar) {
00096 if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
00097 fSolver = gsl_multifit_fdfsolver_alloc(fType, npoints, npar);
00098 }
00099
00100
00101 template<class Func>
00102 int Set(const std::vector<Func> & funcVec, const double * x) {
00103
00104
00105 unsigned int npts = funcVec.size();
00106 if (npts == 0) return -1;
00107
00108 unsigned int npar = funcVec[0].NDim();
00109 typedef typename std::vector<Func> FuncVec;
00110
00111 fFunc.SetFunction(funcVec, npts, npar);
00112
00113 CreateSolver( npts, npar );
00114
00115 if (fVec != 0) gsl_vector_free(fVec);
00116 fVec = gsl_vector_alloc( npar );
00117 std::copy(x,x+npar, fVec->data);
00118 assert(fSolver != 0);
00119 return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
00120 }
00121
00122 std::string Name() const {
00123 if (fSolver == 0) return "undefined";
00124 return std::string(gsl_multifit_fdfsolver_name(fSolver) );
00125 }
00126
00127 int Iterate() {
00128 if (fSolver == 0) return -1;
00129 return gsl_multifit_fdfsolver_iterate(fSolver);
00130 }
00131
00132
00133 const double * X() const {
00134 if (fSolver == 0) return 0;
00135 gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
00136 return x->data;
00137 }
00138
00139
00140 const double * Gradient() const {
00141 if (fSolver == 0) return 0;
00142 gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
00143 return fVec->data;
00144 }
00145
00146
00147 const double * CovarMatrix() const {
00148 if (fSolver == 0) return 0;
00149 if (fCov != 0) gsl_matrix_free(fCov);
00150 unsigned int npar = fSolver->fdf->p;
00151 fCov = gsl_matrix_alloc( npar, npar );
00152 static double kEpsrel = 0.0001;
00153 int ret = gsl_multifit_covar(fSolver->J, kEpsrel, fCov);
00154 if (ret != GSL_SUCCESS) return 0;
00155 return fCov->data;
00156 }
00157
00158
00159 int TestGradient(double absTol) const {
00160 if (fSolver == 0) return -1;
00161 Gradient();
00162 return gsl_multifit_test_gradient( fVec, absTol);
00163 }
00164
00165
00166
00167 int TestDelta(double absTol, double relTol) const {
00168 if (fSolver == 0) return -1;
00169 return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
00170 }
00171
00172
00173 double Edm() const {
00174
00175 double edm = -1;
00176 const double * g = Gradient();
00177 if (g == 0) return edm;
00178 const double * c = CovarMatrix();
00179 if (c == 0) return edm;
00180 gsl_vector * tmp = gsl_vector_alloc( fSolver->fdf->p );
00181 int status = gsl_blas_dgemv(CblasNoTrans, 1.0, fCov, fVec, 0.,tmp);
00182 if (status == 0) status |= gsl_blas_ddot(fVec, tmp, &edm);
00183 gsl_vector_free(tmp);
00184 if (status != 0) return -1;
00185
00186 return 0.5*edm;
00187
00188 }
00189
00190
00191 private:
00192
00193 GSLMultiFitFunctionWrapper fFunc;
00194 gsl_multifit_fdfsolver * fSolver;
00195
00196 mutable gsl_vector * fVec;
00197 mutable gsl_matrix * fCov;
00198 const gsl_multifit_fdfsolver_type * fType;
00199
00200 };
00201
00202 }
00203
00204 }
00205
00206
00207 #endif