Belle II Software  release-06-01-15
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 
32 namespace Belle2 {
38  using Eigen::VectorXd;
39  using Eigen::MatrixXd;
40  using namespace std;
41 
42 
43 //return values of -2*log(p(x)), where p(x) is normalized to 1 over the fitted range
44  VectorXd ChebFitter::getLogFunction(Pars pars) const
45  {
46  VectorXd fVals(m_nodes.size());
47 
48  //calc function values
49  fVals = m_nodes.unaryExpr([&](double x) { return m_myFun(x, pars); });
50 
51  //normalize the function
52  double I = fVals.dot(m_weights);
53 
54  fVals = -2 * log(fVals.array() / I); //normalize by integral
55 
56  return fVals;
57 
58  }
59 
60  void ChebFitter::init(int Size, double xMin, double xMax)
61  {
62  // loading the Cheb nodes
63  m_nodes = (xMax - xMin) * getNodes(Size).array() + xMin;
64 
65  // loding the weights for integration
66  m_weights = (xMax - xMin) * getWeights(Size);
67 
68 
69  // calculate the transformation matrix from pol coefs to grid points
70  m_coefsMat = getCoefsCheb(Size).transpose();
71 
72  m_dataGrid = getDataGrid();
73 
74  }
75 
76 //function assumed to be normalized !!!
77  double ChebFitter::getLogLikelihoodSlow(const Pars& pars) const
78  {
79 
80  double L = 0;
81  for (double d : m_data) {
82  double v = m_myFun(d, pars);
83  L += -2 * log(v);
84  }
85 
86  return L;
87  }
88 
89 //evaluation using cheb pols
90  double ChebFitter::getLogLikelihoodFast(const Pars& pars) const
91  {
92  VectorXd funVals = getLogFunction(pars);
93  double LL = funVals.dot(m_dataGrid);
94 
95  return LL;
96  }
97 
98  double ChebFitter::operator()(const double* par) const
99  {
100  Pars pars;
101  for (unsigned i = 0; i < m_parNames.size(); ++i)
102  pars[m_parNames[i]] = par[i];
103 
104  double LL = m_useCheb ? getLogLikelihoodFast(pars) : getLogLikelihoodSlow(pars);
105 
106  return LL;
107  }
108 
109 // get data transformed into the grid such that (chebFunVals dot m_dataGrid) == logL
110  VectorXd ChebFitter::getDataGrid() const
111  {
112  double a = m_nodes[0];
113  double b = m_nodes[m_nodes.size() - 1];
114 
115 
116  VectorXd polSum = VectorXd::Zero(m_nodes.size());
117  for (double x : m_data) {
118  double xx = (x - a) / (b - a); //normalize between 0 and 1
119  polSum += getPols(m_nodes.size(), xx);
120  }
121 
122 
123  //transform to the basis of the cheb m_nodes
124  VectorXd gridVals = m_coefsMat * polSum;
125 
126  return gridVals;
127  }
128 
129 
130 // get data transformed into the grid such that (chebFunVals dot m_dataGrid) == logL
131  pair<VectorXd, MatrixXd> ChebFitter::getDataGridWithCov() const
132  {
133  double a = m_nodes[0];
134  double b = m_nodes[m_nodes.size() - 1];
135 
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); //normalize between 0 and 1
140  VectorXd pol = getPols(m_nodes.size(), xx);
141  polSum += pol;
142  polSum2 += pol * pol.transpose();
143  }
144 
145 
146  //transform to the basis of the cheb nodes
147  MatrixXd coefs = getCoefsCheb(polSum.size()).transpose();
148  VectorXd gridVals = coefs * polSum;
149  MatrixXd gridValsCov = coefs * polSum2 * coefs.transpose();
150 
151  return make_pair(gridVals, gridValsCov);
152  }
153 
154 
155 //Minimize using ROOT minimizer
156  pair<Pars, MatrixXd> ChebFitter::fitData(Pars pars, Limits limits, bool UseCheb)
157  {
158  m_useCheb = UseCheb;
159 
160  ROOT::Math::Minimizer* minimum =
161  ROOT::Math::Factory::CreateMinimizer("Minuit2", "");
162 
163  // set tolerance , etc...
164  minimum->SetMaxFunctionCalls(10000000); // for Minuit/Minuit2
165  minimum->SetMaxIterations(100000); // for GSL
166  minimum->SetTolerance(10.0);
167  //minimum->SetPrecision(1e-5);
168 
169  //minimum->SetPrintLevel(3); //many outputs
170  minimum->SetPrintLevel(0); //few outputs
171  minimum->SetStrategy(2);
172  minimum->SetErrorDef(1);
173 
174 
175  // Set the free variables to be minimized !
176  m_parNames.clear();
177  int k = 0;
178  for (auto p : pars) {
179  string n = p.first;
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);
186  } else {
187  double step = 1;
188  minimum->SetVariable(k, n, vCnt, step);
189  }
190  m_parNames.push_back(n);
191  ++k;
192  }
193 
194  // create function wrapper for minimizer
195  // a IMultiGenFunction type
196  ROOT::Math::Functor f(*this, pars.size());
197  minimum->SetFunction(f);
198 
199 
200  // do the minimization
201  minimum->Minimize();
202 
203 
204  Pars parsF;
205  for (unsigned i = 0; i < m_parNames.size(); ++i)
206  parsF[m_parNames[i]] = minimum->X()[i];
207 
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);
212 
213  // print pars
214  stringstream log;
215  log << "Minuit status : " << minimum->Status() << ", ";
216  for (auto p : parsF)
217  log << "\"" << p.first << "\" : " << p.second << ", ";
218 
219  B2INFO(log.str());
220 
221  delete minimum;
222 
223  return make_pair(parsF, covMat);
224  }
225 
226 
227  double ChebFitter::getFunctionFast(const Pars& pars, double x)
228  {
229  static VectorXd funVals = getLogFunction(pars);
230 
231 
232  return exp(-0.5 * interpol(m_nodes, funVals, x));
233 
234  }
235 
237 }
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.
Definition: ChebFitter.cc:227
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getDataGridWithCov() const
Calculate Chebyshev coefficients with the covariance (useful for toy studies)
Definition: ChebFitter.cc:131
double operator()(const double *par) const
Evaluate the log likelihood.
Definition: ChebFitter.cc:98
double getLogLikelihoodSlow(const Pars &pars) const
Calculate log likelihood using exact formula.
Definition: ChebFitter.cc:77
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:30
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:156
void init(int Size, double xMin, double xMax)
Initialize the fitter (the Chebyshev coefficients)
Definition: ChebFitter.cc:60
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:78
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:110
Eigen::VectorXd getLogFunction(Pars pars) const
Get the -2*log(p(x)) on the Cheb nodes.
Definition: ChebFitter.cc:44
Eigen::MatrixXd getCoefsCheb(int oldSize)
Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borde...
Definition: nodes.cc:163
double getLogLikelihoodFast(const Pars &pars) const
Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)
Definition: ChebFitter.cc:90
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:66
Abstract base class for different kinds of events.