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

Last change on this file since 9029 was 8946, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 32.8 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 "MParameters.h"
283#include "MPointingPos.h"
284#include "MRawRunHeader.h"
285
286#include "MBinning.h"
287#include "MParList.h"
288
289#include "MLog.h"
290#include "MLogManip.h"
291
292ClassImp(MHEffectiveOnTime);
293
294using namespace std;
295
296// --------------------------------------------------------------------------
297//
298// Default Constructor. It initializes all histograms.
299//
300MHEffectiveOnTime::MHEffectiveOnTime(const char *name, const char *title)
301 : fPointPos(0), fTime(0), fParam(0), fIsFinalized(kFALSE),
302 fNumEvents(200*60), fFirstBin(3), fTotalTime(-1)
303 //fNumEvents(2*60), fFirstBin(1)
304{
305 //
306 // set the name and title of this object
307 //
308 fName = name ? name : "MHEffectiveOnTime";
309 fTitle = title ? title : "Histogram to determin effective On-Time vs Time and Zenith Angle";
310
311 // Main histogram
312 fH2DeltaT.SetName("DeltaT");
313 fH2DeltaT.SetXTitle("\\Delta t [s]");
314 fH2DeltaT.SetYTitle("\\Theta [\\circ]");
315 fH2DeltaT.SetZTitle("Count");
316 fH2DeltaT.UseCurrentStyle();
317 fH2DeltaT.SetDirectory(NULL);
318
319 // Main histogram
320 fH1DeltaT.SetName("DeltaT");
321 fH1DeltaT.SetXTitle("\\Delta t [s]");
322 fH1DeltaT.SetYTitle("Counts");
323 fH1DeltaT.UseCurrentStyle();
324 fH1DeltaT.SetDirectory(NULL);
325
326 // effective on time versus theta
327 fHThetaEffOn.SetName("EffOnTheta");
328 fHThetaEffOn.SetTitle("Effective On Time T_{eff}");
329 fHThetaEffOn.SetXTitle("\\Theta [\\circ]");
330 fHThetaEffOn.SetYTitle("T_{eff} [s]");
331 fHThetaEffOn.UseCurrentStyle();
332 fHThetaEffOn.SetDirectory(NULL);
333 fHThetaEffOn.GetYaxis()->SetTitleOffset(1.2);
334 fHThetaEffOn.GetYaxis()->SetTitleColor(kBlue);
335 fHThetaEffOn.SetLineColor(kBlue);
336 //fHEffOn.Sumw2();
337
338 // effective on time versus time
339 fHTimeEffOn.SetName("EffOnTime");
340 fHTimeEffOn.SetTitle("Effective On Time T_{eff}");
341 fHTimeEffOn.SetXTitle("Time");
342 fHTimeEffOn.SetYTitle("T_{eff} [s]");
343 fHTimeEffOn.UseCurrentStyle();
344 fHTimeEffOn.SetDirectory(NULL);
345 fHTimeEffOn.GetYaxis()->SetTitleOffset(1.2);
346 fHTimeEffOn.GetXaxis()->SetLabelSize(0.033);
347 fHTimeEffOn.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
348 fHTimeEffOn.GetXaxis()->SetTimeDisplay(1);
349 fHTimeEffOn.GetYaxis()->SetTitleColor(kBlue);
350 fHTimeEffOn.SetLineColor(kBlue);
351
352 // chi2 probability
353 fHThetaProb.SetName("ProbTheta");
354 fHThetaProb.SetTitle("\\chi^{2} Probability of Fit");
355 fHThetaProb.SetXTitle("\\Theta [\\circ]");
356 fHThetaProb.SetYTitle("p [%]");
357 fHThetaProb.UseCurrentStyle();
358 fHThetaProb.SetDirectory(NULL);
359 fHThetaProb.GetYaxis()->SetTitleOffset(1.2);
360 fHThetaProb.SetMaximum(101);
361 fHThetaProb.GetYaxis()->SetTitleColor(kBlue);
362 fHThetaProb.SetLineColor(kBlue);
363
364 // chi2 probability
365 fHTimeProb.SetName("ProbTime");
366 fHTimeProb.SetTitle("\\chi^{2} Probability of Fit");
367 fHTimeProb.SetXTitle("Time");
368 fHTimeProb.SetYTitle("p [%]");
369 fHTimeProb.UseCurrentStyle();
370 fHTimeProb.SetDirectory(NULL);
371 fHTimeProb.GetYaxis()->SetTitleOffset(1.2);
372 fHTimeProb.GetXaxis()->SetLabelSize(0.033);
373 fHTimeProb.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
374 fHTimeProb.GetXaxis()->SetTimeDisplay(1);
375 fHTimeProb.SetMaximum(101);
376 fHTimeProb.GetYaxis()->SetTitleColor(kBlue);
377 fHTimeProb.SetLineColor(kBlue);
378
379 // lambda versus theta
380 fHThetaLambda.SetName("LambdaTheta");
381 fHThetaLambda.SetTitle("Slope (Rate) vs Theta");
382 fHThetaLambda.SetXTitle("\\Theta [\\circ]");
383 fHThetaLambda.SetYTitle("\\lambda [s^{-1}]");
384 fHThetaLambda.UseCurrentStyle();
385 fHThetaLambda.SetDirectory(NULL);
386 fHThetaLambda.SetLineColor(kGreen);
387
388 // lambda versus time
389 fHTimeLambda.SetName("LambdaTime");
390 fHTimeLambda.SetTitle("Slope (Rate) vs Time");
391 fHTimeLambda.SetXTitle("\\Time [\\circ]");
392 fHTimeLambda.SetYTitle("\\lambda [s^{-1}]");
393 fHTimeLambda.UseCurrentStyle();
394 fHTimeLambda.SetDirectory(NULL);
395 fHTimeLambda.GetYaxis()->SetTitleOffset(1.2);
396 fHTimeLambda.GetXaxis()->SetLabelSize(0.033);
397 fHTimeLambda.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
398 fHTimeLambda.GetXaxis()->SetTimeDisplay(1);
399 fHTimeLambda.SetLineColor(kGreen);
400
401 // NDoF versus theta
402 fHThetaNDF.SetName("NDofTheta");
403 fHThetaNDF.SetTitle("Number of Degrees of freedom vs Theta");
404 fHThetaNDF.SetXTitle("\\Theta [\\circ]");
405 fHThetaNDF.SetYTitle("NDoF [#]");
406 fHThetaNDF.UseCurrentStyle();
407 fHThetaNDF.SetDirectory(NULL);
408 fHThetaNDF.SetLineColor(kGreen);
409
410 // NDoF versus time
411 /*
412 fHTimeNDF.SetName("NDofTime");
413 fHTimeNDF.SetTitle("Number of Degrees of freedom vs Time");
414 fHTimeNDF.SetXTitle("Time");
415 fHTimeNDF.SetYTitle("NDoF [#]");
416 fHTimeNDF.UseCurrentStyle();
417 fHTimeNDF.SetDirectory(NULL);
418 fHTimeNDF.GetYaxis()->SetTitleOffset(1.2);
419 fHTimeNDF.GetXaxis()->SetLabelSize(0.033);
420 fHTimeNDF.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
421 fHTimeNDF.GetXaxis()->SetTimeDisplay(1);
422 fHTimeNDF.SetLineColor(kBlue);
423 */
424 // setup binning
425 MBinning btheta("BinningTheta");
426 btheta.SetEdgesASin(67, -0.005, 0.665);
427
428 MBinning btime("BinningDeltaT");
429 btime.SetEdges(50, 0, 0.1);
430
431 MH::SetBinning(&fH2DeltaT, &btime, &btheta);
432
433 btime.Apply(fH1DeltaT);
434
435 btheta.Apply(fHThetaEffOn);
436 btheta.Apply(fHThetaLambda);
437 btheta.Apply(fHThetaNDF);
438 btheta.Apply(fHThetaProb);
439 //btheta.Apply(fHChi2);
440}
441
442// --------------------------------------------------------------------------
443//
444// Set the binnings and prepare the filling of the histogram
445//
446Bool_t MHEffectiveOnTime::SetupFill(const MParList *plist)
447{
448 fPointPos = (MPointingPos*)plist->FindObject("MPointingPos");
449 if (!fPointPos)
450 {
451 *fLog << err << dbginf << "MPointingPos not found... aborting." << endl;
452 return kFALSE;
453 }
454
455 // FIXME: Remove const-qualifier from base-class!
456 fTime = (MTime*)const_cast<MParList*>(plist)->FindCreateObj("MTime", "MTimeEffectiveOnTime");
457 if (!fTime)
458 return kFALSE;
459 fParam = (MParameterDerr*)const_cast<MParList*>(plist)->FindCreateObj("MParameterDerr", "MEffectiveOnTime");
460 if (!fParam)
461 return kFALSE;
462
463 const MBinning* binsdtime = (MBinning*)plist->FindObject("BinningDeltaT");
464 const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
465 if (binsdtime)
466 binsdtime->Apply(fH1DeltaT);
467 if (binstheta)
468 {
469 binstheta->Apply(fHThetaEffOn);
470 binstheta->Apply(fHThetaLambda);
471 binstheta->Apply(fHThetaNDF);
472 binstheta->Apply(fHThetaProb);
473 //binstheta->Apply(fHChi2);
474 }
475 if (binstheta && binsdtime)
476 SetBinning(&fH2DeltaT, binsdtime, binstheta);
477
478 fTotalTime = 0;
479
480 return kTRUE;
481}
482
483Bool_t MHEffectiveOnTime::ReInit(MParList *pList)
484{
485 MRawRunHeader *header = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
486 if (!header)
487 {
488 *fLog << err << "MRawRunHeader not found... aborting." << endl;
489 return kFALSE;
490 }
491
492 fTotalTime += header->GetRunLength();
493
494 return kTRUE;
495}
496
497// --------------------------------------------------------------------------
498//
499// Fit a single Delta-T distribution. See source code for more details
500//
501Bool_t MHEffectiveOnTime::FitH(TH1D *h, Double_t *res, Bool_t paint) const
502{
503 // Count also events in under-/overflowbins
504 const Double_t Nm = h->Integral(0, h->GetNbinsX()+1);
505
506 // FIXME: Do fit only if contents of bin has changed
507 if (Nm<=0 || h->GetBinContent(1)<=0)
508 return kFALSE;
509
510 // determine range (yq[0], yq[1]) of time differences
511 // where fit should be performed;
512 // require a fraction >=xq[0] of all entries to lie below yq[0]
513 // and a fraction <=xq[1] of all entries to lie below yq[1];
514 // within the range (yq[0], yq[1]) there must be no empty bin;
515 // choose pedestrian approach as long as GetQuantiles is not available
516 Double_t xq[2] = { 0.6, 0.95 }; // previously 0.99
517 Double_t yq[2];
518 h->GetQuantiles(2, yq, xq);
519
520 //
521 // Determine a good starting value for the slope
522 //
523 const TAxis &axe = *h->GetXaxis();
524 const UInt_t ibin = axe.FindFixBin(yq[1]);
525 const Double_t x1 = axe.GetBinCenter(ibin<=fFirstBin?fFirstBin+1:ibin);
526 const Double_t x0 = axe.GetBinCenter(fFirstBin);
527 const Double_t y1 = h->GetBinContent(ibin)>1 ? TMath::Log(h->GetBinContent(ibin)) : 0;
528 const Double_t y0 = TMath::Log(h->GetBinContent(fFirstBin));
529
530 // Estimated slope
531 const Float_t m = -(y1-y0)/(x1-x0);
532
533 //
534 // Setup exponential function for the fit:
535 //
536 // parameter 0 = rate [Hz]
537 // parameter 1 = normalization
538 //
539 TF1 func("Exp", " exp([1]-[0]*x)");
540
541 func.SetParameter(0, m); // Hz
542 func.SetParameter(1, log(h->GetBinContent(1))); // Hz
543
544 // We set a limit to make sure that almost empty histograms which
545 // are fitted dont't produce hang ups or crashes
546 func.SetParLimits(0, 0, 15000); // Hz
547
548 // options : N do not store the function, do not draw
549 // I use integral of function in bin rather than value at bin center
550 // R use the range specified in the function range
551 // Q quiet mode
552 // L Use log-likelihood (better for low statistics)
553 h->Fit(&func, "NIQEL", "", h->GetBinLowEdge(fFirstBin)/*yq[0]*/, yq[1]);
554
555 const Double_t chi2 = func.GetChisquare();
556 const Double_t prob = func.GetProb();
557 const Int_t NDF = func.GetNDF();
558
559 // was fit successful ?
560 const Bool_t ok = prob>0.001; //NDF>0 && chi2<3*NDF;
561
562 if (paint)
563 {
564 func.SetLineWidth(2);
565 func.SetLineColor(ok ? kGreen : kRed);
566 func.Paint("same");
567 }
568
569 // The effective on time is the "real rate" (slope of the exponential)
570 // divided by the total number of events (histogram integral including
571 // under- and overflows)
572 const Double_t lambda = func.GetParameter(0);
573 const Double_t dldl = func.GetParError(0)*func.GetParError(0);
574 const Double_t teff = lambda==0 ? 0 : Nm / lambda;
575 const Double_t dteff = lambda==0 ? 0 : teff * TMath::Sqrt(dldl/(lambda*lambda) + 1.0/Nm);
576 const Double_t dl = TMath::Sqrt(dldl);
577
578 // the effective on time is Nm/lambda
579 res[0] = teff;
580 res[1] = dteff;
581
582 // plot chi2-probability of fit
583 res[2] = prob*100;
584
585 // lambda of fit
586 res[3] = lambda;
587 res[4] = dl;
588
589 // NDoF of fit
590 res[5] = NDF;
591
592 // Chi2
593 res[6] = chi2;
594
595 return ok;
596}
597
598// --------------------------------------------------------------------------
599//
600// Fit a all bins of the distribution in theta. Store the result in the
601// Theta-Histograms
602//
603void MHEffectiveOnTime::FitThetaBins()
604{
605 fHThetaEffOn.Reset();
606 fHThetaProb.Reset();
607 fHThetaLambda.Reset();
608 fHThetaNDF.Reset();
609
610 // Use a random name to make sure the object is unique
611 const TString name = Form("CalcTheta%d", (UInt_t)gRandom->Uniform(999999999));
612
613 // nbins = number of Theta bins
614 const Int_t nbins = fH2DeltaT.GetNbinsY();
615
616 TH1D *h=0;
617 for (int i=1; i<=nbins; i++)
618 {
619 // TH1D &h = *hist->ProjectionX("Calc-theta", i, i);
620 h = fH2DeltaT.ProjectionX(name, i, i, "E");
621
622 Double_t res[7] = {0, 0, 0, 0, 0, 0, 0};
623 //if (!FitH(h, res))
624 // continue;
625 FitH(h, res);
626
627 if (res[0]==0)
628 continue;
629
630 // the effective on time is Nm/lambda
631 fHThetaEffOn.SetBinContent(i, res[0]);
632 fHThetaEffOn.SetBinError (i, res[1]);
633
634 // plot chi2-probability of fit
635 fHThetaProb.SetBinContent(i, res[2]);
636
637 // plot chi2/NDF of fit
638 //fHChi2.SetBinContent(i, res[3]);
639
640 // lambda of fit
641 fHThetaLambda.SetBinContent(i, res[3]);
642 fHThetaLambda.SetBinError (i, res[4]);
643
644 // NDoF of fit
645 fHThetaNDF.SetBinContent(i, res[5]);
646
647 // Rdead (from fit) is the fraction from real time lost by the dead time
648 //fHRdead.SetBinContent(i, Rdead);
649 //fHRdead.SetBinError (i,dRdead);
650 }
651
652 // Histogram is reused via gROOT->FindObject()
653 // Need to be deleted only once
654 if (h)
655 delete h;
656}
657
658// --------------------------------------------------------------------------
659//
660// Fit the single-time-bin histogram. Store the result in the
661// Time-Histograms
662//
663void MHEffectiveOnTime::FitTimeBin()
664{
665 //
666 // Fit histogram
667 //
668 Double_t res[7];
669 if (!FitH(&fH1DeltaT, res))
670 return;
671
672 // Reset Histogram
673 fH1DeltaT.Reset();
674
675 //
676 // Prepare Histogram
677 //
678
679 // Get number of bins
680 const Int_t n = fHTimeEffOn.GetNbinsX();
681
682 // Enhance binning
683 MBinning bins;
684 bins.SetEdges(fHTimeEffOn, 'x');
685 bins.AddEdge(fLastTime.GetAxisTime());
686 bins.Apply(fHTimeEffOn);
687 bins.Apply(fHTimeProb);
688 bins.Apply(fHTimeLambda);
689 //bins.Apply(fHTimeNDF);
690
691 //
692 // Fill histogram
693 //
694 fHTimeEffOn.SetBinContent(n, res[0]);
695 fHTimeEffOn.SetBinError(n, res[1]);
696
697 fHTimeProb.SetBinContent(n, res[2]);
698
699 fHTimeLambda.SetBinContent(n, res[3]);
700 fHTimeLambda.SetBinError(n, res[4]);
701
702 //fHTimeNDF.SetBinContent(n, res[5]);
703
704 //
705 // Now prepare output
706 //
707 fParam->SetVal(res[0], res[1]);
708 fParam->SetReadyToSave();
709
710 *fTime = fLastTime;
711
712 // Include the current event
713 fTime->Plus1ns();
714
715 *fLog << fLastTime << ": Val=" << res[0] << " Err=" << res[1] << endl;
716}
717
718// --------------------------------------------------------------------------
719//
720// Fill the histogram
721//
722Bool_t MHEffectiveOnTime::Fill(const MParContainer *par, const Stat_t w)
723{
724 const MTime *time = dynamic_cast<const MTime*>(par);
725 if (!time)
726 {
727 *fLog << err << "ERROR - MHEffectiveOnTime::Fill without argument or container doesn't inherit from MTime... abort." << endl;
728 return kFALSE;
729 }
730
731 //
732 // If this is the first call we have to initialize the time-histogram
733 //
734 if (fLastTime==MTime())
735 {
736 MBinning bins;
737 bins.SetEdges(1, time->GetAxisTime()-fNumEvents/200, time->GetAxisTime());
738 bins.Apply(fHTimeEffOn);
739 bins.Apply(fHTimeProb);
740 bins.Apply(fHTimeLambda);
741 //bins.Apply(fHTimeNDF);
742
743 fParam->SetVal(0, 0);
744 fParam->SetReadyToSave();
745
746 *fTime = *time;
747
748 // Make this 1ns before the first event!
749 fTime->Minus1ns();
750 }
751
752 //
753 // Fill time difference into the histograms
754 //
755 const Double_t dt = *time-fLastTime;
756 fLastTime = *time;
757
758 fH2DeltaT.Fill(dt, fPointPos->GetZd(), w);
759 fH1DeltaT.Fill(dt, w);
760
761 //
762 // If we reached the event number limit for the time-bins fit the
763 // histogram - if it fails try again when 1.6% more events available
764 //
765 const UInt_t n = (UInt_t)fH1DeltaT.GetEntries();
766 if (n>=fNumEvents && n%(fNumEvents/60)==0)
767 FitTimeBin();
768
769 return kTRUE;
770}
771
772// --------------------------------------------------------------------------
773//
774// Fit the theta projections of the 2D histogram and the 1D Delta-T
775// distribution
776//
777Bool_t MHEffectiveOnTime::Finalize()
778{
779 FitThetaBins();
780 FitTimeBin();
781
782 fIsFinalized = kTRUE;
783
784 return kTRUE;
785}
786
787// --------------------------------------------------------------------------
788//
789// Paint the integral and the error on top of the histogram
790//
791void MHEffectiveOnTime::PaintText(Double_t val, Double_t error, Double_t range) const
792{
793 TLatex text;
794 text.SetBit(TLatex::kTextNDC);
795 text.SetTextSize(0.04);
796
797 TString txt = Form("T_{eff} = %.1fs \\pm %.1fs", val, error);
798 if (range>0)
799 txt += Form(" T_{axe} = %.1fs", range);
800 if (fTotalTime>0)
801 txt += Form(" T_{max} = %.1fs", fTotalTime);
802
803 text.SetText(0.35, 0.94, txt);
804 text.Paint();
805}
806
807void MHEffectiveOnTime::PaintText(Double_t *res) const
808{
809 TLatex text(0.27, 0.94, Form("T_{eff}=%.1fs\\pm%.1fs \\lambda=%.1f\\pm%.1fHz p=%.1f%% \\chi^{2}/%d=%.1f",
810 res[0], res[1], res[3], res[4], res[2], (int)res[5], res[6]/res[5]));
811 text.SetBit(TLatex::kTextNDC);
812 text.SetTextSize(0.04);
813 text.Paint();
814}
815
816void MHEffectiveOnTime::PaintProb(TH1 &h) const
817{
818 Double_t sum = 0;
819 Int_t n = 0;
820 for (int i=0; i<h.GetNbinsX(); i++)
821 if (h.GetBinContent(i+1)>0)
822 {
823 sum += h.GetBinContent(i+1);
824 n++;
825 }
826
827 if (n==0)
828 return;
829
830 TLatex text(0.47, 0.94, Form("\\bar{p} = %.1f%%", sum/n));
831 text.SetBit(TLatex::kTextNDC);
832 text.SetTextSize(0.04);
833 text.Paint();
834}
835
836void MHEffectiveOnTime::UpdateRightAxis(TH1 &h)
837{
838 const Double_t max = h.GetMaximum()*1.1;
839 if (max==0)
840 return;
841
842 h.SetNormFactor(h.Integral()*gPad->GetUymax()/max);
843
844 TGaxis *axis = (TGaxis*)gPad->FindObject("RightAxis");
845 if (!axis)
846 return;
847
848 axis->SetX1(gPad->GetUxmax());
849 axis->SetX2(gPad->GetUxmax());
850 axis->SetY1(gPad->GetUymin());
851 axis->SetY2(gPad->GetUymax());
852 axis->SetWmax(max);
853}
854
855// --------------------------------------------------------------------------
856//
857// Prepare painting the histograms
858//
859void MHEffectiveOnTime::Paint(Option_t *opt)
860{
861 TH1D *h=0;
862 TPaveStats *st=0;
863
864 TString o(opt);
865 if (o==(TString)"fit")
866 {
867 TVirtualPad *pad = gPad;
868
869 for (int x=0; x<2; x++)
870 for (int y=0; y<3; y++)
871 {
872 TVirtualPad *p=gPad->GetPad(x+1)->GetPad(y+1);
873 if (!(st = dynamic_cast<TPaveStats*>(p->GetPrimitive("stats"))))
874 continue;
875
876 if (st->GetOptStat()==11)
877 continue;
878
879 const Double_t y1 = st->GetY1NDC();
880 const Double_t y2 = st->GetY2NDC();
881 const Double_t x1 = st->GetX1NDC();
882 const Double_t x2 = st->GetX2NDC();
883
884 st->SetY1NDC((y2-y1)/3+y1);
885 st->SetX1NDC((x2-x1)/3+x1);
886 st->SetOptStat(11);
887 }
888
889 pad->GetPad(1)->cd(1);
890 if ((h = (TH1D*)gPad->FindObject("ProjDeltaT"/*fNameProjDeltaT*/)))
891 {
892 h = fH2DeltaT.ProjectionX("ProjDeltaT"/*fNameProjDeltaT*/, -1, -1, "E");
893 if (h->GetEntries()>0)
894 gPad->SetLogy();
895 }
896
897 pad->GetPad(2)->cd(1);
898 if ((h = (TH1D*)gPad->FindObject("ProjTheta"/*fNameProjTheta*/)))
899 fH2DeltaT.ProjectionY("ProjTheta"/*fNameProjTheta*/, -1, -1, "E");
900
901 if (!fIsFinalized)
902 FitThetaBins();
903 return;
904 }
905 if (o==(TString)"paint")
906 {
907 if ((h = (TH1D*)gPad->FindObject("ProjDeltaT"/*fNameProjDeltaT*/)))
908 {
909 Double_t res[7];
910 FitH(h, res, kTRUE);
911 PaintText(res);
912 }
913 return;
914 }
915
916 if (o==(TString)"timendf")
917 {
918 // UpdateRightAxis(fHTimeNDF);
919 // FIXME: first bin?
920 PaintProb(fHTimeProb);
921 }
922
923 if (o==(TString)"thetandf")
924 {
925 UpdateRightAxis(fHThetaNDF);
926 // FIXME: first bin?
927 PaintProb(fHThetaProb);
928 }
929
930 h=0;
931
932 Double_t range=-1;
933 if (o==(TString)"theta")
934 {
935 h = &fHThetaEffOn;
936 UpdateRightAxis(fHThetaLambda);
937 }
938 if (o==(TString)"time")
939 {
940 h = &fHTimeEffOn;
941 UpdateRightAxis(fHTimeLambda);
942 range = h->GetXaxis()->GetXmax() - h->GetXaxis()->GetXmin();
943 }
944
945 if (!h)
946 return;
947
948 Double_t error = 0;
949 for (int i=0; i<h->GetXaxis()->GetNbins(); i++)
950 error += h->GetBinError(i);
951
952 PaintText(h->Integral(), error, range);
953}
954
955void MHEffectiveOnTime::DrawRightAxis(const char *title)
956{
957 TGaxis *axis = new TGaxis(gPad->GetUxmax(), gPad->GetUymin(),
958 gPad->GetUxmax(), gPad->GetUymax(),
959 0, 1, 510, "+L");
960 axis->SetName("RightAxis");
961 axis->SetTitle(title);
962 axis->SetTitleOffset(0.9);
963 axis->SetTextColor(kGreen);
964 axis->CenterTitle();
965 axis->SetBit(kCanDelete);
966 axis->Draw();
967}
968
969// --------------------------------------------------------------------------
970//
971// Draw the histogram
972//
973void MHEffectiveOnTime::Draw(Option_t *opt)
974{
975 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
976 pad->SetBorderMode(0);
977
978 AppendPad("fit");
979
980 pad->Divide(2, 1, 1e-10, 1e-10);
981
982 TH1 *h;
983
984 pad->cd(1);
985 gPad->SetBorderMode(0);
986 gPad->Divide(1, 3, 1e-10, 1e-10);
987 pad->GetPad(1)->cd(1);
988 gPad->SetBorderMode(0);
989 h = fH2DeltaT.ProjectionX("ProjDeltaT"/*fNameProjDeltaT*/, -1, -1, "E");
990 h->SetTitle("Distribution of \\Delta t [s]");
991 h->SetXTitle("\\Delta t [s]");
992 h->SetYTitle("Counts");
993 h->SetDirectory(NULL);
994 h->SetMarkerStyle(kFullDotMedium);
995 h->SetBit(kCanDelete);
996 h->Draw();
997 AppendPad("paint");
998
999 pad->GetPad(1)->cd(2);
1000 gPad->SetBorderMode(0);
1001 fHTimeProb.Draw();
1002 AppendPad("timendf");
1003 //fHTimeNDF.Draw("same");
1004 //DrawRightAxis("NDF");
1005
1006 pad->GetPad(1)->cd(3);
1007 gPad->SetBorderMode(0);
1008 fHTimeEffOn.Draw();
1009 AppendPad("time");
1010 fHTimeLambda.Draw("same");
1011 DrawRightAxis("\\lambda [s^{-1}]");
1012
1013 pad->cd(2);
1014 gPad->SetBorderMode(0);
1015 gPad->Divide(1, 3, 1e-10, 1e-10);
1016
1017 pad->GetPad(2)->cd(1);
1018 gPad->SetBorderMode(0);
1019 h = fH2DeltaT.ProjectionY("ProjTheta"/*fNameProjTheta*/, -1, -1, "E");
1020 h->SetTitle("Distribution of \\Theta [\\circ]");
1021 h->SetXTitle("\\Theta [\\circ]");
1022 h->SetYTitle("Counts");
1023 h->SetDirectory(NULL);
1024 h->SetMarkerStyle(kFullDotMedium);
1025 h->SetBit(kCanDelete);
1026 h->GetYaxis()->SetTitleOffset(1.1);
1027 h->Draw();
1028
1029 pad->GetPad(2)->cd(2);
1030 gPad->SetBorderMode(0);
1031 fHThetaProb.Draw();
1032 AppendPad("thetandf");
1033 fHThetaNDF.Draw("same");
1034 DrawRightAxis("NDF");
1035
1036 pad->GetPad(2)->cd(3);
1037 gPad->SetBorderMode(0);
1038 fHThetaEffOn.Draw();
1039 AppendPad("theta");
1040 fHThetaLambda.Draw("same");
1041 DrawRightAxis("\\lambda [s^{-1}]");
1042}
1043
1044// --------------------------------------------------------------------------
1045//
1046// The following resources are available:
1047//
1048// MHEffectiveOnTime.FistBin: 3
1049// MHEffectiveOnTime.NumEvents: 12000
1050//
1051Int_t MHEffectiveOnTime::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1052{
1053 Bool_t rc = kFALSE;
1054 if (IsEnvDefined(env, prefix, "FirstBin", print))
1055 {
1056 rc = kTRUE;
1057 SetFirstBin(GetEnvValue(env, prefix, "FirstBin", (Int_t)fFirstBin));
1058 }
1059 if (IsEnvDefined(env, prefix, "NumEvents", print))
1060 {
1061 rc = kTRUE;
1062 SetNumEvents(GetEnvValue(env, prefix, "NumEvents", (Int_t)fNumEvents));
1063 }
1064 return rc;
1065}
Note: See TracBrowser for help on using the repository browser.