Belle II Software development
ChebFitter.h
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#pragma once
11
12#include <vector>
13#include <map>
14#include <Eigen/Dense>
15#include <functional>
16
17
18namespace Belle2 {
25 typedef std::map<std::string, double> Pars;
26
28 typedef std::map<std::string, std::pair<double, double>> Limits;
29
30
32 class ChebFitter {
33
34 public:
35
37 void setDataAndFunction(std::function<double(double, Pars)> fun, const std::vector<double>& data)
38 {
39 m_myFun = fun;
40 m_data = data;
41 }
42
44 void init(int Size, double xMin, double xMax);
45
47 double operator()(const double* par) const;
48
49
51 std::pair<Pars, Eigen::MatrixXd> fitData(Pars pars, Limits limits, bool UseCheb = true);
52
53
54 private:
55
57 Eigen::VectorXd getLogFunction(Pars pars) const;
58
60 double getLogLikelihoodSlow(const Pars& pars) const;
61
63 double getLogLikelihoodFast(const Pars& pars) const;
64
65
67 Eigen::VectorXd getDataGrid() const;
68
70 std::pair<Eigen::VectorXd, Eigen::MatrixXd> getDataGridWithCov() const;
71
73 double getFunctionFast(const Pars& pars, double x);
74
75
76
77
78 std::vector<double> m_data;
79 Eigen::VectorXd m_dataGrid;
80 Eigen::MatrixXd m_dataGridCov;
81 Eigen::MatrixXd m_coefsMat;
82
83 Eigen::VectorXd m_nodes;
84 Eigen::VectorXd m_weights;
85
86 std::vector<std::string> m_parNames;
87 bool m_useCheb = true;
88
89 std::function<double(double, Pars)> m_myFun;
90
91 };
92
94}
Unbinned Maximum Likelihood fitter with a possibility to use Chebyshev interpolation.
Definition: ChebFitter.h:32
std::vector< double > m_data
vector with the data points to be fitted
Definition: ChebFitter.h:78
void setDataAndFunction(std::function< double(double, Pars)> fun, const std::vector< double > &data)
Set the fitted data and the fitted function, should be called before init.
Definition: ChebFitter.h:37
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::MatrixXd m_dataGridCov
vector with the data points covariances (useful for toy studies)
Definition: ChebFitter.h:80
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
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
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
std::map< std::string, std::pair< double, double > > Limits
limits of parameters in ML fit
Definition: ChebFitter.h:28
Abstract base class for different kinds of events.