source: trunk/MagicSoft/Mars/mhflux/MHEffectiveOnTime.cc@ 8802

Last change on this file since 8802 was 8695, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 32.2 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): Thomas Bretz, 8/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-2004
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27//
28// MHEffectiveOnTime
29//
30// Filling this you will get the effective on-time versus theta and
31// observation time.
32//
33// From this histogram the effective on-time is determined by a fit.
34// The result of the fit (see Fit()) and the fit-parameters (like chi^2)
35// are stored in corresponding histograms
36//
37// To determin the efective on time a poisson fit is done. For more details
38// please have a look into the source code of FitH() it should be simple
39// to understand. In this function a Delta-T distribution is fitted, while
40// Delta-T is the time between two consecutive events.
41//
42// The fit is done for projections of a 2D histogram in Theta and Delta T.
43// So you get the effective on time versus theta.
44//
45// To get the effective on-time versus time a histogram is filled with
46// the Delta-T distribution of a number of events set by SetNumEvents().
47// The default is 12000 (roughly 1min at 200Hz)
48//
49// For each "time-bin" the histogram is fitted and the resulting effective
50// on-time is stored in the fHTimeEffOn histogram. Each entry in this
51// histogram is the effective observation time between the upper and
52// lower edges of the bins.
53//
54// In addition the calculated effective on time is stored in a
55// "MEffectiveOnTime [MParameterDerr]" and the corresponding time-stamp
56// (the upper edge of the bin) "MTimeEffectiveOnTime [MTime]"
57//
58// The class takes two binnings from the Parameter list; if these binnings
59// are not available the defaultbinning is used:
60// MBinning("BinningDeltaT"); // Units of seconds
61// MBinning("BinningTheta"); // Units of degrees
62//
63//
64// Usage:
65// ------
66// MFillH fill("MHEffectiveOnTime", "MTime");
67// tlist.AddToList(&fill);
68//
69//
70// Input Container:
71// MPointingPos
72//
73// Output Container:
74// MEffectiveOnTime [MParameterDerr]
75// MTimeEffectiveOnTime [MTime]
76//
77//
78// Class version 2:
79// ----------------
80// + UInt_t fFirstBin;
81// + UInt_t fNumEvents;
82// - Int_t fNumEvents;
83//
84//
85// ==========================================================================
86// Dear Colleagues,
87//
88// for the case that we are taking calibration events interleaved with
89// cosmics events the calculation of the effective observation time has to
90// be modified. I have summarized the proposed procedures in the note at the
91// end of this message. The formulas have been checked by a simulation.
92//
93// Comments are welcome.
94//
95// Regards, Wolfgang
96// --------------------------------------------------------------------------
97// Wolfgang Wittek
98// 2 Dec. 2004
99//
100// Calculation of the effective observation time when cosmics and calibration
101// events are taken simultaneously.
102// --------------------------------
103//
104// I. Introduction
105// ---------------
106// It is planned to take light calibration events (at a certain fixed frequency
107// lambda_calib) interlaced with cosmics events. The advantages of this
108// procedure are :
109//
110// - the pedestals, which would be determined from the cosmics, could be
111// used for both the calibration and the cosmics events
112//
113// - because calibration and cosmics events are taken quasi simultaneously,
114// rapid variations (in the order of a few minutes) of base lines and of the
115// photon/ADC conversion factors could be recognized and taken into account
116//
117// The effective observation time T_eff is defined as that time range, within
118// which the recorded number of events N_cosmics would be obtained under ideal
119// conditions (only cosmics, no dead time, no calibration events, ...).
120//
121// In the absence of calibration events the effective observation time can
122// be determined from the distribution of time differences 'dt' between
123// successive cosmics events (see first figure in the attached ps file).
124// The exponential slope 'lambda' of this distribution is the ideal cosmics
125// event rate. If 'N_cosmics' is the total number of recorded cosmics events,
126// T_eff is obtained by
127//
128// T_eff = N_cosmics / lambda
129//
130// In the case of a finite dead time 'dead', the distribution (for dt > dead) is
131// still exponential with the same slope 'lambda'. 'lambda' should be determined
132// in a region of 'dt' which is not affected by the dead time, i.e. at not too
133// low 'dt'.
134//
135//
136//
137// II. Problems in the presence of calibration events
138// --------------------------------------------------
139// If calibration events are taken interlaced with cosmics, and if the dead time
140// is negligible, the distribution of time differences 'dt' between cosmics can
141// be used for calculating the effective observation time, as if the calibration
142// events were not present.
143//
144// In the case of a non-negligible dead time 'dead', however, the distribution of
145// time differences between cosmics is distorted, because a cosmics event may be
146// lost due to the dead time after a calibration event. Even if the time
147// intervals are ignored which contain a calibration event,
148//
149//
150// ---|---------o--------|---------> t
151//
152// cosmics calib cosmics
153//
154// <----------------> <==== time interval to be ignored
155//
156//
157// the distribution of 'dt' is still distorted, because there would be no
158// 'dt' with dt > tau_calib = 1/lambda_calib. The distribution would also be
159// distorted in the region dt < tau_calib, due to calibration events occuring
160// shortly after cosmics events. As a result, the slope of the distribution of
161// 'dt' would not reflect the ideal cosmics event rate (see second figure; the
162// values assumed in the simulation are lambda = 200 Hz, lambda_calib = 50
163// Hz, dead = 0.001 sec, total time = 500 sec, number of generated cosmics
164// events = 100 000).
165//
166//
167// Note also that some calibration events will not be recorded due to the dead
168// time after a cosmics event.
169//
170//
171// III. Proposed procedures
172// ------------------------
173//
174// A) The ideal event rate 'lambda' may be calculated from the distribution of
175// the time difference 'dt_first' between a calibration event and the first
176// recorded cosmics event after the calibration event. In the region
177//
178// dead < dt_first < tau_calib
179//
180// the probability distribution of dt_first is given by
181//
182// p(dt_first) = c * exp(-lambda*dt_first)
183//
184// where c is a normalization constant. 'lambda' can be obtained by a simple
185// exponential fit to the experimental distribution of dt_first (see third
186// figure). The fit range should start well above the average value of the dead
187// time 'dead'.
188//
189//
190// B) One may consider those time intervals between recorded cosmics events, which
191// are completely contained in the region
192//
193// t_calib < t < t_calib + tau_calib
194//
195// where t_calib is the time of a recorded calibration event.
196//
197//
198// <--------------- tau_calib ----------->
199//
200//
201// 0 1 2 3 4 5 6 7 8 9 10
202// --|-o---|-|---|--|-|----|--|---|---|-|----o-|---|-|---------> t
203// ^ ^
204// | |
205// t_calib t_calib + tau_calib
206//
207//
208// In this example, of the time intervals 0 to 10 only the intervals 1 to 9
209// should be retained and plotted. The distribution of the length 'dt' of these
210// intervals in the region
211//
212// dead < dt < tau_calib
213//
214// is given by
215//
216// p(dt) = c * (tau_calib-dt-dead) * exp(-lambda*dt)
217//
218// A fit of this expression to the experimental distribution of 'dt' yields
219// 'lambda' (see fourth figure). For 'dead' an average value of the dead time
220// should be chosen, and the fit range should end well before dt = tau_calib-dead.
221//
222//
223// Method A has the advantage that the p(dt_first) does not depend on 'dead'.
224// 'dead' has to be considered when defining the fit range, both in method A and
225// in method B. In method B the event statistics is larger leading to a smaller
226// fitted error of 'lambda' than method A (see the figures).
227//
228//
229// The effective observation time is again obtained by
230//
231// T_eff = N_cosmics / lambda
232//
233// where N_cosmics is the total number of recorded cosmics events. Note that
234// N_cosmics is equal to
235//
236// N_cosmics = N_tot - N_calib
237//
238// where N_tot is the total number of recorded events (including the calibration
239// events) and N_calib is the number of recorded calibration events.
240//
241// Note that if time intervals are discarded for the determination of lambda,
242// the corresponding cosmics events need not and should not be discarded.
243//
244//
245// IV. Procedure if the calibration events are taken in bunches
246// ------------------------------------------------------------
247// In November 2004 the rate of calibration events is not constant. The events
248// are taken in 200 Hz bunches every second, such that the rate is 200 Hz for
249// 0.25 sec, followed by a gap of 0.75 sec. Then follows the next 200 Hz bunch.
250//
251// In this case it is proposed to consider for the calculation of 'lambda' only
252// the cosmics events within the gaps of 0.75 sec. For these cosmics events one
253// of the methods described in III. can be applied.
254//
255//
256// V. Alternative pocedure
257// -----------------------
258// The effective observation time can also be determined from the total
259// observation time and the total dead time. The latter is written out by the DAQ.
260// In this case it has to be made sure that the dead time is available in Mars
261// when the effective observation time is calculated.
262//
263//////////////////////////////////////////////////////////////////////////////
264#include "MHEffectiveOnTime.h"
265
266#include <TF1.h>
267#include <TMinuit.h>
268#include <TRandom.h>
269
270#include <TLatex.h>
271#include <TGaxis.h>
272#include <TCanvas.h>
273#include <TPaveStats.h>
274
275#include "MTime.h"
276#include "MParameters.h"
277#include "MPointingPos.h"
278
279#include "MBinning.h"
280#include "MParList.h"
281
282#include "MLog.h"
283#include "MLogManip.h"
284
285ClassImp(MHEffectiveOnTime);
286
287using namespace std;
288
289// --------------------------------------------------------------------------
290//
291// Default Constructor. It initializes all histograms.
292//
293MHEffectiveOnTime::MHEffectiveOnTime(const char *name, const char *title)
294 : fPointPos(0), fTime(0), fParam(0), fIsFinalized(kFALSE),
295 fNumEvents(200*60), fFirstBin(3)
296 //fNumEvents(2*60), fFirstBin(1)
297{
298 //
299 // set the name and title of this object
300 //
301 fName = name ? name : "MHEffectiveOnTime";
302 fTitle = title ? title : "Histogram to determin effective On-Time vs Time and Zenith Angle";
303
304 // Main histogram
305 fH2DeltaT.SetName("DeltaT");
306 fH2DeltaT.SetXTitle("\\Delta t [s]");
307 fH2DeltaT.SetYTitle("\\Theta [\\circ]");
308 fH2DeltaT.SetZTitle("Count");
309 fH2DeltaT.UseCurrentStyle();
310 fH2DeltaT.SetDirectory(NULL);
311
312 // Main histogram
313 fH1DeltaT.SetName("DeltaT");
314 fH1DeltaT.SetXTitle("\\Delta t [s]");
315 fH1DeltaT.SetYTitle("Counts");
316 fH1DeltaT.UseCurrentStyle();
317 fH1DeltaT.SetDirectory(NULL);
318
319 // effective on time versus theta
320 fHThetaEffOn.SetName("EffOnTheta");
321 fHThetaEffOn.SetTitle("Effective On Time T_{eff}");
322 fHThetaEffOn.SetXTitle("\\Theta [\\circ]");
323 fHThetaEffOn.SetYTitle("T_{eff} [s]");
324 fHThetaEffOn.UseCurrentStyle();
325 fHThetaEffOn.SetDirectory(NULL);
326 fHThetaEffOn.GetYaxis()->SetTitleOffset(1.2);
327 fHThetaEffOn.GetYaxis()->SetTitleColor(kBlue);
328 fHThetaEffOn.SetLineColor(kBlue);
329 //fHEffOn.Sumw2();
330
331 // effective on time versus time
332 fHTimeEffOn.SetName("EffOnTime");
333 fHTimeEffOn.SetTitle("Effective On Time T_{eff}");
334 fHTimeEffOn.SetXTitle("Time");
335 fHTimeEffOn.SetYTitle("T_{eff} [s]");
336 fHTimeEffOn.UseCurrentStyle();
337 fHTimeEffOn.SetDirectory(NULL);
338 fHTimeEffOn.GetYaxis()->SetTitleOffset(1.2);
339 fHTimeEffOn.GetXaxis()->SetLabelSize(0.033);
340 fHTimeEffOn.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
341 fHTimeEffOn.GetXaxis()->SetTimeDisplay(1);
342 fHTimeEffOn.GetYaxis()->SetTitleColor(kBlue);
343 fHTimeEffOn.SetLineColor(kBlue);
344
345 // chi2 probability
346 fHThetaProb.SetName("ProbTheta");
347 fHThetaProb.SetTitle("\\chi^{2} Probability of Fit");
348 fHThetaProb.SetXTitle("\\Theta [\\circ]");
349 fHThetaProb.SetYTitle("p [%]");
350 fHThetaProb.UseCurrentStyle();
351 fHThetaProb.SetDirectory(NULL);
352 fHThetaProb.GetYaxis()->SetTitleOffset(1.2);
353 fHThetaProb.SetMaximum(101);
354 fHThetaProb.GetYaxis()->SetTitleColor(kBlue);
355 fHThetaProb.SetLineColor(kBlue);
356
357 // chi2 probability
358 fHTimeProb.SetName("ProbTime");
359 fHTimeProb.SetTitle("\\chi^{2} Probability of Fit");
360 fHTimeProb.SetXTitle("Time");
361 fHTimeProb.SetYTitle("p [%]");
362 fHTimeProb.UseCurrentStyle();
363 fHTimeProb.SetDirectory(NULL);
364 fHTimeProb.GetYaxis()->SetTitleOffset(1.2);
365 fHTimeProb.GetXaxis()->SetLabelSize(0.033);
366 fHTimeProb.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
367 fHTimeProb.GetXaxis()->SetTimeDisplay(1);
368 fHTimeProb.SetMaximum(101);
369 fHTimeProb.GetYaxis()->SetTitleColor(kBlue);
370 fHTimeProb.SetLineColor(kBlue);
371
372 // lambda versus theta
373 fHThetaLambda.SetName("LambdaTheta");
374 fHThetaLambda.SetTitle("Slope (Rate) vs Theta");
375 fHThetaLambda.SetXTitle("\\Theta [\\circ]");
376 fHThetaLambda.SetYTitle("\\lambda [s^{-1}]");
377 fHThetaLambda.UseCurrentStyle();
378 fHThetaLambda.SetDirectory(NULL);
379 fHThetaLambda.SetLineColor(kGreen);
380
381 // lambda versus time
382 fHTimeLambda.SetName("LambdaTime");
383 fHTimeLambda.SetTitle("Slope (Rate) vs Time");
384 fHTimeLambda.SetXTitle("\\Time [\\circ]");
385 fHTimeLambda.SetYTitle("\\lambda [s^{-1}]");
386 fHTimeLambda.UseCurrentStyle();
387 fHTimeLambda.SetDirectory(NULL);
388 fHTimeLambda.GetYaxis()->SetTitleOffset(1.2);
389 fHTimeLambda.GetXaxis()->SetLabelSize(0.033);
390 fHTimeLambda.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
391 fHTimeLambda.GetXaxis()->SetTimeDisplay(1);
392 fHTimeLambda.SetLineColor(kGreen);
393
394 // NDoF versus theta
395 fHThetaNDF.SetName("NDofTheta");
396 fHThetaNDF.SetTitle("Number of Degrees of freedom vs Theta");
397 fHThetaNDF.SetXTitle("\\Theta [\\circ]");
398 fHThetaNDF.SetYTitle("NDoF [#]");
399 fHThetaNDF.UseCurrentStyle();
400 fHThetaNDF.SetDirectory(NULL);
401 fHThetaNDF.SetLineColor(kGreen);
402
403 // NDoF versus time
404 /*
405 fHTimeNDF.SetName("NDofTime");
406 fHTimeNDF.SetTitle("Number of Degrees of freedom vs Time");
407 fHTimeNDF.SetXTitle("Time");
408 fHTimeNDF.SetYTitle("NDoF [#]");
409 fHTimeNDF.UseCurrentStyle();
410 fHTimeNDF.SetDirectory(NULL);
411 fHTimeNDF.GetYaxis()->SetTitleOffset(1.2);
412 fHTimeNDF.GetXaxis()->SetLabelSize(0.033);
413 fHTimeNDF.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
414 fHTimeNDF.GetXaxis()->SetTimeDisplay(1);
415 fHTimeNDF.SetLineColor(kBlue);
416 */
417 // setup binning
418 MBinning btheta("BinningTheta");
419 btheta.SetEdgesASin(67, -0.005, 0.665);
420
421 MBinning btime("BinningDeltaT");
422 btime.SetEdges(50, 0, 0.1);
423
424 MH::SetBinning(&fH2DeltaT, &btime, &btheta);
425
426 btime.Apply(fH1DeltaT);
427
428 btheta.Apply(fHThetaEffOn);
429 btheta.Apply(fHThetaLambda);
430 btheta.Apply(fHThetaNDF);
431 btheta.Apply(fHThetaProb);
432 //btheta.Apply(fHChi2);
433}
434
435// --------------------------------------------------------------------------
436//
437// Set the binnings and prepare the filling of the histogram
438//
439Bool_t MHEffectiveOnTime::SetupFill(const MParList *plist)
440{
441 fPointPos = (MPointingPos*)plist->FindObject("MPointingPos");
442 if (!fPointPos)
443 {
444 *fLog << err << dbginf << "MPointingPos not found... aborting." << endl;
445 return kFALSE;
446 }
447
448 // FIXME: Remove const-qualifier from base-class!
449 fTime = (MTime*)const_cast<MParList*>(plist)->FindCreateObj("MTime", "MTimeEffectiveOnTime");
450 if (!fTime)
451 return kFALSE;
452 fParam = (MParameterDerr*)const_cast<MParList*>(plist)->FindCreateObj("MParameterDerr", "MEffectiveOnTime");
453 if (!fParam)
454 return kFALSE;
455
456 const MBinning* binsdtime = (MBinning*)plist->FindObject("BinningDeltaT");
457 const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
458 if (binsdtime)
459 binsdtime->Apply(fH1DeltaT);
460 if (binstheta)
461 {
462 binstheta->Apply(fHThetaEffOn);
463 binstheta->Apply(fHThetaLambda);
464 binstheta->Apply(fHThetaNDF);
465 binstheta->Apply(fHThetaProb);
466 //binstheta->Apply(fHChi2);
467 }
468 if (binstheta && binsdtime)
469 SetBinning(&fH2DeltaT, binsdtime, binstheta);
470
471 return kTRUE;
472}
473
474// --------------------------------------------------------------------------
475//
476// Fit a single Delta-T distribution. See source code for more details
477//
478Bool_t MHEffectiveOnTime::FitH(TH1D *h, Double_t *res, Bool_t paint) const
479{
480 const Double_t Nm = h->Integral();
481
482 // FIXME: Do fit only if contents of bin has changed
483 if (Nm<=0 || h->GetBinContent(1)<=0)
484 return kFALSE;
485
486 // determine range (yq[0], yq[1]) of time differences
487 // where fit should be performed;
488 // require a fraction >=xq[0] of all entries to lie below yq[0]
489 // and a fraction <=xq[1] of all entries to lie below yq[1];
490 // within the range (yq[0], yq[1]) there must be no empty bin;
491 // choose pedestrian approach as long as GetQuantiles is not available
492 Double_t xq[2] = { 0.6, 0.95 }; // previously 0.99
493 Double_t yq[2];
494 h->GetQuantiles(2, yq, xq);
495
496 //
497 // Determine a good starting value for the slope
498 //
499 const TAxis &axe = *h->GetXaxis();
500 const UInt_t ibin = axe.FindFixBin(yq[1]);
501 const Double_t x1 = axe.GetBinCenter(ibin<=fFirstBin?fFirstBin+1:ibin);
502 const Double_t x0 = axe.GetBinCenter(fFirstBin);
503 const Double_t y1 = h->GetBinContent(ibin)>1 ? TMath::Log(h->GetBinContent(ibin)) : 0;
504 const Double_t y0 = TMath::Log(h->GetBinContent(fFirstBin));
505
506 // Estimated slope
507 const Float_t m = -(y1-y0)/(x1-x0);
508
509 //
510 // Setup exponential function for the fit:
511 //
512 // parameter 0 = rate [Hz]
513 // parameter 1 = normalization
514 //
515 TF1 func("Exp", " exp([1]-[0]*x)");
516
517 func.SetParameter(0, m); // Hz
518 func.SetParameter(1, log(h->GetBinContent(1))); // Hz
519
520 // We set a limit to make sure that almost empty histograms which
521 // are fitted dont't produce hang ups or crashes
522 func.SetParLimits(0, 0, 15000); // Hz
523
524 // options : N do not store the function, do not draw
525 // I use integral of function in bin rather than value at bin center
526 // R use the range specified in the function range
527 // Q quiet mode
528 // L Use log-likelihood (better for low statistics)
529 h->Fit(&func, "NIQEL", "", h->GetBinLowEdge(fFirstBin)/*yq[0]*/, yq[1]);
530
531 const Double_t chi2 = func.GetChisquare();
532 const Double_t prob = func.GetProb();
533 const Int_t NDF = func.GetNDF();
534
535 // was fit successful ?
536 const Bool_t ok = prob>0.001; //NDF>0 && chi2<3*NDF;
537
538 if (paint)
539 {
540 func.SetLineWidth(2);
541 func.SetLineColor(ok ? kGreen : kRed);
542 func.Paint("same");
543 }
544
545 // The effective on time is the "real rate" (slope of the exponential)
546 // divided by the total number of events (histogram integral including
547 // under- and overflows)
548 const Double_t lambda = func.GetParameter(0);
549 const Double_t dldl = func.GetParError(0)*func.GetParError(0);
550 const Double_t teff = lambda==0 ? 0 : Nm / lambda;
551 const Double_t dteff = lambda==0 ? 0 : teff * TMath::Sqrt(dldl/(lambda*lambda) + 1.0/Nm);
552 const Double_t dl = TMath::Sqrt(dldl);
553
554 // the effective on time is Nm/lambda
555 res[0] = teff;
556 res[1] = dteff;
557
558 // plot chi2-probability of fit
559 res[2] = prob*100;
560
561 // lambda of fit
562 res[3] = lambda;
563 res[4] = dl;
564
565 // NDoF of fit
566 res[5] = NDF;
567
568 // Chi2
569 res[6] = chi2;
570
571 return ok;
572}
573
574// --------------------------------------------------------------------------
575//
576// Fit a all bins of the distribution in theta. Store the result in the
577// Theta-Histograms
578//
579void MHEffectiveOnTime::FitThetaBins()
580{
581 fHThetaEffOn.Reset();
582 fHThetaProb.Reset();
583 fHThetaLambda.Reset();
584 fHThetaNDF.Reset();
585
586 // Use a random name to make sure the object is unique
587 const TString name = Form("CalcTheta%d", (UInt_t)gRandom->Uniform(999999999));
588
589 // nbins = number of Theta bins
590 const Int_t nbins = fH2DeltaT.GetNbinsY();
591
592 TH1D *h=0;
593 for (int i=1; i<=nbins; i++)
594 {
595 // TH1D &h = *hist->ProjectionX("Calc-theta", i, i);
596 h = fH2DeltaT.ProjectionX(name, i, i, "E");
597
598 Double_t res[7] = {0, 0, 0, 0, 0, 0, 0};
599 //if (!FitH(h, res))
600 // continue;
601 FitH(h, res);
602
603 if (res[0]==0)
604 continue;
605
606 // the effective on time is Nm/lambda
607 fHThetaEffOn.SetBinContent(i, res[0]);
608 fHThetaEffOn.SetBinError (i, res[1]);
609
610 // plot chi2-probability of fit
611 fHThetaProb.SetBinContent(i, res[2]);
612
613 // plot chi2/NDF of fit
614 //fHChi2.SetBinContent(i, res[3]);
615
616 // lambda of fit
617 fHThetaLambda.SetBinContent(i, res[3]);
618 fHThetaLambda.SetBinError (i, res[4]);
619
620 // NDoF of fit
621 fHThetaNDF.SetBinContent(i, res[5]);
622
623 // Rdead (from fit) is the fraction from real time lost by the dead time
624 //fHRdead.SetBinContent(i, Rdead);
625 //fHRdead.SetBinError (i,dRdead);
626 }
627
628 // Histogram is reused via gROOT->FindObject()
629 // Need to be deleted only once
630 if (h)
631 delete h;
632}
633
634// --------------------------------------------------------------------------
635//
636// Fit the single-time-bin histogram. Store the result in the
637// Time-Histograms
638//
639void MHEffectiveOnTime::FitTimeBin()
640{
641 //
642 // Fit histogram
643 //
644 Double_t res[7];
645 if (!FitH(&fH1DeltaT, res))
646 return;
647
648 // Reset Histogram
649 fH1DeltaT.Reset();
650
651 //
652 // Prepare Histogram
653 //
654
655 // Get number of bins
656 const Int_t n = fHTimeEffOn.GetNbinsX();
657
658 // Enhance binning
659 MBinning bins;
660 bins.SetEdges(fHTimeEffOn, 'x');
661 bins.AddEdge(fLastTime.GetAxisTime());
662 bins.Apply(fHTimeEffOn);
663 bins.Apply(fHTimeProb);
664 bins.Apply(fHTimeLambda);
665 //bins.Apply(fHTimeNDF);
666
667 //
668 // Fill histogram
669 //
670 fHTimeEffOn.SetBinContent(n, res[0]);
671 fHTimeEffOn.SetBinError(n, res[1]);
672
673 fHTimeProb.SetBinContent(n, res[2]);
674
675 fHTimeLambda.SetBinContent(n, res[3]);
676 fHTimeLambda.SetBinError(n, res[4]);
677
678 //fHTimeNDF.SetBinContent(n, res[5]);
679
680 //
681 // Now prepare output
682 //
683 fParam->SetVal(res[0], res[1]);
684 fParam->SetReadyToSave();
685
686 *fTime = fLastTime;
687
688 // Include the current event
689 fTime->Plus1ns();
690
691 *fLog << fLastTime << ": Val=" << res[0] << " Err=" << res[1] << endl;
692}
693
694// --------------------------------------------------------------------------
695//
696// Fill the histogram
697//
698Bool_t MHEffectiveOnTime::Fill(const MParContainer *par, const Stat_t w)
699{
700 const MTime *time = dynamic_cast<const MTime*>(par);
701 if (!time)
702 {
703 *fLog << err << "ERROR - MHEffectiveOnTime::Fill without argument or container doesn't inherit from MTime... abort." << endl;
704 return kFALSE;
705 }
706
707 //
708 // If this is the first call we have to initialize the time-histogram
709 //
710 if (fLastTime==MTime())
711 {
712 MBinning bins;
713 bins.SetEdges(1, time->GetAxisTime()-fNumEvents/200, time->GetAxisTime());
714 bins.Apply(fHTimeEffOn);
715 bins.Apply(fHTimeProb);
716 bins.Apply(fHTimeLambda);
717 //bins.Apply(fHTimeNDF);
718
719 fParam->SetVal(0, 0);
720 fParam->SetReadyToSave();
721
722 *fTime = *time;
723
724 // Make this 1ns before the first event!
725 fTime->Minus1ns();
726 }
727
728 //
729 // Fill time difference into the histograms
730 //
731 const Double_t dt = *time-fLastTime;
732 fLastTime = *time;
733
734 fH2DeltaT.Fill(dt, fPointPos->GetZd(), w);
735 fH1DeltaT.Fill(dt, w);
736
737 //
738 // If we reached the event number limit for the time-bins fit the
739 // histogram - if it fails try again when 1.6% more events available
740 //
741 const UInt_t n = (UInt_t)fH1DeltaT.GetEntries();
742 if (n>=fNumEvents && n%(fNumEvents/60)==0)
743 FitTimeBin();
744
745 return kTRUE;
746}
747
748// --------------------------------------------------------------------------
749//
750// Fit the theta projections of the 2D histogram and the 1D Delta-T
751// distribution
752//
753Bool_t MHEffectiveOnTime::Finalize()
754{
755 FitThetaBins();
756 FitTimeBin();
757
758 fIsFinalized = kTRUE;
759
760 return kTRUE;
761}
762
763// --------------------------------------------------------------------------
764//
765// Paint the integral and the error on top of the histogram
766//
767void MHEffectiveOnTime::PaintText(Double_t val, Double_t error, Double_t range) const
768{
769 TLatex text;
770 text.SetBit(TLatex::kTextNDC);
771 text.SetTextSize(0.04);
772
773 text.SetText(0.45, 0.94, Form("T_{eff} = %.1fs \\pm %.1fs", val, error));
774 text.Paint();
775
776 if (range<0)
777 return;
778
779 text.SetText(0.66, 0.94, Form("T_{axe} = %.1fs", range));
780 text.Paint();
781}
782
783void MHEffectiveOnTime::PaintText(Double_t *res) const
784{
785 TLatex text(0.27, 0.94, Form("T_{eff}=%.1fs\\pm%.1fs \\lambda=%.1f\\pm%.1fHz p=%.1f%% \\chi^{2}/%d=%.1f",
786 res[0], res[1], res[3], res[4], res[2], (int)res[5], res[6]/res[5]));
787 text.SetBit(TLatex::kTextNDC);
788 text.SetTextSize(0.04);
789 text.Paint();
790}
791
792void MHEffectiveOnTime::PaintProb(TH1 &h) const
793{
794 Double_t sum = 0;
795 Int_t n = 0;
796 for (int i=0; i<h.GetNbinsX(); i++)
797 if (h.GetBinContent(i+1)>0)
798 {
799 sum += h.GetBinContent(i+1);
800 n++;
801 }
802
803 if (n==0)
804 return;
805
806 TLatex text(0.47, 0.94, Form("\\bar{p} = %.1f%%", sum/n));
807 text.SetBit(TLatex::kTextNDC);
808 text.SetTextSize(0.04);
809 text.Paint();
810}
811
812void MHEffectiveOnTime::UpdateRightAxis(TH1 &h)
813{
814 const Double_t max = h.GetMaximum()*1.1;
815 if (max==0)
816 return;
817
818 h.SetNormFactor(h.Integral()*gPad->GetUymax()/max);
819
820 TGaxis *axis = (TGaxis*)gPad->FindObject("RightAxis");
821 if (!axis)
822 return;
823
824 axis->SetX1(gPad->GetUxmax());
825 axis->SetX2(gPad->GetUxmax());
826 axis->SetY1(gPad->GetUymin());
827 axis->SetY2(gPad->GetUymax());
828 axis->SetWmax(max);
829}
830
831// --------------------------------------------------------------------------
832//
833// Prepare painting the histograms
834//
835void MHEffectiveOnTime::Paint(Option_t *opt)
836{
837 TH1D *h=0;
838 TPaveStats *st=0;
839
840 TString o(opt);
841 if (o==(TString)"fit")
842 {
843 TVirtualPad *pad = gPad;
844
845 for (int x=0; x<2; x++)
846 for (int y=0; y<3; y++)
847 {
848 TVirtualPad *p=gPad->GetPad(x+1)->GetPad(y+1);
849 if (!(st = dynamic_cast<TPaveStats*>(p->GetPrimitive("stats"))))
850 continue;
851
852 if (st->GetOptStat()==11)
853 continue;
854
855 const Double_t y1 = st->GetY1NDC();
856 const Double_t y2 = st->GetY2NDC();
857 const Double_t x1 = st->GetX1NDC();
858 const Double_t x2 = st->GetX2NDC();
859
860 st->SetY1NDC((y2-y1)/3+y1);
861 st->SetX1NDC((x2-x1)/3+x1);
862 st->SetOptStat(11);
863 }
864
865 pad->GetPad(1)->cd(1);
866 if ((h = (TH1D*)gPad->FindObject("ProjDeltaT"/*fNameProjDeltaT*/)))
867 {
868 h = fH2DeltaT.ProjectionX("ProjDeltaT"/*fNameProjDeltaT*/, -1, -1, "E");
869 if (h->GetEntries()>0)
870 gPad->SetLogy();
871 }
872
873 pad->GetPad(2)->cd(1);
874 if ((h = (TH1D*)gPad->FindObject("ProjTheta"/*fNameProjTheta*/)))
875 fH2DeltaT.ProjectionY("ProjTheta"/*fNameProjTheta*/, -1, -1, "E");
876
877 if (!fIsFinalized)
878 FitThetaBins();
879 return;
880 }
881 if (o==(TString)"paint")
882 {
883 if ((h = (TH1D*)gPad->FindObject("ProjDeltaT"/*fNameProjDeltaT*/)))
884 {
885 Double_t res[7];
886 FitH(h, res, kTRUE);
887 PaintText(res);
888 }
889 return;
890 }
891
892 if (o==(TString)"timendf")
893 {
894 // UpdateRightAxis(fHTimeNDF);
895 // FIXME: first bin?
896 PaintProb(fHTimeProb);
897 }
898
899 if (o==(TString)"thetandf")
900 {
901 UpdateRightAxis(fHThetaNDF);
902 // FIXME: first bin?
903 PaintProb(fHThetaProb);
904 }
905
906 h=0;
907
908 Double_t range=-1;
909 if (o==(TString)"theta")
910 {
911 h = &fHThetaEffOn;
912 UpdateRightAxis(fHThetaLambda);
913 }
914 if (o==(TString)"time")
915 {
916 h = &fHTimeEffOn;
917 UpdateRightAxis(fHTimeLambda);
918 range = h->GetXaxis()->GetXmax() - h->GetXaxis()->GetXmin();
919 }
920
921 if (!h)
922 return;
923
924 Double_t error = 0;
925 for (int i=0; i<h->GetXaxis()->GetNbins(); i++)
926 error += h->GetBinError(i);
927
928 PaintText(h->Integral(), error, range);
929}
930
931void MHEffectiveOnTime::DrawRightAxis(const char *title)
932{
933 TGaxis *axis = new TGaxis(gPad->GetUxmax(), gPad->GetUymin(),
934 gPad->GetUxmax(), gPad->GetUymax(),
935 0, 1, 510, "+L");
936 axis->SetName("RightAxis");
937 axis->SetTitle(title);
938 axis->SetTitleOffset(0.9);
939 axis->SetTextColor(kGreen);
940 axis->CenterTitle();
941 axis->SetBit(kCanDelete);
942 axis->Draw();
943}
944
945// --------------------------------------------------------------------------
946//
947// Draw the histogram
948//
949void MHEffectiveOnTime::Draw(Option_t *opt)
950{
951 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
952 pad->SetBorderMode(0);
953
954 AppendPad("fit");
955
956 pad->Divide(2, 1, 1e-10, 1e-10);
957
958 TH1 *h;
959
960 pad->cd(1);
961 gPad->SetBorderMode(0);
962 gPad->Divide(1, 3, 1e-10, 1e-10);
963 pad->GetPad(1)->cd(1);
964 gPad->SetBorderMode(0);
965 h = fH2DeltaT.ProjectionX("ProjDeltaT"/*fNameProjDeltaT*/, -1, -1, "E");
966 h->SetTitle("Distribution of \\Delta t [s]");
967 h->SetXTitle("\\Delta t [s]");
968 h->SetYTitle("Counts");
969 h->SetDirectory(NULL);
970 h->SetMarkerStyle(kFullDotMedium);
971 h->SetBit(kCanDelete);
972 h->Draw();
973 AppendPad("paint");
974
975 pad->GetPad(1)->cd(2);
976 gPad->SetBorderMode(0);
977 fHTimeProb.Draw();
978 AppendPad("timendf");
979 //fHTimeNDF.Draw("same");
980 //DrawRightAxis("NDF");
981
982 pad->GetPad(1)->cd(3);
983 gPad->SetBorderMode(0);
984 fHTimeEffOn.Draw();
985 AppendPad("time");
986 fHTimeLambda.Draw("same");
987 DrawRightAxis("\\lambda [s^{-1}]");
988
989 pad->cd(2);
990 gPad->SetBorderMode(0);
991 gPad->Divide(1, 3, 1e-10, 1e-10);
992
993 pad->GetPad(2)->cd(1);
994 gPad->SetBorderMode(0);
995 h = fH2DeltaT.ProjectionY("ProjTheta"/*fNameProjTheta*/, -1, -1, "E");
996 h->SetTitle("Distribution of \\Theta [\\circ]");
997 h->SetXTitle("\\Theta [\\circ]");
998 h->SetYTitle("Counts");
999 h->SetDirectory(NULL);
1000 h->SetMarkerStyle(kFullDotMedium);
1001 h->SetBit(kCanDelete);
1002 h->GetYaxis()->SetTitleOffset(1.1);
1003 h->Draw();
1004
1005 pad->GetPad(2)->cd(2);
1006 gPad->SetBorderMode(0);
1007 fHThetaProb.Draw();
1008 AppendPad("thetandf");
1009 fHThetaNDF.Draw("same");
1010 DrawRightAxis("NDF");
1011
1012 pad->GetPad(2)->cd(3);
1013 gPad->SetBorderMode(0);
1014 fHThetaEffOn.Draw();
1015 AppendPad("theta");
1016 fHThetaLambda.Draw("same");
1017 DrawRightAxis("\\lambda [s^{-1}]");
1018}
1019
1020// --------------------------------------------------------------------------
1021//
1022// The following resources are available:
1023//
1024// MHEffectiveOnTime.FistBin: 3
1025// MHEffectiveOnTime.NumEvents: 12000
1026//
1027Int_t MHEffectiveOnTime::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1028{
1029 Bool_t rc = kFALSE;
1030 if (IsEnvDefined(env, prefix, "FirstBin", print))
1031 {
1032 rc = kTRUE;
1033 SetFirstBin(GetEnvValue(env, prefix, "FirstBin", (Int_t)fFirstBin));
1034 }
1035 if (IsEnvDefined(env, prefix, "NumEvents", print))
1036 {
1037 rc = kTRUE;
1038 SetNumEvents(GetEnvValue(env, prefix, "NumEvents", (Int_t)fNumEvents));
1039 }
1040 return rc;
1041}
Note: See TracBrowser for help on using the repository browser.