source: trunk/MagicSoft/Mars/mtemp/mucm/classes/MHMyFindSignificanceONOFF.h@ 6976

Last change on this file since 6976 was 6634, checked in by marcos, 20 years ago
*** empty log message ***
File size: 3.7 KB
Line 
1#ifndef MARS_MHMyFindSignificanceONOFF
2#define MARS_MHMyFindSignificanceONOFF
3
4#ifndef MARS_MH
5#include "MH.h"
6#endif
7#ifndef MARS_MHFindSignificance
8#include "MHFindSignificance.h"
9#endif
10
11#ifndef ROOT_TArrayD
12#include <TArrayD.h>
13#endif
14
15class TF1;
16class TH1;
17class TCanvas;
18
19class MHMyFindSignificanceONOFF : public MH
20{
21 private:
22
23 TH1 *fHistOrigON; // original plot of |alpha| (0.0 to 90.0 degrees)
24 TH1 *fHistON; // copy of fHistOrig or rebinned histogram
25 TH1 *fHistOrigOFF; // original plot of |alpha| (0.0 to 90.0 degrees)
26 TH1 *fHistOFF; // copy of fHistOrig or rebinned histogram
27 TH1 *fHistOFFNormalized; // fHistOFF normalized (contents and errors) with
28
29 MHFindSignificance fFindsigON; // Will contain the fit of the on data
30 MHFindSignificance fFindsigOFF;
31
32
33 Double_t fAlphasig;
34 Double_t fAlphaminON;
35 Double_t fAlphamaxON;
36 Double_t fAlphaminOFF;
37 Double_t fAlphamaxOFF;
38 Int_t fDegree;
39
40 Double_t fNormFactor; // Normalization factor between ON and OFF hists.
41 Double_t fNormFactorError;
42
43 // Couted quatities
44 Double_t fNon;
45 Double_t fdNon;
46 Double_t fNbg; // Counts in the signal region of OFF hist, NO normalized !!!
47 Double_t fdNbg;
48 Double_t fNex;
49 Double_t fdNex;
50 Double_t fSigLiMa; // significance of gamma signal according to Li & Ma
51
52 // from the polynomial fit of the OFF data
53 Double_t fNbgFit;
54 Double_t fdNbgFit;
55 Double_t fNexFit;
56 Double_t fdNexFit;
57 Double_t fGamma;
58 Double_t fNoff;
59 Double_t fSigLiMaFit; // significance of gamma signal according to Li & Ma
60
61 Bool_t fRebin; // if true : allow rebinning of the alpha plot
62 Bool_t fReduceDegree; // if true : allow reducing of the order of the polynomial
63
64 const static Double_t fEps = 1.e-4; // tolerance for floating point
65 // comparisons
66
67
68 Double_t fBestSigma;
69 Double_t fBestAlphaCut;
70 Double_t fBestExcess;
71
72 void NormalizedHistOFF();
73 void DrawResults(Option_t* option="all");
74
75public:
76
77 MHMyFindSignificanceONOFF(const char *name=NULL, const char *title=NULL);
78
79 Bool_t FindSigmaONOFF(TH1 *histON, TH1 *histOFF, Double_t alphasig, Double_t alphamin, Double_t alphamax, Int_t degree);
80
81 Int_t CountEventsInRange(TH1* hist, Double_t alphalo, Double_t alphaup, Double_t* numevents, Double_t* numeventserror);
82
83 Int_t CalcNormalizationFactor(TH1* hon, TH1* hoff, Double_t alphalo, Double_t alphaup);
84
85 Bool_t SigmaLiMa(Double_t non, Double_t noff, Double_t gamma,
86 Double_t *siglima);
87
88
89 Double_t GetNormFactor() const { return fNormFactor; }
90 Double_t GetNon() const { return fNon; }
91 Double_t GetNbg() const { return fNbg*fNormFactor; } // <------!!!
92 Double_t GetNex() const { return fNex; }
93 Double_t GetdNex() const { return fdNex; }
94
95 Double_t GetSigLiMa() const { return fSigLiMa; }
96 Double_t GetNbgFit() const { return fNbgFit*fNormFactor; } //<------!!!
97 Double_t GetNexFit() const { return fNexFit; }
98 Double_t GetSigLiMaFit() const { return fSigLiMaFit; }
99
100
101 void SetRebin(Bool_t b=kTRUE);
102 void SetReduceDegree(Bool_t b=kTRUE);
103
104
105 void Draw(Option_t* option="all");
106
107 Bool_t SigmaVsAlpha(TH1 *histON, TH1 *histOFF, Double_t alphamin, Double_t alphamax, Int_t degree, Bool_t draw=kTRUE);
108
109
110
111 Double_t GetBestSigma() const {return fBestSigma;}
112 Double_t GetBestAlphaCut() const {return fBestAlphaCut;}
113 Double_t GetBestExcess() const {return fBestExcess;}
114
115 ClassDef(MHMyFindSignificanceONOFF, 1) // Determine significance from ON and OFF data
116};
117
118#endif
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
Note: See TracBrowser for help on using the repository browser.