BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
ROOT::Math::GSLMultiFit Class Reference

Description

GSLMultiFit, internal class for implementing GSL non linear least square GSL fitting

Definition at line 52 of file GSLMultiFit.h.

Collaboration diagram for ROOT::Math::GSLMultiFit:
[legend]

Public Member Functions

 GSLMultiFit (const gsl_multifit_fdfsolver_type *type=0)
 
 ~GSLMultiFit ()
 
const double * CovarMatrix () const
 return covariance matrix of the parameters More...
 
void CreateSolver (unsigned int npoints, unsigned int npar)
 create the minimizer from the type and size of number of fitting points and number of parameters More...
 
double Edm () const
 
const double * Gradient () const
 gradient value at the minimum More...
 
int Iterate ()
 
std::string Name () const
 
template<class Func >
int Set (const std::vector< Func > &funcVec, const double *x)
 set the solver parameters More...
 
int TestDelta (double absTol, double relTol) const
 test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component More...
 
int TestGradient (double absTol) const
 test gradient (ask from solver gradient vector) More...
 
const double * X () const
 parameter values at the minimum More...
 

Private Member Functions

 GSLMultiFit (const GSLMultiFit &)
 
GSLMultiFitoperator= (const GSLMultiFit &rhs)
 

Private Attributes

gsl_matrix * fCov
 
GSLMultiFitFunctionWrapper fFunc
 
gsl_multifit_fdfsolver * fSolver
 
gsl_vector * fTmp
 
const gsl_multifit_fdfsolver_type * fType
 
gsl_vector * fVec
 

Constructor & Destructor Documentation

◆ GSLMultiFit() [1/2]

ROOT::Math::GSLMultiFit::GSLMultiFit ( const gsl_multifit_fdfsolver_type *  type = 0)
inline

Default constructor No need to specify the type so far since only one solver exists so far

Definition at line 60 of file GSLMultiFit.h.

60  :
61  fSolver(0),
62  fVec(0),
63  fTmp(0),
64  fCov(0),
65 #if GSL_MAJOR_VERSION > 1
66  fJac(0),
67 #endif
68  fType(type)
69  {
70  if (fType == 0) fType = gsl_multifit_fdfsolver_lmsder; // default value
71  }
gsl_multifit_fdfsolver * fSolver
Definition: GSLMultiFit.h:229
const gsl_multifit_fdfsolver_type * fType
Definition: GSLMultiFit.h:237

References fType.

◆ ~GSLMultiFit()

ROOT::Math::GSLMultiFit::~GSLMultiFit ( )
inline

Destructor (no operations)

Definition at line 76 of file GSLMultiFit.h.

76  {
77  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
78  if (fVec != 0) gsl_vector_free(fVec);
79  if (fTmp != 0) gsl_vector_free(fTmp);
80  if (fCov != 0) gsl_matrix_free(fCov);
81 #if GSL_MAJOR_VERSION > 1
82  if (fJac != 0) gsl_matrix_free(fJac);
83 #endif
84  }

References fCov, fSolver, fTmp, and fVec.

◆ GSLMultiFit() [2/2]

ROOT::Math::GSLMultiFit::GSLMultiFit ( const GSLMultiFit )
inlineprivate

Copy constructor

Definition at line 92 of file GSLMultiFit.h.

92 {}

Member Function Documentation

◆ CovarMatrix()

const double* ROOT::Math::GSLMultiFit::CovarMatrix ( ) const
inline

return covariance matrix of the parameters

Definition at line 181 of file GSLMultiFit.h.

181  {
182  if (fSolver == 0) return 0;
183  static double kEpsrel = 0.0001;
184 #if GSL_MAJOR_VERSION > 1
185  gsl_multifit_fdfsolver_jac (fSolver, fJac);
186  int ret = gsl_multifit_covar(fJac, kEpsrel, fCov);
187 #else
188  int ret = gsl_multifit_covar(fSolver->J, kEpsrel, fCov);
189 #endif
190  if (ret != GSL_SUCCESS) return 0;
191  return fCov->data;
192  }

References fCov, and fSolver.

Referenced by Edm().

◆ CreateSolver()

void ROOT::Math::GSLMultiFit::CreateSolver ( unsigned int  npoints,
unsigned int  npar 
)
inline

create the minimizer from the type and size of number of fitting points and number of parameters

Definition at line 106 of file GSLMultiFit.h.

106  {
107  if (fSolver) gsl_multifit_fdfsolver_free(fSolver);
108  fSolver = gsl_multifit_fdfsolver_alloc(fType, npoints, npar);
109  if (fVec != 0) gsl_vector_free(fVec);
110  fVec = gsl_vector_alloc( npar );
111  if (fTmp != 0) gsl_vector_free(fTmp);
112  fTmp = gsl_vector_alloc( npar );
113  if (fCov != 0) gsl_matrix_free(fCov);
114  fCov = gsl_matrix_alloc( npar, npar );
115 #if GSL_MAJOR_VERSION > 1
116  if (fJac != 0) gsl_matrix_free(fJac);
117  fJac = gsl_matrix_alloc( npoints, npar );
118 #endif
119  }

References fCov, fSolver, fTmp, fType, and fVec.

Referenced by Set().

◆ Edm()

double ROOT::Math::GSLMultiFit::Edm ( ) const
inline

Definition at line 209 of file GSLMultiFit.h.

209  {
210  // product C * g
211  double edm = -1;
212  const double * g = Gradient();
213  if (g == 0) return edm;
214  const double * c = CovarMatrix();
215  if (c == 0) return edm;
216  if (fTmp == 0) return edm;
217  int status = gsl_blas_dgemv(CblasNoTrans, 1.0, fCov, fVec, 0.,fTmp);
218  if (status == 0) status |= gsl_blas_ddot(fVec, fTmp, &edm);
219  if (status != 0) return -1;
220  // need to divide by 2 ??
221  return 0.5*edm;
222 
223  }
const double * Gradient() const
gradient value at the minimum
Definition: GSLMultiFit.h:170
const double * CovarMatrix() const
return covariance matrix of the parameters
Definition: GSLMultiFit.h:181

References CovarMatrix(), fCov, fTmp, fVec, and Gradient().

Here is the call graph for this function:

◆ Gradient()

const double* ROOT::Math::GSLMultiFit::Gradient ( ) const
inline

gradient value at the minimum

Definition at line 170 of file GSLMultiFit.h.

170  {
171  if (fSolver == 0) return 0;
172 #if GSL_MAJOR_VERSION > 1
173  fType->gradient(fSolver->state, fVec);
174 #else
175  gsl_multifit_gradient(fSolver->J, fSolver->f,fVec);
176 #endif
177  return fVec->data;
178  }

References fSolver, fType, and fVec.

Referenced by Edm(), and TestGradient().

◆ Iterate()

int ROOT::Math::GSLMultiFit::Iterate ( )
inline

Definition at line 157 of file GSLMultiFit.h.

157  {
158  if (fSolver == 0) return -1;
159  return gsl_multifit_fdfsolver_iterate(fSolver);
160  }

References fSolver.

◆ Name()

std::string ROOT::Math::GSLMultiFit::Name ( ) const
inline

Definition at line 152 of file GSLMultiFit.h.

152  {
153  if (fSolver == 0) return "undefined";
154  return std::string(gsl_multifit_fdfsolver_name(fSolver) );
155  }

References fSolver.

◆ operator=()

GSLMultiFit& ROOT::Math::GSLMultiFit::operator= ( const GSLMultiFit rhs)
inlineprivate

Assignment operator

Definition at line 97 of file GSLMultiFit.h.

97  {
98  if (this == &rhs) return *this; // time saving self-test
99  return *this;
100  }

◆ Set()

template<class Func >
int ROOT::Math::GSLMultiFit::Set ( const std::vector< Func > &  funcVec,
const double *  x 
)
inline

set the solver parameters

Definition at line 123 of file GSLMultiFit.h.

123  {
124  // create a vector of the fit contributions
125  // create function wrapper from an iterator of functions
126  unsigned int npts = funcVec.size();
127  if (npts == 0) return -1;
128 
129  unsigned int npar = funcVec[0].NDim();
130  // Remove unused typedef to remove warning in GCC48
131  // http://gcc.gnu.org/gcc-4.8/porting_to.html
132  // typedef typename std::vector<Func> FuncVec;
133  //FuncIt funcIter = funcVec.begin();
134  fFunc.SetFunction(funcVec, npts, npar);
135  // create solver object
136  CreateSolver( npts, npar );
137  assert(fSolver != 0);
138  // set initial values
139  assert(fVec != 0);
140  std::copy(x,x+npar, fVec->data);
141  assert(fTmp != 0);
142  gsl_vector_set_zero(fTmp);
143  assert(fCov != 0);
144  gsl_matrix_set_zero(fCov);
145 #if GSL_MAJOR_VERSION > 1
146  assert(fJac != 0);
147  gsl_matrix_set_zero(fJac);
148 #endif
149  return gsl_multifit_fdfsolver_set(fSolver, fFunc.GetFunc(), fVec);
150  }
void SetFunction(const FuncVector &f, unsigned int nres, unsigned int npar)
Fill gsl function structure from a C++ function iterator and size and number of residuals.
void CreateSolver(unsigned int npoints, unsigned int npar)
create the minimizer from the type and size of number of fitting points and number of parameters
Definition: GSLMultiFit.h:106
GSLMultiFitFunctionWrapper fFunc
Definition: GSLMultiFit.h:228

References CreateSolver(), fCov, fFunc, fSolver, fTmp, fVec, ROOT::Math::GSLMultiFitFunctionWrapper::GetFunc(), and ROOT::Math::GSLMultiFitFunctionWrapper::SetFunction().

Here is the call graph for this function:

◆ TestDelta()

int ROOT::Math::GSLMultiFit::TestDelta ( double  absTol,
double  relTol 
) const
inline

test using abs and relative tolerance |dx| < absTol + relTol*|x| for every component

Definition at line 203 of file GSLMultiFit.h.

203  {
204  if (fSolver == 0) return -1;
205  return gsl_multifit_test_delta(fSolver->dx, fSolver->x, absTol, relTol);
206  }

References fSolver.

◆ TestGradient()

int ROOT::Math::GSLMultiFit::TestGradient ( double  absTol) const
inline

test gradient (ask from solver gradient vector)

Definition at line 195 of file GSLMultiFit.h.

195  {
196  if (fSolver == 0) return -1;
197  Gradient();
198  return gsl_multifit_test_gradient( fVec, absTol);
199  }

References fSolver, fVec, and Gradient().

Here is the call graph for this function:

◆ X()

const double* ROOT::Math::GSLMultiFit::X ( ) const
inline

parameter values at the minimum

Definition at line 163 of file GSLMultiFit.h.

163  {
164  if (fSolver == 0) return 0;
165  gsl_vector * x = gsl_multifit_fdfsolver_position(fSolver);
166  return x->data;
167  }

References fSolver.

Member Data Documentation

◆ fCov

gsl_matrix* ROOT::Math::GSLMultiFit::fCov
mutableprivate

Definition at line 233 of file GSLMultiFit.h.

Referenced by ~GSLMultiFit(), CovarMatrix(), CreateSolver(), Edm(), and Set().

◆ fFunc

GSLMultiFitFunctionWrapper ROOT::Math::GSLMultiFit::fFunc
private

Definition at line 228 of file GSLMultiFit.h.

Referenced by Set().

◆ fSolver

gsl_multifit_fdfsolver* ROOT::Math::GSLMultiFit::fSolver
private

◆ fTmp

gsl_vector* ROOT::Math::GSLMultiFit::fTmp
mutableprivate

Definition at line 232 of file GSLMultiFit.h.

Referenced by ~GSLMultiFit(), CreateSolver(), Edm(), and Set().

◆ fType

const gsl_multifit_fdfsolver_type* ROOT::Math::GSLMultiFit::fType
private

Definition at line 237 of file GSLMultiFit.h.

Referenced by GSLMultiFit(), CreateSolver(), and Gradient().

◆ fVec

gsl_vector* ROOT::Math::GSLMultiFit::fVec
mutableprivate

Definition at line 231 of file GSLMultiFit.h.

Referenced by ~GSLMultiFit(), CreateSolver(), Edm(), Gradient(), Set(), and TestGradient().


The documentation for this class was generated from the following file: