1 | /* ======================================================================== *\
|
---|
2 | !
|
---|
3 | ! *
|
---|
4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
8 | ! *
|
---|
9 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
10 | ! * documentation for any purpose is hereby granted without fee,
|
---|
11 | ! * provided that the above copyright notice appear in all copies and
|
---|
12 | ! * that both that copyright notice and this permission notice appear
|
---|
13 | ! * in supporting documentation. It is provided "as is" without express
|
---|
14 | ! * or implied warranty.
|
---|
15 | ! *
|
---|
16 | !
|
---|
17 | !
|
---|
18 | ! Author(s): Thomas Bretz 1/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | ! Author(s): Wolfgang Wittek 1/2002 <mailto:wittek@mppmu.mpg.de>
|
---|
20 | !
|
---|
21 | ! Copyright: MAGIC Software Development, 2000-2002
|
---|
22 | !
|
---|
23 | !
|
---|
24 | \* ======================================================================== */
|
---|
25 |
|
---|
26 | //////////////////////////////////////////////////////////////////////////////
|
---|
27 | // //
|
---|
28 | // MHEffOnTimeTime //
|
---|
29 | // //
|
---|
30 | // calculates the effective on time for each bin in time //
|
---|
31 | // //
|
---|
32 | //////////////////////////////////////////////////////////////////////////////
|
---|
33 |
|
---|
34 | #include "MHEffOnTimeTime.h"
|
---|
35 |
|
---|
36 | #include <TStyle.h>
|
---|
37 |
|
---|
38 | #include <TF1.h>
|
---|
39 | #include <TH2.h>
|
---|
40 | #include <TCanvas.h>
|
---|
41 |
|
---|
42 | #include "MTime.h"
|
---|
43 |
|
---|
44 | #include "MBinning.h"
|
---|
45 | #include "MParList.h"
|
---|
46 |
|
---|
47 | #include "MLog.h"
|
---|
48 | #include "MLogManip.h"
|
---|
49 |
|
---|
50 | ClassImp(MHEffOnTimeTime);
|
---|
51 |
|
---|
52 |
|
---|
53 | // --------------------------------------------------------------------------
|
---|
54 | //
|
---|
55 | // Default Constructor. It sets name and title of the histograms.
|
---|
56 | //
|
---|
57 | MHEffOnTimeTime::MHEffOnTimeTime(const char *name, const char *title)
|
---|
58 | : fHEffOn()
|
---|
59 | {
|
---|
60 | //
|
---|
61 | // set the name and title of this object
|
---|
62 | //
|
---|
63 | fName = name ? name : "MHEffOnTimeTime";
|
---|
64 | fTitle = title ? title : "1-D histogram of Eff On Time";
|
---|
65 |
|
---|
66 | // effective on time versus time
|
---|
67 | fHEffOn.SetName("EffOn");
|
---|
68 | fHEffOn.SetTitle("Effective On Time vs. Time");
|
---|
69 |
|
---|
70 | fHEffOn.SetDirectory(NULL);
|
---|
71 |
|
---|
72 | fHEffOn.SetXTitle("time [s]");
|
---|
73 | fHEffOn.SetYTitle("t_{eff} [s]");
|
---|
74 |
|
---|
75 | // chi2/NDF versus time
|
---|
76 | fHChi2.SetName("Chi2/NDF");
|
---|
77 | fHChi2.SetTitle("Chi2/NDF of OnTimeFit vs. Time");
|
---|
78 |
|
---|
79 | fHChi2.SetDirectory(NULL);
|
---|
80 |
|
---|
81 | fHChi2.SetXTitle("time [s]");
|
---|
82 | fHChi2.SetYTitle("chi^{2}/NDF");
|
---|
83 |
|
---|
84 | // lambda versus time
|
---|
85 | fHLambda.SetName("lambda");
|
---|
86 | fHLambda.SetTitle("lambda of OnTimeFit vs. Time");
|
---|
87 |
|
---|
88 | fHLambda.SetDirectory(NULL);
|
---|
89 |
|
---|
90 | fHLambda.SetXTitle("time [s]");
|
---|
91 | fHLambda.SetYTitle("\\lambda [Hz]");
|
---|
92 |
|
---|
93 | // N0del versus time
|
---|
94 | fHN0del.SetName("N0del");
|
---|
95 | fHN0del.SetTitle("N0del of OnTimeFit vs. Time");
|
---|
96 |
|
---|
97 | fHN0del.SetDirectory(NULL);
|
---|
98 |
|
---|
99 | fHN0del.SetXTitle("time [s]");
|
---|
100 | fHN0del.SetYTitle("N0del");
|
---|
101 | }
|
---|
102 |
|
---|
103 | // -----------------------------------------------------------------------
|
---|
104 | //
|
---|
105 | // Calculate the effective on time by fitting the distribution of
|
---|
106 | // time differences
|
---|
107 | //
|
---|
108 | void MHEffOnTimeTime::Calc(TH2D *hist)
|
---|
109 | {
|
---|
110 | // nbins = number of time bins
|
---|
111 | const Int_t nbins = hist->GetNbinsY();
|
---|
112 |
|
---|
113 | for (int i=1; i<=nbins; i++)
|
---|
114 | {
|
---|
115 |
|
---|
116 | // TH1D &h = *hist->ProjectionX("Calc-time", i, i, "E");
|
---|
117 | TH1D &h = *hist->ProjectionX("Calc-time", i, i, "E");
|
---|
118 |
|
---|
119 |
|
---|
120 | // Nmdel = Nm * binwidth, with Nm = number of observed events
|
---|
121 | const Double_t Nmdel = h.Integral("width");
|
---|
122 | const Double_t Nm = h.Integral();
|
---|
123 | // Double_t mean = h->GetMean();
|
---|
124 |
|
---|
125 | //...................................................
|
---|
126 | // determine range (yq[0], yq[1]) of time differences
|
---|
127 | // where fit should be performed;
|
---|
128 | // require a fraction >=xq[0] of all entries to lie below yq[0]
|
---|
129 | // and a fraction <=xq[1] of all entries to lie below yq[1];
|
---|
130 | // within the range (yq[0], yq[1]) there must be no empty bin;
|
---|
131 | // choose pedestrian approach as long as GetQuantiles is not available
|
---|
132 |
|
---|
133 | Double_t xq[2] = { 0.15, 0.95 };
|
---|
134 | Double_t yq[2];
|
---|
135 |
|
---|
136 | // GetQuantiles doesn't seem to be available in root 3.01/06
|
---|
137 | // h->GetQuantiles(2,yq,xq);
|
---|
138 |
|
---|
139 | const Double_t sumtot = h.Integral();
|
---|
140 | const Int_t jbins = h.GetNbinsX();
|
---|
141 |
|
---|
142 | if (sumtot > 0.0)
|
---|
143 | {
|
---|
144 | // char txt[100];
|
---|
145 | // sprintf(txt, "time_bin:%d", i);
|
---|
146 | // new TCanvas(txt, txt);
|
---|
147 |
|
---|
148 | Double_t sum1 = 0.0;
|
---|
149 | yq[0] = h.GetBinLowEdge(jbins+1);
|
---|
150 | for (int j=1; j<=jbins; j++)
|
---|
151 | {
|
---|
152 | if (sum1 >= xq[0]*sumtot)
|
---|
153 | {
|
---|
154 | yq[0] = h.GetBinLowEdge(j);
|
---|
155 | break;
|
---|
156 | }
|
---|
157 | sum1 += h.GetBinContent(j);
|
---|
158 | }
|
---|
159 |
|
---|
160 | Double_t sum2 = 0.0;
|
---|
161 | yq[1] = h.GetBinLowEdge(jbins+1);
|
---|
162 | for (int j=1; j<=jbins; j++)
|
---|
163 | {
|
---|
164 | const Double_t content = h.GetBinContent(j);
|
---|
165 | sum2 += content;
|
---|
166 | if (sum2 >= xq[1]*sumtot || content == 0.0)
|
---|
167 | {
|
---|
168 | yq[1] = h.GetBinLowEdge(j);
|
---|
169 | break;
|
---|
170 | }
|
---|
171 | }
|
---|
172 |
|
---|
173 | //...................................................
|
---|
174 |
|
---|
175 | // parameter 0 = lambda
|
---|
176 | // parameter 1 = N0*del with N0 = ideal number of events
|
---|
177 | // and del = bin width of time difference
|
---|
178 | TF1 func("Poisson", "[1] * [0] * exp(-[0] *x)", yq[0], yq[1]);
|
---|
179 |
|
---|
180 | func.SetParameter(0, 100); // [Hz]
|
---|
181 | func.SetParameter(1, Nmdel);
|
---|
182 |
|
---|
183 | func.SetParLimits(0, 0, 1000); // [Hz]
|
---|
184 | func.SetParLimits(1, 0, 10*Nmdel);
|
---|
185 |
|
---|
186 | func.SetParName(0, "lambda");
|
---|
187 | func.SetParName(1, "Nmdel");
|
---|
188 |
|
---|
189 | // options : 0 (=zero) do not plot the function
|
---|
190 | // I use integral of function in bin rather than value at bin center
|
---|
191 | // R use the range specified in the function range
|
---|
192 | // Q quiet mode
|
---|
193 | h.Fit("Poisson", "0IRQ");
|
---|
194 |
|
---|
195 | // gPad->SetLogy();
|
---|
196 | // gStyle->SetOptFit(1011);
|
---|
197 | // h->GetXaxis()->SetTitle("time difference [s]");
|
---|
198 | // h->GetYaxis()->SetTitle("Counts");
|
---|
199 | // h->DrawCopy();
|
---|
200 |
|
---|
201 | // func.SetRange(yq[0], yq[1]); // Range of Drawing
|
---|
202 | // func.DrawCopy("same");
|
---|
203 |
|
---|
204 | const Double_t lambda = func.GetParameter(0);
|
---|
205 | const Double_t N0del = func.GetParameter(1);
|
---|
206 | const Double_t chi2 = func.GetChisquare();
|
---|
207 | const Int_t NDF = func.GetNDF();
|
---|
208 |
|
---|
209 | // was fit successful ?
|
---|
210 | if (NDF>0 && chi2<2.5*NDF)
|
---|
211 | {
|
---|
212 | // the effective on time is Nm/lambda
|
---|
213 | fHEffOn.SetBinContent(i, Nm/lambda);
|
---|
214 |
|
---|
215 | // plot chi2/NDF of fit
|
---|
216 | fHChi2.SetBinContent(i, NDF ? chi2/NDF : 0.0);
|
---|
217 |
|
---|
218 | // lambda of fit
|
---|
219 | fHLambda.SetBinContent(i, lambda);
|
---|
220 |
|
---|
221 | // N0del of fit
|
---|
222 | fHN0del.SetBinContent(i, N0del);
|
---|
223 |
|
---|
224 | delete &h;
|
---|
225 | continue;
|
---|
226 | }
|
---|
227 | }
|
---|
228 |
|
---|
229 | fHEffOn.SetBinContent (i, 1.e-20);
|
---|
230 | fHChi2.SetBinContent (i, 1.e-20);
|
---|
231 | fHLambda.SetBinContent(i, 1.e-20);
|
---|
232 | fHN0del.SetBinContent (i, 1.e-20);
|
---|
233 |
|
---|
234 | delete &h;
|
---|
235 | }
|
---|
236 | }
|
---|
237 |
|
---|
238 | // -------------------------------------------------------------------------
|
---|
239 | //
|
---|
240 | // Set the binnings and prepare the filling of the histograms
|
---|
241 | //
|
---|
242 | Bool_t MHEffOnTimeTime::SetupFill(const MParList *plist)
|
---|
243 | {
|
---|
244 | const MBinning* binstime = (MBinning*)plist->FindObject("BinningTime");
|
---|
245 | if (!binstime)
|
---|
246 | {
|
---|
247 | *fLog << err << dbginf << "BinningTime [MBinning] not found... aborting." << endl;
|
---|
248 | return kFALSE;
|
---|
249 | }
|
---|
250 |
|
---|
251 | SetBinning(&fHEffOn, binstime);
|
---|
252 | SetBinning(&fHChi2, binstime);
|
---|
253 | SetBinning(&fHLambda, binstime);
|
---|
254 | SetBinning(&fHN0del, binstime);
|
---|
255 |
|
---|
256 | fHEffOn.Sumw2();
|
---|
257 | fHChi2.Sumw2();
|
---|
258 | fHLambda.Sumw2();
|
---|
259 | fHN0del.Sumw2();
|
---|
260 |
|
---|
261 | return kTRUE;
|
---|
262 | }
|
---|
263 |
|
---|
264 | // -------------------------------------------------------------------------
|
---|
265 | //
|
---|
266 | // Dummy Fill
|
---|
267 | // without it get error message :
|
---|
268 | // Error: MHEffOnTimeTime() no default constructor FILE:macros/wowflux.C LINE:359
|
---|
269 | //*** Interpreter error recovered ***
|
---|
270 | Bool_t MHEffOnTimeTime::Fill(const MParContainer *par)
|
---|
271 | {
|
---|
272 | return kTRUE;
|
---|
273 | }
|
---|
274 |
|
---|
275 | // -------------------------------------------------------------------------
|
---|
276 | //
|
---|
277 | // Draw a copy of the histogram
|
---|
278 | //
|
---|
279 | TObject *MHEffOnTimeTime::DrawClone(Option_t *opt) const
|
---|
280 | {
|
---|
281 | TCanvas &c = *MakeDefCanvas("EffOnTimeTime", "Results from on time fit vs. time");
|
---|
282 | c.Divide(2, 2);
|
---|
283 |
|
---|
284 | gROOT->SetSelectedPad(NULL);
|
---|
285 |
|
---|
286 | c.cd(1);
|
---|
287 | ((TH2*)&fHEffOn)->DrawCopy(opt);
|
---|
288 |
|
---|
289 | c.cd(2);
|
---|
290 | ((TH2*)&fHChi2)->DrawCopy(opt);
|
---|
291 |
|
---|
292 | c.cd(3);
|
---|
293 | ((TH2*)&fHLambda)->DrawCopy(opt);
|
---|
294 |
|
---|
295 | c.cd(4);
|
---|
296 | ((TH2*)&fHN0del)->DrawCopy(opt);
|
---|
297 |
|
---|
298 | c.Modified();
|
---|
299 | c.Update();
|
---|
300 |
|
---|
301 | return &c;
|
---|
302 | }
|
---|
303 |
|
---|
304 | // -------------------------------------------------------------------------
|
---|
305 | //
|
---|
306 | // Draw the histogram
|
---|
307 | //
|
---|
308 | void MHEffOnTimeTime::Draw(Option_t *opt)
|
---|
309 | {
|
---|
310 | if (!gPad)
|
---|
311 | MakeDefCanvas("EffOnTimeTime", "Results from on time fit vs. time");
|
---|
312 |
|
---|
313 | gPad->Divide(2,2);
|
---|
314 |
|
---|
315 | gPad->cd(1);
|
---|
316 | fHEffOn.Draw(opt);
|
---|
317 |
|
---|
318 | gPad->cd(2);
|
---|
319 | fHChi2.Draw(opt);
|
---|
320 |
|
---|
321 | gPad->cd(3);
|
---|
322 | fHLambda.Draw(opt);
|
---|
323 |
|
---|
324 | gPad->cd(4);
|
---|
325 | fHN0del.Draw(opt);
|
---|
326 |
|
---|
327 | gPad->Modified();
|
---|
328 | gPad->Update();
|
---|
329 | }
|
---|
330 |
|
---|
331 |
|
---|
332 |
|
---|
333 |
|
---|
334 |
|
---|