source: trunk/MagicSoft/Mars/mhist/MHEffOnTime.cc@ 2777

Last change on this file since 2777 was 2173, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 11.0 KB
Line 
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): Wolfgang Wittek 5/2002 <mailto:wittek@mppmu.mpg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MHEffOnTime //
28// //
29// calculates the effective on time for each bin in the variable Var //
30// (Var may be time or Theta) //
31// //
32// Important : The sample of events used for this should be as close //
33// to the sample of recorded events as possible. //
34// This means that NO additional cuts (be it quality cuts or //
35// gamma/hadron cuts) should be applied (see MAGIC-TDAS 02-02).//
36// //
37//////////////////////////////////////////////////////////////////////////////
38
39#include "MHEffOnTime.h"
40
41#include <TStyle.h>
42
43#include <TMinuit.h>
44#include <TFitter.h>
45
46#include <TF1.h>
47#include <TH2.h>
48#include <TCanvas.h>
49
50#include "MBinning.h"
51#include "MParList.h"
52
53#include "MLog.h"
54#include "MLogManip.h"
55
56#include "MHTimeDiffTheta.h"
57
58ClassImp(MHEffOnTime);
59
60using namespace std;
61
62// --------------------------------------------------------------------------
63//
64// Default Constructor. It sets the name of the variable ("Time" or "Theta")
65// and the units of the variable ("[s]" or "[\\circ])")
66//
67MHEffOnTime::MHEffOnTime(const char *varname, const char *unit)
68 : fHEffOn(), fHProb(), fHLambda(), fHRdead()
69{
70 fVarname = varname;
71 fUnit = unit;
72
73 TString strg = fVarname + " " + fUnit;
74
75 //
76 // set the name and title of this object
77 //
78 fName = TString("MHEffOnTime-")+fVarname;
79 fTitle = "1-D histogram of Eff On Time";
80
81 // effective on time versus Var
82 fHEffOn.SetName("EffOn");
83 fHEffOn.SetTitle(TString("T_{on, eff} vs. ")+fVarname);
84
85 fHEffOn.SetDirectory(NULL);
86
87 fHEffOn.SetXTitle(strg);
88 fHEffOn.SetYTitle("T_{on, eff} [s]");
89
90 // chi2-probability versus Var
91 fHProb.SetName("Chi2-prob");
92 TString strg3("\\chi^{2}-prob of OnTimeFit vs. ");
93 strg3 += fVarname;
94 fHProb.SetTitle(strg3);
95
96 fHProb.SetDirectory(NULL);
97
98 fHProb.SetXTitle(strg);
99 fHProb.SetYTitle("\\chi^{2}-probability");
100
101 // lambda versus Var
102 fHLambda.SetName("lambda");
103 fHLambda.SetTitle(TString("\\lambda from OnTimeFit vs. ")+fVarname);
104
105 fHLambda.SetDirectory(NULL);
106
107 fHLambda.SetXTitle(strg);
108 fHLambda.SetYTitle("\\lambda [Hz]");
109
110 // Rdead versus Var
111 // Rdead is the fraction from real time lost by the dead time
112 fHRdead.SetName("Rdead");
113
114 fHRdead.SetTitle(TString("Rdead of OnTimeFit vs. ")+fVarname);
115
116 fHRdead.SetDirectory(NULL);
117
118 fHRdead.SetXTitle(strg);
119 fHRdead.SetYTitle("Rdead");
120
121}
122
123// --------------------------------------------------------------------
124//
125// determine range (yq[0], yq[1]) of time differences
126// where fit should be performed;
127// require a fraction >=min of all entries to lie below yq[0]
128// and a fraction <=max of all entries to lie below yq[1];
129// within the range (yq[0], yq[1]) there must be no empty bin;
130//
131void MHEffOnTime::GetQuantiles(const TH1 &h, Double_t min, Double_t max, Double_t yq[2]) const
132{
133#if ROOT_VERSION_CODE < ROOT_VERSION(3,02,07)
134 // WOrkaround for missing const qualifier of TH1::Integral
135 TH1 &h1 = (TH1&)h;
136
137 // choose pedestrian approach as long as GetQuantiles is
138 // not available
139 const Int_t jbins = h1.GetNbinsX();
140 const Double_t Nm = h1.Integral();
141
142 const Double_t xq[2] = { 0.15*Nm, 0.98*Nm };
143
144 yq[0] = yq[1] = h1.GetBinLowEdge(jbins+1);
145
146 for (int j=1; j<=jbins; j++)
147 if (h1.Integral(2, j) >= xq[0])
148 {
149 yq[0] = h1.GetBinLowEdge(j);
150 break;
151 }
152
153 for (int j=1; j<=jbins; j++)
154 if (h1.Integral(1, j) >= xq[1] || h1.GetBinContent(j)==0)
155 {
156 yq[1] = h1.GetBinLowEdge(j);
157 break;
158 }
159#else
160 // GetQuantiles doesn't seem to be available in root 3.01/06
161 Double_t xq[2] = { min, max };
162 ((TH1&)h).GetQuantiles(2, yq, xq);
163#endif
164}
165
166void MHEffOnTime::DrawBin(TH1 &h, Int_t i) const
167{
168 TString strg1 = fVarname+"-bin #";
169 strg1 += i;
170
171 new TCanvas(strg1, strg1);
172
173 gPad->SetLogy();
174 gStyle->SetOptFit(1011);
175
176 TString name="Bin_";
177 name += i;
178
179 h.SetName(name);
180 h.SetXTitle("\\Delta t [s]");
181 h.SetYTitle("Counts");
182 h.DrawCopy();
183}
184
185Bool_t MHEffOnTime::CalcResults(const TF1 &func, Double_t Nm, Int_t i)
186{
187 const Double_t lambda = func.GetParameter(0);
188 const Double_t N0 = func.GetParameter(1);
189 const Double_t prob = func.GetProb();
190 const Int_t NDF = func.GetNDF();
191
192 Double_t xmin, xmax;
193 ((TF1&)func).GetRange(xmin, xmax);
194
195 *fLog << inf;
196 *fLog << "Fitted bin #" << i << " from " << xmin << " to " << xmax;
197 *fLog << ", got: lambda=" << lambda << "Hz N0=" << N0 << endl;
198
199 if (prob<=0.01)
200 {
201 *fLog << warn << "WARN - Fit bin#" << i << " gives:";
202 *fLog << " Chi^2-Probab(" << prob << ")<0.01";
203 *fLog << " NoFitPts=" << func.GetNumberFitPoints();
204 *fLog << " Chi^2=" << func.GetChisquare();
205 }
206
207 // was fit successful ?
208 if (NDF<=0 || /*prob<=0.001 ||*/ lambda<=0 || N0<=0)
209 {
210 *fLog << warn << dbginf << "Fit failed bin #" << i << ": ";
211 if (NDF<=0)
212 *fLog << " NDF(Number of Degrees of Freedom)=0";
213 if (lambda<=0)
214 *fLog << " Parameter#0(lambda)=0";
215 if (N0<=0)
216 *fLog << " Parameter#1(N0)=0";
217 *fLog << endl;
218
219 return kFALSE;
220 }
221
222 //
223 // -------------- start error calculation ----------------
224 //
225 Double_t emat[2][2];
226 gMinuit->mnemat(&emat[0][0], 2);
227
228 //
229 // Rdead : fraction of real time lost by the dead time
230 // kappa = number of observed events (Nm) divided by
231 // the number of genuine events (N0)
232 // Teff : effective on-time
233 //
234 const Double_t Teff = Nm/lambda;
235 const Double_t kappa = Nm/N0;
236 const Double_t Rdead = 1.0 - kappa;
237
238 const Double_t dldl = emat[0][0];
239 const Double_t dN0dN0 = emat[1][1];
240
241 const Double_t dTeff = Teff * sqrt(dldl/(lambda*lambda) + 1.0/Nm);
242 const Double_t dl = sqrt(dldl);
243 const Double_t dRdead = kappa * sqrt(dN0dN0/(N0*N0) + 1.0/Nm);
244 //
245 // -------------- end error calculation ----------------
246 //
247
248 // the effective on time is Nm/lambda
249 fHEffOn.SetBinContent(i, Teff);
250 fHEffOn.SetBinError (i, dTeff);
251
252 // plot chi2-probability of fit
253 fHProb.SetBinContent(i, prob);
254
255 // lambda from fit
256 fHLambda.SetBinContent(i, lambda);
257 fHLambda.SetBinError (i, dl);
258
259 // Rdead from fit
260 fHRdead.SetBinContent(i, Rdead);
261 fHRdead.SetBinError (i,dRdead);
262
263 return kTRUE;
264}
265
266void MHEffOnTime::ResetBin(Int_t i)
267{
268 fHEffOn.SetBinContent (i, 1.e-20);
269 fHProb.SetBinContent (i, 1.e-20);
270 fHLambda.SetBinContent(i, 1.e-20);
271 fHRdead.SetBinContent (i, 1.e-20);
272}
273
274// -----------------------------------------------------------------------
275//
276// Calculate the effective on time by fitting the distribution of
277// time differences
278//
279void MHEffOnTime::Calc(const TH2D *hist, const Bool_t draw)
280{
281 // nbins = number of Var bins
282 const Int_t nbins = hist->GetNbinsY();
283
284 for (int i=1; i<=nbins; i++)
285 {
286 TH1D &h = *((TH2D*)hist)->ProjectionX(TString("Calc-")+fVarname,
287 i, i, "E");
288 if (draw)
289 DrawBin(h, i);
290
291 ResetBin(i);
292
293 // Nmdel = Nm * binwidth, with Nm = number of observed events
294 const Double_t Nm = h.Integral();
295 const Double_t Nmdel = h.Integral("width");
296 if (Nm <= 0)
297 {
298 *fLog << warn << dbginf << "Nm<0 for bin #" << i << endl;
299 delete &h;
300 continue;
301 }
302
303 Double_t yq[2];
304 GetQuantiles(h, 0.15, 0.95, yq);
305
306 //
307 // Setup Poisson function for the fit:
308 // lambda [Hz], N0 = ideal no of evts, del = bin width of dt
309 //
310 TF1 func("Poisson", " [1]*[2] * [0] * exp(-[0] *x)", yq[0], yq[1]);
311 func.SetParNames("lambda", "N0", "del");
312
313 func.SetParameter(0, 1);
314 func.SetParameter(1, Nm);
315 func.FixParameter(2, Nmdel/Nm);
316
317 // For me (TB) it seems that setting parameter limits isn't necessary
318
319 // For a description of the options see TH1::Fit
320 h.Fit("Poisson", "0IRQ");
321
322 // Calc and fill results of fit into the histograms
323 const Bool_t rc = CalcResults(func, Nm, i);
324
325 // draw the distribution of time differences if requested
326 if (rc && draw)
327 func.DrawCopy("same");
328
329 delete &h;
330 }
331}
332
333// -------------------------------------------------------------------------
334//
335// Set the binnings and prepare the filling of the histograms
336//
337Bool_t MHEffOnTime::SetupFill(const MParList *plist)
338{
339 TString bn = "Binning"+fVarname;
340
341 const MBinning* binsVar = (MBinning*)plist->FindObject(bn);
342
343 if (!binsVar)
344 {
345 *fLog << err << dbginf << bn << " [MBinning] not found... aborting." << endl;
346 return kFALSE;
347 }
348
349 SetBinning(&fHEffOn, binsVar);
350 SetBinning(&fHProb, binsVar);
351 SetBinning(&fHLambda, binsVar);
352 SetBinning(&fHRdead, binsVar);
353
354 fHEffOn.Sumw2();
355 fHProb.Sumw2();
356 fHLambda.Sumw2();
357 fHRdead.Sumw2();
358
359 return kTRUE;
360}
361
362// -------------------------------------------------------------------------
363//
364// Draw the histogram
365//
366void MHEffOnTime::Draw(Option_t *opt)
367{
368 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
369 pad->SetBorderMode(0);
370
371 pad->Divide(2,2);
372
373 pad->cd(1);
374 gPad->SetBorderMode(0);
375 fHEffOn.Draw(opt);
376
377 pad->cd(2);
378 gPad->SetBorderMode(0);
379 fHProb.Draw(opt);
380
381 pad->cd(3);
382 gPad->SetBorderMode(0);
383 fHLambda.Draw(opt);
384
385 pad->cd(4);
386 gPad->SetBorderMode(0);
387 fHRdead.Draw(opt);
388
389 pad->Modified();
390 pad->Update();
391}
392
393void MHEffOnTime::Calc(const MHTimeDiffTheta &hist, const Bool_t Draw)
394{
395 Calc(hist.GetHist(), Draw);
396}
397
Note: See TracBrowser for help on using the repository browser.