Belle II Software development
ChebFitter.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9
10#include <cmath>
11#include <iostream>
12#include <sstream>
13
14#include <Math/Minimizer.h>
15#include <Math/Factory.h>
16#include <Math/Functor.h>
17
18#include <Eigen/Core>
19
20//If compiled within BASF2
21#ifdef _PACKAGE_
22#include <tracking/calibration/ChebFitter.h>
23#include <tracking/calibration/nodes.h>
24#include <framework/logging/Logger.h>
25#else
26#include <ChebFitter.h>
27#include <nodes.h>
28#define B2INFO(arg) { std::cout << arg << std::endl;}
29#endif
30
31
32namespace Belle2 {
38 using Eigen::VectorXd;
39 using Eigen::MatrixXd;
40
41
42//return values of -2*log(p(x)), where p(x) is normalized to 1 over the fitted range
43 VectorXd ChebFitter::getLogFunction(Pars pars) const
44 {
45 VectorXd fVals(m_nodes.size());
46
47 //calc function values
48 fVals = m_nodes.unaryExpr([&](double x) { return m_myFun(x, pars); });
49
50 //normalize the function
51 double I = fVals.dot(m_weights);
52
53 fVals = -2 * log(fVals.array() / I); //normalize by integral
54
55 return fVals;
56
57 }
58
59 void ChebFitter::init(int Size, double xMin, double xMax)
60 {
61 // loading the Cheb nodes
62 m_nodes = (xMax - xMin) * getNodes(Size).array() + xMin;
63
64 // loding the weights for integration
65 m_weights = (xMax - xMin) * getWeights(Size);
66
67
68 // calculate the transformation matrix from pol coefs to grid points
69 m_coefsMat = getCoefsCheb(Size).transpose();
70
72
73 }
74
75//function assumed to be normalized !!!
76 double ChebFitter::getLogLikelihoodSlow(const Pars& pars) const
77 {
78
79 double L = 0;
80 for (double d : m_data) {
81 double v = m_myFun(d, pars);
82 L += -2 * log(v);
83 }
84
85 return L;
86 }
87
88//evaluation using cheb pols
89 double ChebFitter::getLogLikelihoodFast(const Pars& pars) const
90 {
91 VectorXd funVals = getLogFunction(pars);
92 double LL = funVals.dot(m_dataGrid);
93
94 return LL;
95 }
96
97 double ChebFitter::operator()(const double* par) const
98 {
99 Pars pars;
100 for (unsigned i = 0; i < m_parNames.size(); ++i)
101 pars[m_parNames[i]] = par[i];
102
103 double LL = m_useCheb ? getLogLikelihoodFast(pars) : getLogLikelihoodSlow(pars);
104
105 return LL;
106 }
107
108// get data transformed into the grid such that (chebFunVals dot m_dataGrid) == logL
109 VectorXd ChebFitter::getDataGrid() const
110 {
111 double a = m_nodes[0];
112 double b = m_nodes[m_nodes.size() - 1];
113
114
115 VectorXd polSum = VectorXd::Zero(m_nodes.size());
116 for (double x : m_data) {
117 double xx = (x - a) / (b - a); //normalize between 0 and 1
118 polSum += getPols(m_nodes.size(), xx);
119 }
120
121
122 //transform to the basis of the cheb m_nodes
123 VectorXd gridVals = m_coefsMat * polSum;
124
125 return gridVals;
126 }
127
128
129// get data transformed into the grid such that (chebFunVals dot m_dataGrid) == logL
130 std::pair<VectorXd, MatrixXd> ChebFitter::getDataGridWithCov() const
131 {
132 double a = m_nodes[0];
133 double b = m_nodes[m_nodes.size() - 1];
134
135 MatrixXd polSum2 = MatrixXd::Zero(m_nodes.size(), m_nodes.size());
136 VectorXd polSum = VectorXd::Zero(m_nodes.size());
137 for (double x : m_data) {
138 double xx = (x - a) / (b - a); //normalize between 0 and 1
139 VectorXd pol = getPols(m_nodes.size(), xx);
140 polSum += pol;
141 polSum2 += pol * pol.transpose();
142 }
143
144
145 //transform to the basis of the cheb nodes
146 MatrixXd coefs = getCoefsCheb(polSum.size()).transpose();
147 VectorXd gridVals = coefs * polSum;
148 MatrixXd gridValsCov = coefs * polSum2 * coefs.transpose();
149
150 return std::make_pair(gridVals, gridValsCov);
151 }
152
153
154//Minimize using ROOT minimizer
155 std::pair<Pars, MatrixXd> ChebFitter::fitData(Pars pars, Limits limits, bool UseCheb)
156 {
157 m_useCheb = UseCheb;
158
159 ROOT::Math::Minimizer* minimum =
160 ROOT::Math::Factory::CreateMinimizer("Minuit2", "");
161
162 // set tolerance , etc...
163 minimum->SetMaxFunctionCalls(10000000); // for Minuit/Minuit2
164 minimum->SetMaxIterations(100000); // for GSL
165 minimum->SetTolerance(10.0);
166 //minimum->SetPrecision(1e-5);
167
168 //minimum->SetPrintLevel(3); //many outputs
169 minimum->SetPrintLevel(0); //few outputs
170 minimum->SetStrategy(2);
171 minimum->SetErrorDef(1);
172
173
174 // Set the free variables to be minimized !
175 m_parNames.clear();
176 int k = 0;
177 for (auto p : pars) {
178 std::string n = p.first;
179 double vCnt = p.second;
180 if (limits.count(n) == 1) {
181 double vMin = limits.at(n).first;
182 double vMax = limits.at(n).second;
183 double step = (vMax - vMin) / 100;
184 minimum->SetLimitedVariable(k, n, vCnt, step, vMin, vMax);
185 } else {
186 double step = 1;
187 minimum->SetVariable(k, n, vCnt, step);
188 }
189 m_parNames.push_back(n);
190 ++k;
191 }
192
193 // create function wrapper for minimizer
194 // a IMultiGenFunction type
195 ROOT::Math::Functor f(*this, pars.size());
196 minimum->SetFunction(f);
197
198
199 // do the minimization
200 minimum->Minimize();
201
202
203 Pars parsF;
204 for (unsigned i = 0; i < m_parNames.size(); ++i)
205 parsF[m_parNames[i]] = minimum->X()[i];
206
207 MatrixXd covMat(parsF.size(), parsF.size());
208 for (unsigned i = 0; i < parsF.size(); ++i)
209 for (unsigned j = 0; j < parsF.size(); ++j)
210 covMat(i, j) = minimum->CovMatrix(i, j);
211
212 // print pars
213 std::stringstream log;
214 log << "Minuit status : " << minimum->Status() << ", ";
215 for (auto p : parsF)
216 log << "\"" << p.first << "\" : " << p.second << ", ";
217
218 B2INFO(log.str());
219
220 delete minimum;
221
222 return std::make_pair(parsF, covMat);
223 }
224
225
226 double ChebFitter::getFunctionFast(const Pars& pars, double x)
227 {
228 static VectorXd funVals = getLogFunction(pars);
229
230
231 return exp(-0.5 * interpol(m_nodes, funVals, x));
232
233 }
234
236}
std::vector< double > m_data
vector with the data points to be fitted
Definition: ChebFitter.h:78
bool m_useCheb
flag to use approximation based on Chebyshev polynomials
Definition: ChebFitter.h:87
Eigen::VectorXd m_weights
vector with cheb weights for integration
Definition: ChebFitter.h:84
std::vector< std::string > m_parNames
vector with names of the parameters
Definition: ChebFitter.h:86
std::function< double(double, Pars)> m_myFun
function to fit
Definition: ChebFitter.h:89
Eigen::VectorXd m_dataGrid
vector with the data points related to the cheb nodes (m_dataGrid.size = nodes.size)
Definition: ChebFitter.h:79
Eigen::VectorXd m_nodes
vector with cheb nodes
Definition: ChebFitter.h:83
Eigen::MatrixXd m_coefsMat
transformation matrix from chebPol to gridPoints
Definition: ChebFitter.h:81
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.
Definition: ChebFitter.cc:155
double getFunctionFast(const Pars &pars, double x)
Evaluate the fitted function approximated with the Chebyshev polynomial, typically runs faster.
Definition: ChebFitter.cc:226
double operator()(const double *par) const
Evaluate the log likelihood.
Definition: ChebFitter.cc:97
double getLogLikelihoodSlow(const Pars &pars) const
Calculate log likelihood using exact formula.
Definition: ChebFitter.cc:76
Eigen::VectorXd getWeights(int Size)
Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definit...
Definition: nodes.cc:29
void init(int Size, double xMin, double xMax)
Initialize the fitter (the Chebyshev coefficients)
Definition: ChebFitter.cc:59
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....
Definition: nodes.cc:77
std::map< std::string, double > Pars
values of parameters in ML fit
Definition: ChebFitter.h:25
Eigen::VectorXd getDataGrid() const
Calculate Chebyshev coefficients for the data set.
Definition: ChebFitter.cc:109
Eigen::VectorXd getLogFunction(Pars pars) const
Get the -2*log(p(x)) on the Cheb nodes.
Definition: ChebFitter.cc:43
Eigen::MatrixXd getCoefsCheb(int oldSize)
Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borde...
Definition: nodes.cc:162
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getDataGridWithCov() const
Calculate Chebyshev coefficients with the covariance (useful for toy studies)
Definition: ChebFitter.cc:130
double getLogLikelihoodFast(const Pars &pars) const
Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)
Definition: ChebFitter.cc:89
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...
std::map< std::string, std::pair< double, double > > Limits
limits of parameters in ML fit
Definition: ChebFitter.h:28
Eigen::VectorXd getNodes(int Size)
Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1,...
Definition: nodes.cc:65
Abstract base class for different kinds of events.