source: branches/MarsMoreSimulationTruth/mhflux/MHEffectiveOnTime.cc@ 19085

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