14 #include <Math/Minimizer.h>
15 #include <Math/Factory.h>
16 #include <Math/Functor.h>
22 #include <tracking/calibration/ChebFitter.h>
23 #include <tracking/calibration/nodes.h>
24 #include <framework/logging/Logger.h>
26 #include <ChebFitter.h>
28 #define B2INFO(arg) { std::cout << arg << std::endl;}
38 using Eigen::VectorXd;
39 using Eigen::MatrixXd;
46 VectorXd fVals(m_nodes.size());
49 fVals = m_nodes.unaryExpr([&](
double x) {
return m_myFun(x, pars); });
52 double I = fVals.dot(m_weights);
54 fVals = -2 * log(fVals.array() / I);
63 m_nodes = (xMax - xMin) *
getNodes(Size).array() + xMin;
72 m_dataGrid = getDataGrid();
81 for (
double d : m_data) {
82 double v = m_myFun(d, pars);
92 VectorXd funVals = getLogFunction(pars);
93 double LL = funVals.dot(m_dataGrid);
101 for (
unsigned i = 0; i < m_parNames.size(); ++i)
102 pars[m_parNames[i]] = par[i];
104 double LL = m_useCheb ? getLogLikelihoodFast(pars) : getLogLikelihoodSlow(pars);
112 double a = m_nodes[0];
113 double b = m_nodes[m_nodes.size() - 1];
116 VectorXd polSum = VectorXd::Zero(m_nodes.size());
117 for (
double x : m_data) {
118 double xx = (x - a) / (b - a);
119 polSum +=
getPols(m_nodes.size(), xx);
124 VectorXd gridVals = m_coefsMat * polSum;
133 double a = m_nodes[0];
134 double b = m_nodes[m_nodes.size() - 1];
136 MatrixXd polSum2 = MatrixXd::Zero(m_nodes.size(), m_nodes.size());
137 VectorXd polSum = VectorXd::Zero(m_nodes.size());
138 for (
double x : m_data) {
139 double xx = (x - a) / (b - a);
140 VectorXd pol =
getPols(m_nodes.size(), xx);
142 polSum2 += pol * pol.transpose();
147 MatrixXd coefs =
getCoefsCheb(polSum.size()).transpose();
148 VectorXd gridVals = coefs * polSum;
149 MatrixXd gridValsCov = coefs * polSum2 * coefs.transpose();
151 return make_pair(gridVals, gridValsCov);
160 ROOT::Math::Minimizer* minimum =
161 ROOT::Math::Factory::CreateMinimizer(
"Minuit2",
"");
164 minimum->SetMaxFunctionCalls(10000000);
165 minimum->SetMaxIterations(100000);
166 minimum->SetTolerance(10.0);
170 minimum->SetPrintLevel(0);
171 minimum->SetStrategy(2);
172 minimum->SetErrorDef(1);
178 for (
auto p : pars) {
180 double vCnt = p.second;
181 if (limits.count(n) == 1) {
182 double vMin = limits.at(n).first;
183 double vMax = limits.at(n).second;
184 double step = (vMax - vMin) / 100;
185 minimum->SetLimitedVariable(k, n, vCnt, step, vMin, vMax);
188 minimum->SetVariable(k, n, vCnt, step);
190 m_parNames.push_back(n);
196 ROOT::Math::Functor f(*
this, pars.size());
197 minimum->SetFunction(f);
205 for (
unsigned i = 0; i < m_parNames.size(); ++i)
206 parsF[m_parNames[i]] = minimum->X()[i];
208 MatrixXd covMat(parsF.size(), parsF.size());
209 for (
unsigned i = 0; i < parsF.size(); ++i)
210 for (
unsigned j = 0; j < parsF.size(); ++j)
211 covMat(i, j) = minimum->CovMatrix(i, j);
215 log <<
"Minuit status : " << minimum->Status() <<
", ";
217 log <<
"\"" << p.first <<
"\" : " << p.second <<
", ";
223 return make_pair(parsF, covMat);
229 static VectorXd funVals = getLogFunction(pars);
232 return exp(-0.5 *
interpol(m_nodes, funVals, x));
Eigen::VectorXd interpol(const Eigen::VectorXd &xi, double x)
Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolat...
double getFunctionFast(const Pars &pars, double x)
Evaluate the fitted function approximated with the Chebyshev polynomial, typically runs faster.
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getDataGridWithCov() const
Calculate Chebyshev coefficients with the covariance (useful for toy studies)
double operator()(const double *par) const
Evaluate the log likelihood.
double getLogLikelihoodSlow(const Pars &pars) const
Calculate log likelihood using exact formula.
Eigen::VectorXd getWeights(int Size)
Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definit...
std::pair< Pars, Eigen::MatrixXd > fitData(Pars pars, Limits limits, bool UseCheb=true)
Fit the data with specified initial values of parameters and limits on them.
void init(int Size, double xMin, double xMax)
Initialize the fitter (the Chebyshev coefficients)
Eigen::VectorXd getPols(int Size, double x)
Evaluate Chebyshev polynomials up to Size at point x It returns a vector of the P_i(x) for i=0....
std::map< std::string, double > Pars
values of parameters in ML fit
Eigen::VectorXd getDataGrid() const
Calculate Chebyshev coefficients for the data set.
Eigen::VectorXd getLogFunction(Pars pars) const
Get the -2*log(p(x)) on the Cheb nodes.
Eigen::MatrixXd getCoefsCheb(int oldSize)
Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borde...
double getLogLikelihoodFast(const Pars &pars) const
Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)
std::map< std::string, std::pair< double, double > > Limits
limits of parameters in ML fit
Eigen::VectorXd getNodes(int Size)
Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1,...
Abstract base class for different kinds of events.