/usr/include/root/Math/GSLNLSMinimizer.h is in libroot-math-mathmore-dev 5.34.30-0ubuntu8.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 | // @(#)root/mathmore:$Id$
// Author: L. Moneta Wed Dec 20 17:16:32 2006
/**********************************************************************
* *
* Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
* *
* This library is free software; you can redistribute it and/or *
* modify it under the terms of the GNU General Public License *
* as published by the Free Software Foundation; either version 2 *
* of the License, or (at your option) any later version. *
* *
* This library is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
* General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this library (see file COPYING); if not, write *
* to the Free Software Foundation, Inc., 59 Temple Place, Suite *
* 330, Boston, MA 02111-1307 USA, or contact the author. *
* *
**********************************************************************/
// Header file for class GSLNLSMinimizer
#ifndef ROOT_Math_GSLNLSMinimizer
#define ROOT_Math_GSLNLSMinimizer
#ifndef ROOT_Math_BasicMinimizer
#include "Math/BasicMinimizer.h"
#endif
#ifndef ROOT_Math_IFunctionfwd
#include "Math/IFunctionfwd.h"
#endif
#ifndef ROOT_Math_IParamFunctionfwd
#include "Math/IParamFunctionfwd.h"
#endif
#ifndef ROOT_Math_FitMethodFunction
#include "Math/FitMethodFunction.h"
#endif
#ifndef ROOT_Math_MinimTransformVariable
#include "Math/MinimTransformVariable.h"
#endif
#include <vector>
#include <map>
#include <string>
namespace ROOT {
namespace Math {
class GSLMultiFit;
//________________________________________________________________________________
/**
LSResidualFunc class description.
Internal class used for accessing the residuals of the Least Square function
and their derivates which are estimated numerically using GSL numerical derivation.
The class contains a pointer to the fit method function and an index specifying
the i-th residual and wraps it in a multi-dim gradient function interface
ROOT::Math::IGradientFunctionMultiDim.
The class is used by ROOT::Math::GSLNLSMinimizer (GSL non linear least square fitter)
@ingroup MultiMin
*/
class LSResidualFunc : public IMultiGradFunction {
public:
//default ctor (required by CINT)
LSResidualFunc() : fIndex(0), fChi2(0)
{}
LSResidualFunc(const ROOT::Math::FitMethodFunction & func, unsigned int i) :
fIndex(i),
fChi2(&func),
fX2(std::vector<double>(func.NDim() ) )
{}
// copy ctor
LSResidualFunc(const LSResidualFunc & rhs) :
IMultiGenFunction(),
IMultiGradFunction()
{
operator=(rhs);
}
// assignment
LSResidualFunc & operator= (const LSResidualFunc & rhs)
{
fIndex = rhs.fIndex;
fChi2 = rhs.fChi2;
fX2 = rhs.fX2;
return *this;
}
IMultiGenFunction * Clone() const {
return new LSResidualFunc(*fChi2,fIndex);
}
unsigned int NDim() const { return fChi2->NDim(); }
void Gradient( const double * x, double * g) const {
double f0 = 0;
FdF(x,f0,g);
}
void FdF (const double * x, double & f, double * g) const {
unsigned int n = NDim();
std::copy(x,x+n,fX2.begin());
const double kEps = 1.0E-4;
f = DoEval(x);
for (unsigned int i = 0; i < n; ++i) {
fX2[i] += kEps;
g[i] = ( DoEval(&fX2.front()) - f )/kEps;
fX2[i] = x[i];
}
}
private:
double DoEval (const double * x) const {
return fChi2->DataElement(x, fIndex);
}
double DoDerivative(const double * x, unsigned int icoord) const {
//return ROOT::Math::Derivator::Eval(*this, x, icoord, 1E-8);
std::copy(x,x+NDim(),fX2.begin());
const double kEps = 1.0E-4;
fX2[icoord] += kEps;
return ( DoEval(&fX2.front()) - DoEval(x) )/kEps;
}
unsigned int fIndex;
const ROOT::Math::FitMethodFunction * fChi2;
mutable std::vector<double> fX2; // cached vector
};
//_____________________________________________________________________________________________________
/**
GSLNLSMinimizer class for Non Linear Least Square fitting
It Uses the Levemberg-Marquardt algorithm from
<A HREF="http://www.gnu.org/software/gsl/manual/html_node/Nonlinear-Least_002dSquares-Fitting.html">
GSL Non Linear Least Square fitting</A>.
@ingroup MultiMin
*/
class GSLNLSMinimizer : public ROOT::Math::BasicMinimizer {
public:
/**
Default constructor
*/
GSLNLSMinimizer (int type = 0);
/**
Destructor (no operations)
*/
~GSLNLSMinimizer ();
private:
// usually copying is non trivial, so we make this unaccessible
/**
Copy constructor
*/
GSLNLSMinimizer(const GSLNLSMinimizer &) : ROOT::Math::BasicMinimizer() {}
/**
Assignment operator
*/
GSLNLSMinimizer & operator = (const GSLNLSMinimizer & rhs) {
if (this == &rhs) return *this; // time saving self-test
return *this;
}
public:
/// set the function to minimize
virtual void SetFunction(const ROOT::Math::IMultiGenFunction & func);
/// set gradient the function to minimize
virtual void SetFunction(const ROOT::Math::IMultiGradFunction & func);
/// method to perform the minimization
virtual bool Minimize();
/// return expected distance reached from the minimum
virtual double Edm() const { return fEdm; } // not impl. }
/// return pointer to gradient values at the minimum
virtual const double * MinGradient() const;
/// number of function calls to reach the minimum
virtual unsigned int NCalls() const { return (fChi2Func) ? fChi2Func->NCalls() : 0; }
/// number of free variables (real dimension of the problem)
/// this is <= Function().NDim() which is the total
// virtual unsigned int NFree() const { return fNFree; }
/// minimizer provides error and error matrix
virtual bool ProvidesError() const { return true; }
/// return errors at the minimum
virtual const double * Errors() const { return (fErrors.size() > 0) ? &fErrors.front() : 0; }
// {
// static std::vector<double> err;
// err.resize(fDim);
// return &err.front();
// }
/** return covariance matrices elements
if the variable is fixed the matrix is zero
The ordering of the variables is the same as in errors
*/
virtual double CovMatrix(unsigned int , unsigned int ) const;
/// return covariance matrix status
virtual int CovMatrixStatus() const;
protected:
private:
unsigned int fNFree; // dimension of the internal function to be minimized
unsigned int fSize; // number of fit points (residuals)
ROOT::Math::GSLMultiFit * fGSLMultiFit; // pointer to GSL multi fit solver
const ROOT::Math::FitMethodFunction * fChi2Func; // pointer to Least square function
double fEdm; // edm value
double fLSTolerance; // Line Search Tolerance
std::vector<double> fErrors;
std::vector<double> fCovMatrix; // cov matrix (stored as cov[ i * dim + j]
std::vector<LSResidualFunc> fResiduals; //! transient Vector of the residual functions
};
} // end namespace Math
} // end namespace ROOT
#endif /* ROOT_Math_GSLNLSMinimizer */
|