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

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