PhysicsTools
Utilities
interface
HistoChiSquare.h
Go to the documentation of this file.
1
#ifndef PhysicsTools_Utilities_HistoChiSquare_h
2
#define PhysicsTools_Utilities_HistoChiSquare_h
3
#include "
PhysicsTools/Utilities/interface/RootMinuitResultPrinter.h
"
4
#include <vector>
5
#include "TH1.h"
6
#include "TMath.h"
7
8
namespace
fit
{
9
template
<
typename
T>
10
class
HistoChiSquare
{
11
public
:
12
HistoChiSquare
() {}
13
HistoChiSquare
(
T
&
t
, TH1 *
histo
,
double
rangeMin,
double
rangeMax)
14
:
t_
(&
t
),
rangeMin_
(rangeMin),
rangeMax_
(rangeMax) {
15
nBins_
=
histo
->GetNbinsX();
16
xMin_
=
histo
->GetXaxis()->GetXmin();
17
xMax_
=
histo
->GetXaxis()->GetXmax();
18
deltaX_
= (
xMax_
-
xMin_
) /
nBins_
;
19
for
(
size_t
i
= 0;
i
<
nBins_
; ++
i
) {
20
cont_
.push_back(
histo
->GetBinContent(
i
+ 1));
21
err_
.push_back(
histo
->GetBinError(
i
+ 1));
22
}
23
}
24
double
operator()
()
const
{
25
double
chi2
= 0;
26
for
(
size_t
i
= 0;
i
<
nBins_
; ++
i
) {
27
double
x =
xMin_
+ (
i
+ .5) *
deltaX_
;
28
if
((x >
rangeMin_
) && (x <
rangeMax_
) && (
err_
[
i
] > 0)) {
29
double
r
= (
cont_
[
i
] - (*t_)(x)) /
err_
[
i
];
30
chi2
+= (
r
*
r
);
31
}
32
}
33
return
chi2
;
34
}
35
void
setHistos
(TH1 *
histo
) {
36
nBins_
=
histo
->GetNbinsX();
37
xMin_
=
histo
->GetXaxis()->GetXmin();
38
xMax_
=
histo
->GetXaxis()->GetXmax();
39
deltaX_
= (
xMax_
-
xMin_
) /
nBins_
;
40
}
41
size_t
numberOfBins
()
const
{
42
size_t
fullBins = 0;
43
for
(
size_t
i
= 0;
i
<
nBins_
; ++
i
) {
44
double
x =
xMin_
+ (
i
+ .5) *
deltaX_
;
45
if
((x >
rangeMin_
) && (x <
rangeMax_
) && (
err_
[
i
] > 0))
46
fullBins++;
47
}
48
return
fullBins;
49
}
50
T
&
function
() {
return
*
t_
; }
51
const
T
&
function
()
const
{
return
*
t_
; }
52
53
private
:
54
T
*
t_
;
55
double
rangeMin_
,
rangeMax_
;
56
size_t
nBins_
;
57
double
xMin_
,
xMax_
,
deltaX_
;
58
std::vector<double>
cont_
;
59
std::vector<double>
err_
;
60
};
61
62
template
<
typename
T>
63
struct
RootMinuitResultPrinter
<
HistoChiSquare
<
T
> > {
64
static
void
print
(
double
amin,
unsigned
int
numberOfFreeParameters,
const
HistoChiSquare<T>
&
f
) {
65
unsigned
int
ndof
=
f
.numberOfBins() - numberOfFreeParameters;
66
std::cout
<<
"chi-squared/n.d.o.f. = "
<< amin <<
"/"
<<
ndof
<<
" = "
<< amin /
ndof
67
<<
"; prob: "
<< TMath::Prob(amin,
ndof
) << std::endl;
68
}
69
};
70
}
// namespace fit
71
72
#endif
fit::HistoChiSquare::setHistos
void setHistos(TH1 *histo)
Definition:
HistoChiSquare.h:35
fit::HistoChiSquare::deltaX_
double deltaX_
Definition:
HistoChiSquare.h:57
mps_fire.i
i
Definition:
mps_fire.py:355
f
double f[11][100]
Definition:
MuScleFitUtils.cc:78
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
timingPdfMaker.histo
histo
Definition:
timingPdfMaker.py:279
fit::HistoChiSquare::xMax_
double xMax_
Definition:
HistoChiSquare.h:57
fit::HistoChiSquare::HistoChiSquare
HistoChiSquare(T &t, TH1 *histo, double rangeMin, double rangeMax)
Definition:
HistoChiSquare.h:13
hltPixelTracks_cff.chi2
chi2
Definition:
hltPixelTracks_cff.py:25
fit::HistoChiSquare::nBins_
size_t nBins_
Definition:
HistoChiSquare.h:56
fit::HistoChiSquare::cont_
std::vector< double > cont_
Definition:
HistoChiSquare.h:58
fit::HistoChiSquare::rangeMax_
double rangeMax_
Definition:
HistoChiSquare.h:55
ndof
Definition:
HIMultiTrackSelector.h:49
fit::RootMinuitResultPrinter< HistoChiSquare< T > >::print
static void print(double amin, unsigned int numberOfFreeParameters, const HistoChiSquare< T > &f)
Definition:
HistoChiSquare.h:64
OrderedSet.t
t
Definition:
OrderedSet.py:90
fit::HistoChiSquare::rangeMin_
double rangeMin_
Definition:
HistoChiSquare.h:55
fit::HistoChiSquare
Definition:
HistoChiSquare.h:10
fit::HistoChiSquare::t_
T * t_
Definition:
HistoChiSquare.h:54
fit::HistoChiSquare::numberOfBins
size_t numberOfBins() const
Definition:
HistoChiSquare.h:41
fit::HistoChiSquare::xMin_
double xMin_
Definition:
HistoChiSquare.h:57
fit::HistoChiSquare::HistoChiSquare
HistoChiSquare()
Definition:
HistoChiSquare.h:12
RootMinuitResultPrinter.h
alignCSCRings.r
r
Definition:
alignCSCRings.py:93
T
long double T
Definition:
Basic3DVectorLD.h:48
fit::RootMinuitResultPrinter
Definition:
RootMinuitResultPrinter.h:8
fit
Definition:
CombinedChiSquaredLikelihood.h:6
fit::HistoChiSquare::operator()
double operator()() const
Definition:
HistoChiSquare.h:24
fit::HistoChiSquare::err_
std::vector< double > err_
Definition:
HistoChiSquare.h:59
Generated for CMSSW Reference Manual by
1.8.16