source: trunk/MagicSoft/Mars/mcalib/MHGausEvents.cc@ 3700

Last change on this file since 3700 was 3677, checked in by gaug, 21 years ago
*** empty log message ***
File size: 26.1 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): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MHGausEvents
28//
29// A base class for events which are believed to follow a Gaussian distribution
30// with time, e.g. calibration events, observables containing white noise, ...
31//
32// MHGausEvents derives from MH, thus all features of MH can be used by a class
33// deriving from MHGausEvents, especially the filling functions.
34//
35// The central objects are:
36//
37// 1) The TH1F fHGausHist:
38// ====================
39//
40// It is created with a default name and title and resides in directory NULL.
41// - Any class deriving from MHGausEvents needs to apply a binning to fHGausHist
42// (e.g. by setting the variables fNbins, fFirst, fLast and calling the function
43// InitBins() or by directly calling fHGausHist.SetBins(..) )
44// - The histogram is filled with the functions FillHist() or FillHistAndArray().
45// - The histogram can be fitted with the function FitGaus(). This involves stripping
46// of all zeros at the lower and upper end of the histogram and re-binning to
47// a new number of bins, specified in the variable fBinsAfterStripping.
48// - The fit result's probability is compared to a reference probability fProbLimit
49// The NDF is compared to fNDFLimit and a check is made whether results are NaNs.
50// Accordingly, a flag IsGausFitOK() is set.
51// - One can repeat the fit within a given amount of sigmas from the previous mean
52// (specified by the variables fPickupLimit and fBlackoutLimit) with the function RepeatFit()
53// - One can completely skip the fitting to set mean, sigma and its errors directly
54// from the histograms with the function BypassFit()
55//
56// 2) The TArrayF fEvents:
57// ==========================
58//
59// It is created with 0 entries and not expanded unless FillArray() or FillHistAndArray()
60// are called.
61// - A first call to FillArray() or FillHistAndArray() initializes fEvents by default
62// to 512 entries.
63// - Any later call to FillArray() or FillHistAndArray() fills up the array.
64// Reaching the limit, the array is expanded by a factor 2.
65// - The array can be fourier-transformed into the array fPowerSpectrum.
66// Note that any FFT accepts only number of events which are a power of 2.
67// Thus, fEvents is cut to the next power of 2 smaller than its actual number of entries.
68// Be aware that you might lose information at the end of your analysis.
69// - Calling the function CreateFourierSpectrum() creates the array fPowerSpectrum
70// and its projection fHPowerProbability which in turn is fit to an exponential.
71// - The fit result's probability is compared to a referenc probability fProbLimit
72// and accordingly the flag IsExpFitOK() is set.
73// - The flag IsFourierSpectrumOK() is set accordingly to IsExpFitOK().
74// Later, a closer check will be installed.
75// - You can display all arrays by calls to: CreateGraphEvents() and
76// CreateGraphPowerSpectrum() and successive calls to the corresponding Getters.
77//
78//////////////////////////////////////////////////////////////////////////////
79#include "MHGausEvents.h"
80
81#include <TH1.h>
82#include <TF1.h>
83#include <TGraph.h>
84#include <TPad.h>
85#include <TVirtualPad.h>
86#include <TCanvas.h>
87#include <TStyle.h>
88
89#include "MFFT.h"
90#include "MArray.h"
91
92#include "MH.h"
93
94#include "MLog.h"
95#include "MLogManip.h"
96
97ClassImp(MHGausEvents);
98
99using namespace std;
100
101const Int_t MHGausEvents::fgBinsAfterStripping = 40;
102const Float_t MHGausEvents::fgBlackoutLimit = 5.;
103const Int_t MHGausEvents::fgNDFLimit = 2;
104const Float_t MHGausEvents::fgPickupLimit = 5.;
105const Float_t MHGausEvents::fgProbLimit = 0.001;
106const Int_t MHGausEvents::fgPowerProbabilityBins = 20;
107// --------------------------------------------------------------------------
108//
109// Default Constructor.
110// Sets:
111// - the default Probability Bins for fPowerProbabilityBins (fgPowerProbabilityBins)
112// - the default Probability Limit for fProbLimit (fgProbLimit)
113// - the default NDF Limit for fNDFLimit (fgNDFLimit)
114// - the default number for fPickupLimit (fgPickupLimit)
115// - the default number for fBlackoutLimit (fgBlackoutLimit)
116// - the default number of bins after stripping for fBinsAfterStipping (fgBinsAfterStipping)
117// - the default name of the fHGausHist ("HGausHist")
118// - the default title of the fHGausHist ("Histogram of Events with Gaussian Distribution")
119// - the default directory of the fHGausHist (NULL)
120// - the default for fNbins (100)
121// - the default for fFirst (0.)
122// - the default for fLast (100.)
123//
124// Initializes:
125// - fEvents to 0 entries
126// - fHGausHist()
127// - all other pointers to NULL
128// - all variables to 0.
129// - all flags to kFALSE
130//
131MHGausEvents::MHGausEvents(const char *name, const char *title)
132 : fEventFrequency(0), fHPowerProbability(NULL),
133 fPowerSpectrum(NULL),
134 fGraphEvents(NULL), fGraphPowerSpectrum(NULL),
135 fEvents(0), fFGausFit(NULL), fFExpFit(NULL),
136 fFirst(0.), fHGausHist(), fLast(100.),
137 fNbins(100), fPixId(-1)
138{
139
140 fName = name ? name : "MHGausEvents";
141 fTitle = title ? title : "Events with expected Gaussian distributions";
142
143 Clear();
144
145 SetBinsAfterStripping();
146 SetBlackoutLimit();
147 SetNDFLimit();
148 SetPickupLimit();
149 SetPowerProbabilityBins();
150 SetProbLimit();
151
152 fHGausHist.SetName("HGausHist");
153 fHGausHist.SetTitle("Histogram of Events with Gaussian Distribution");
154 // important, other ROOT will not draw the axes:
155 fHGausHist.UseCurrentStyle();
156 fHGausHist.SetDirectory(NULL);
157}
158
159
160
161// --------------------------------------------------------------------------
162//
163// Default Destructor.
164//
165// Deletes (if Pointer is not NULL):
166//
167// - fHPowerProbability
168// - fFGausFit
169// - fFExpFit
170// - fPowerSpectrum
171// - fGraphEvents
172// - fGraphPowerSpectrum
173//
174MHGausEvents::~MHGausEvents()
175{
176
177 // delete histograms
178 if (fHPowerProbability)
179 delete fHPowerProbability;
180
181 // delete fits
182 if (fFGausFit)
183 delete fFGausFit;
184 if (fFExpFit)
185 delete fFExpFit;
186
187 // delete arrays
188 if (fPowerSpectrum)
189 delete fPowerSpectrum;
190
191 // delete graphs
192 if (fGraphEvents)
193 delete fGraphEvents;
194 if (fGraphPowerSpectrum)
195 delete fGraphPowerSpectrum;
196}
197
198// --------------------------------------------------------------------------
199//
200// Default Clear(), can be overloaded.
201//
202// Sets:
203// - all other pointers to NULL
204// - all variables to 0., except fPixId to -1 and keep fEventFrequency
205// - all flags to kFALSE
206//
207// Deletes (if not NULL):
208// - all pointers
209//
210void MHGausEvents::Clear(Option_t *o)
211{
212
213 SetGausFitOK ( kFALSE );
214 SetExpFitOK ( kFALSE );
215 SetFourierSpectrumOK( kFALSE );
216 SetExcluded ( kFALSE );
217
218 fMean = 0.;
219 fSigma = 0.;
220 fMeanErr = 0.;
221 fSigmaErr = 0.;
222 fProb = 0.;
223
224 fCurrentSize = 0;
225 fPixId = -1;
226
227 if (fHPowerProbability)
228 {
229 delete fHPowerProbability;
230 fHPowerProbability = NULL;
231 }
232
233 // delete fits
234 if (fFGausFit)
235 {
236 delete fFGausFit;
237 fFGausFit = NULL;
238 }
239
240 if (fFExpFit)
241 {
242 delete fFExpFit;
243 fFExpFit = NULL;
244 }
245
246 // delete arrays
247 if (fPowerSpectrum)
248 {
249 delete fPowerSpectrum;
250 fPowerSpectrum = NULL;
251 }
252
253 // delete graphs
254 if (fGraphEvents)
255 {
256 delete fGraphEvents;
257 fGraphEvents = NULL;
258 }
259
260 if (fGraphPowerSpectrum)
261 {
262 delete fGraphPowerSpectrum;
263 fGraphPowerSpectrum = NULL;
264 }
265}
266
267// --------------------------------------------------------------------------
268//
269// Default Reset(), can be overloaded.
270//
271// Executes:
272// - Clear()
273// - fHGausHist.Reset()
274// - fEvents.Set(0)
275//
276void MHGausEvents::Reset()
277{
278
279 Clear();
280 fHGausHist.Reset();
281 fEvents.Set(0);
282
283}
284
285// --------------------------------------------------------------------------
286//
287// Default InitBins, can be overloaded.
288//
289// Executes:
290// - fHGausHist.SetBins(fNbins,fFirst,fLast)
291//
292void MHGausEvents::InitBins()
293{
294 fHGausHist.SetBins(fNbins,fFirst,fLast);
295}
296
297// --------------------------------------------------------------------------
298//
299// Executes:
300// - FillArray()
301// - FillHist()
302//
303Bool_t MHGausEvents::FillHistAndArray(const Float_t f)
304{
305
306 FillArray(f);
307 return FillHist(f);
308}
309
310// --------------------------------------------------------------------------
311//
312// Fills fHGausHist with f
313// Returns kFALSE, if overflow or underflow occurred, else kTRUE
314//
315Bool_t MHGausEvents::FillHist(const Float_t f)
316{
317 return fHGausHist.Fill(f) > -1;
318}
319
320// --------------------------------------------------------------------------
321//
322// Fill fEvents with f
323// If size of fEvents is 0, initializes it to 512
324// If size of fEvents is smaller than fCurrentSize, double the size
325// Increase fCurrentSize by 1
326//
327void MHGausEvents::FillArray(const Float_t f)
328{
329 if (fEvents.GetSize() == 0)
330 fEvents.Set(512);
331
332 if (fCurrentSize >= fEvents.GetSize())
333 fEvents.Set(fEvents.GetSize()*2);
334
335 fEvents.AddAt(f,fCurrentSize++);
336}
337
338// --------------------------------------------------------------------------
339//
340// Set Excluded bit from outside
341//
342void MHGausEvents::SetExcluded(const Bool_t b)
343{
344 b ? SETBIT(fFlags,kExcluded) : CLRBIT(fFlags,kExcluded);
345}
346
347
348const Double_t MHGausEvents::GetChiSquare() const
349{
350 return ( fFGausFit ? fFGausFit->GetChisquare() : 0.);
351}
352
353const Int_t MHGausEvents::GetNdf() const
354{
355 return ( fFGausFit ? fFGausFit->GetNDF() : 0);
356}
357
358
359const Double_t MHGausEvents::GetExpChiSquare() const
360{
361 return ( fFExpFit ? fFExpFit->GetChisquare() : 0.);
362}
363
364
365const Int_t MHGausEvents::GetExpNdf() const
366{
367 return ( fFExpFit ? fFExpFit->GetNDF() : 0);
368}
369
370
371const Double_t MHGausEvents::GetExpProb() const
372{
373 return ( fFExpFit ? fFExpFit->GetProb() : 0.);
374}
375
376
377const Double_t MHGausEvents::GetOffset() const
378{
379 return ( fFExpFit ? fFExpFit->GetParameter(0) : 0.);
380}
381
382
383const Double_t MHGausEvents::GetSlope() const
384{
385 return ( fFExpFit ? fFExpFit->GetParameter(1) : 0.);
386}
387
388
389const Bool_t MHGausEvents::IsEmpty() const
390{
391 return !(fHGausHist.GetEntries());
392}
393
394
395const Bool_t MHGausEvents::IsFourierSpectrumOK() const
396{
397 return TESTBIT(fFlags,kFourierSpectrumOK);
398}
399
400
401const Bool_t MHGausEvents::IsGausFitOK() const
402{
403 return TESTBIT(fFlags,kGausFitOK);
404}
405
406
407const Bool_t MHGausEvents::IsExpFitOK() const
408{
409 return TESTBIT(fFlags,kExpFitOK);
410}
411
412const Bool_t MHGausEvents::IsExcluded() const
413{
414 return TESTBIT(fFlags,kExcluded);
415}
416
417
418// -------------------------------------------------------------------
419//
420// The flag setters are to be used ONLY for Monte-Carlo!!
421//
422void MHGausEvents::SetGausFitOK(const Bool_t b)
423{
424 b ? SETBIT(fFlags,kGausFitOK) : CLRBIT(fFlags,kGausFitOK);
425
426}
427// -------------------------------------------------------------------
428//
429// The flag setters are to be used ONLY for Monte-Carlo!!
430//
431void MHGausEvents::SetExpFitOK(const Bool_t b)
432{
433
434 b ? SETBIT(fFlags,kExpFitOK) : CLRBIT(fFlags,kExpFitOK);
435}
436
437// -------------------------------------------------------------------
438//
439// The flag setters are to be used ONLY for Monte-Carlo!!
440//
441void MHGausEvents::SetFourierSpectrumOK(const Bool_t b)
442{
443
444 b ? SETBIT(fFlags,kFourierSpectrumOK) : CLRBIT(fFlags,kFourierSpectrumOK);
445}
446
447
448// -------------------------------------------------------------------
449//
450// Create the fourier spectrum using the class MFFT.
451// The result is projected into a histogram and fitted by an exponential
452//
453void MHGausEvents::CreateFourierSpectrum()
454{
455
456 if (fFExpFit)
457 return;
458
459 if (fEvents.GetSize() < 8)
460 {
461 *fLog << warn << "Cannot create Fourier spectrum in pixel: " << fPixId
462 << ". Number of events smaller than 8 " << endl;
463 return;
464 }
465
466
467 //
468 // The number of entries HAS to be a potence of 2,
469 // so we can only cut out from the last potence of 2 to the rest.
470 // Another possibility would be to fill everything with
471 // zeros, but that gives a low frequency peak, which we would
472 // have to cut out later again.
473 //
474 // So, we have to live with the possibility that at the end
475 // of the calibration run, something has happened without noticing
476 // it...
477 //
478
479 // This cuts only the non-used zero's, but MFFT will later cut the rest
480 MArray::StripZeros(fEvents);
481
482 if (fEvents.GetSize() < 8)
483 {
484 *fLog << warn << "Cannot create Fourier spectrum. " << endl;
485 *fLog << warn << "Number of events (after stripping of zeros) is smaller than 8 "
486 << "in pixel: " << fPixId << endl;
487 return;
488 }
489
490 MFFT fourier;
491
492 fPowerSpectrum = fourier.PowerSpectrumDensity(&fEvents);
493 fHPowerProbability = ProjectArray(*fPowerSpectrum, fPowerProbabilityBins,
494 "PowerProbability",
495 "Probability of Power occurrance");
496 fHPowerProbability->SetXTitle("P(f)");
497 fHPowerProbability->SetDirectory(NULL);
498 //
499 // First guesses for the fit (should be as close to reality as possible,
500 //
501 const Double_t xmax = fHPowerProbability->GetXaxis()->GetXmax();
502
503 fFExpFit = new TF1("FExpFit","exp([0]-[1]*x)",0.,xmax);
504
505 const Double_t slope_guess = (TMath::Log(fHPowerProbability->GetEntries())+1.)/xmax;
506 const Double_t offset_guess = slope_guess*xmax;
507
508 fFExpFit->SetParameters(offset_guess, slope_guess);
509 fFExpFit->SetParNames("Offset","Slope");
510 fFExpFit->SetParLimits(0,offset_guess/2.,2.*offset_guess);
511 fFExpFit->SetParLimits(1,slope_guess/1.5,1.5*slope_guess);
512 fFExpFit->SetRange(0.,xmax);
513
514 fHPowerProbability->Fit(fFExpFit,"RQL0");
515
516 if (GetExpProb() > fProbLimit)
517 SetExpFitOK(kTRUE);
518
519 //
520 // For the moment, this is the only check, later we can add more...
521 //
522 SetFourierSpectrumOK(IsExpFitOK());
523
524 return;
525}
526
527// -------------------------------------------------------------------
528//
529// Fit fGausHist with a Gaussian after stripping zeros from both ends
530// and rebinned to the number of bins specified in fBinsAfterStripping
531//
532// The fit results are retrieved and stored in class-own variables.
533//
534// A flag IsGausFitOK() is set according to whether the fit probability
535// is smaller or bigger than fProbLimit, whether the NDF is bigger than
536// fNDFLimit and whether results are NaNs.
537//
538Bool_t MHGausEvents::FitGaus(Option_t *option, const Double_t xmin, const Double_t xmax)
539{
540
541 if (IsGausFitOK())
542 return kTRUE;
543
544 //
545 // First, strip the zeros from the edges which contain only zeros and rebin
546 // to about fBinsAfterStripping bins.
547 //
548 // (ATTENTION: The Chisquare method is more sensitive,
549 // the _less_ bins, you have!)
550 //
551 StripZeros(&fHGausHist,fBinsAfterStripping);
552
553 //
554 // Get the fitting ranges
555 //
556 Axis_t rmin = (xmin==0.) && (xmax==0.) ? fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst()) : xmin;
557 Axis_t rmax = (xmin==0.) && (xmax==0.) ? fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast()) : xmax;
558
559 //
560 // First guesses for the fit (should be as close to reality as possible,
561 //
562 const Stat_t entries = fHGausHist.Integral("width");
563 const Double_t mu_guess = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
564 const Double_t sigma_guess = fHGausHist.GetRMS();
565 const Double_t area_guess = entries/TMath::Sqrt(TMath::TwoPi())/sigma_guess;
566
567 fFGausFit = new TF1("GausFit","gaus",rmin,rmax);
568
569 if (!fFGausFit)
570 {
571 *fLog << warn << dbginf << "WARNING: Could not create fit function for Gauss fit "
572 << "in pixel: " << fPixId << endl;
573 return kFALSE;
574 }
575
576 fFGausFit->SetParameters(area_guess,mu_guess,sigma_guess);
577 fFGausFit->SetParNames("Area","#mu","#sigma");
578 fFGausFit->SetParLimits(0,0.,area_guess*1.5);
579 fFGausFit->SetParLimits(1,rmin,rmax);
580 fFGausFit->SetParLimits(2,0.,rmax-rmin);
581 fFGausFit->SetRange(rmin,rmax);
582
583 fHGausHist.Fit(fFGausFit,option);
584
585
586 fMean = fFGausFit->GetParameter(1);
587 fSigma = fFGausFit->GetParameter(2);
588 fMeanErr = fFGausFit->GetParError(1);
589 fSigmaErr = fFGausFit->GetParError(2);
590 fProb = fFGausFit->GetProb();
591 //
592 // The fit result is accepted under condition:
593 // 1) The results are not nan's
594 // 2) The NDF is not smaller than fNDFLimit (default: fgNDFLimit)
595 // 3) The Probability is greater than fProbLimit (default: fgProbLimit)
596 //
597 if ( TMath::IsNaN(fMean)
598 || TMath::IsNaN(fMeanErr)
599 || TMath::IsNaN(fProb)
600 || TMath::IsNaN(fSigma)
601 || TMath::IsNaN(fSigmaErr)
602 || fFGausFit->GetNDF() < fNDFLimit
603 || fProb < fProbLimit )
604 return kFALSE;
605
606 SetGausFitOK(kTRUE);
607 return kTRUE;
608}
609
610// -----------------------------------------------------------------------------
611//
612// If flag IsGausFitOK() is set (histogram already successfully fitted),
613// returns kTRUE
614//
615// If both fMean and fSigma are still zero, call FitGaus()
616//
617// Repeats the Gauss fit in a smaller range, defined by:
618//
619// min = GetMean() - fBlackoutLimit * GetSigma();
620// max = GetMean() + fPickupLimit * GetSigma();
621//
622// The fit results are retrieved and stored in class-own variables.
623//
624// A flag IsGausFitOK() is set according to whether the fit probability
625// is smaller or bigger than fProbLimit, whether the NDF is bigger than
626// fNDFLimit and whether results are NaNs.
627//
628Bool_t MHGausEvents::RepeatFit(const Option_t *option)
629{
630
631 if (IsGausFitOK())
632 return kTRUE;
633
634 if ((fMean == 0.) && (fSigma == 0.))
635 return FitGaus();
636
637 //
638 // Get new fitting ranges
639 //
640 Axis_t rmin = fMean - fBlackoutLimit * fSigma;
641 Axis_t rmax = fMean + fPickupLimit * fSigma;
642
643 Axis_t hmin = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst());
644 Axis_t hmax = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast()) ;
645
646 fFGausFit->SetRange(hmin < rmin ? rmin : hmin , hmax > rmax ? rmax : hmax);
647
648 fHGausHist.Fit(fFGausFit,option);
649
650 fMean = fFGausFit->GetParameter(1);
651 fMeanErr = fFGausFit->GetParameter(2);
652 fSigma = fFGausFit->GetParError(1) ;
653 fSigmaErr = fFGausFit->GetParError(2) ;
654 fProb = fFGausFit->GetProb() ;
655
656 //
657 // The fit result is accepted under condition:
658 // 1) The results are not nan's
659 // 2) The NDF is not smaller than fNDFLimit (default: fgNDFLimit)
660 // 3) The Probability is greater than fProbLimit (default: fgProbLimit)
661 //
662 if ( TMath::IsNaN ( fMean )
663 || TMath::IsNaN ( fMeanErr )
664 || TMath::IsNaN ( fProb )
665 || TMath::IsNaN ( fSigma )
666 || TMath::IsNaN ( fSigmaErr )
667 || fFGausFit->GetNDF() < fNDFLimit
668 || fProb < fProbLimit )
669 return kFALSE;
670
671 SetGausFitOK(kTRUE);
672 return kTRUE;
673
674}
675
676// -----------------------------------------------------------------------------
677//
678// Bypasses the Gauss fit by taking mean and RMS from the histogram
679//
680// Errors are determined in the following way:
681// MeanErr = RMS / Sqrt(entries)
682// SigmaErr = RMS / (2.*Sqrt(entries) )
683//
684void MHGausEvents::BypassFit()
685{
686
687 const Stat_t entries = fHGausHist.GetEntries();
688
689 if (entries <= 0.)
690 {
691 *fLog << warn << GetDescriptor()
692 << ": Cannot bypass fit. Number of entries smaller or equal 0 in pixel: " << fPixId << endl;
693 return;
694 }
695
696 fMean = fHGausHist.GetMean();
697 fMeanErr = fHGausHist.GetRMS() / TMath::Sqrt(entries);
698 fSigma = fHGausHist.GetRMS() ;
699 fSigmaErr = fHGausHist.GetRMS() / TMath::Sqrt(entries) / 2.;
700}
701
702// -----------------------------------------------------------------------------------
703//
704// A default print
705//
706void MHGausEvents::Print(const Option_t *o) const
707{
708
709 *fLog << all << endl;
710 *fLog << all << "Results of the Gauss Fit in pixel: " << fPixId << endl;
711 *fLog << all << "Mean: " << GetMean() << endl;
712 *fLog << all << "Sigma: " << GetSigma() << endl;
713 *fLog << all << "Chisquare: " << GetChiSquare() << endl;
714 *fLog << all << "DoF: " << GetNdf() << endl;
715 *fLog << all << "Probability: " << GetProb() << endl;
716 *fLog << all << endl;
717
718}
719
720// --------------------------------------------------------------------------
721//
722// - Set fPixId to id
723//
724// Add id to names and titles of:
725// - fHGausHist
726//
727void MHGausEvents::ChangeHistId(const Int_t id)
728{
729
730 fPixId = id;
731
732 fHGausHist.SetName( Form("%s%d", fHGausHist.GetName(), id));
733 fHGausHist.SetTitle( Form("%s%d", fHGausHist.GetTitle(), id));
734
735 fName = Form("%s%d", fName.Data(), id);
736 fTitle = Form("%s%d", fTitle.Data(), id);
737
738}
739
740// --------------------------------------------------------------------------
741//
742// Re-normalize the results, has to be overloaded
743//
744void MHGausEvents::Renorm()
745{
746}
747
748// ----------------------------------------------------------------------------------
749//
750// Create a graph to display the array fEvents
751// If the variable fEventFrequency is set, the x-axis is transformed into real time.
752//
753void MHGausEvents::CreateGraphEvents()
754{
755
756 MArray::StripZeros(fEvents);
757
758 const Int_t n = fEvents.GetSize();
759
760 fGraphEvents = new TGraph(n,CreateEventXaxis(n),fEvents.GetArray());
761 fGraphEvents->SetTitle("Evolution of Events with time");
762 fGraphEvents->GetXaxis()->SetTitle((fEventFrequency) ? "Time [s]" : "Event Nr.");
763}
764
765// ----------------------------------------------------------------------------------
766//
767// Create a graph to display the array fPowerSpectrum
768// If the variable fEventFrequency is set, the x-axis is transformed into real frequency.
769//
770void MHGausEvents::CreateGraphPowerSpectrum()
771{
772
773 MArray::StripZeros(*fPowerSpectrum);
774
775 const Int_t n = fPowerSpectrum->GetSize();
776
777 fGraphPowerSpectrum = new TGraph(n,CreatePSDXaxis(n),fPowerSpectrum->GetArray());
778 fGraphPowerSpectrum->SetTitle("Power Spectrum Density");
779 fGraphPowerSpectrum->GetXaxis()->SetTitle((fEventFrequency) ? "Frequency [Hz]" : "Frequency");
780 fGraphPowerSpectrum->GetYaxis()->SetTitle("P(f)");
781}
782
783// -----------------------------------------------------------------------------
784//
785// Create the x-axis for the event graph
786//
787Float_t *MHGausEvents::CreateEventXaxis(Int_t n)
788{
789
790 Float_t *xaxis = new Float_t[n];
791
792 if (fEventFrequency)
793 for (Int_t i=0;i<n;i++)
794 xaxis[i] = (Float_t)i/fEventFrequency;
795 else
796 for (Int_t i=0;i<n;i++)
797 xaxis[i] = (Float_t)i;
798
799 return xaxis;
800
801}
802
803// -----------------------------------------------------------------------------
804//
805// Create the x-axis for the event graph
806//
807Float_t *MHGausEvents::CreatePSDXaxis(Int_t n)
808{
809
810 Float_t *xaxis = new Float_t[n];
811
812 if (fEventFrequency)
813 for (Int_t i=0;i<n;i++)
814 xaxis[i] = 0.5*(Float_t)i*fEventFrequency/n;
815 else
816 for (Int_t i=0;i<n;i++)
817 xaxis[i] = 0.5*(Float_t)i/n;
818
819 return xaxis;
820
821}
822
823// -----------------------------------------------------------------------------
824//
825// Default draw:
826//
827// The following options can be chosen:
828//
829// "EVENTS": displays a TGraph of the array fEvents
830// "FOURIER": display a TGraph of the fourier transform of fEvents
831// displays the projection of the fourier transform with the fit
832//
833void MHGausEvents::Draw(const Option_t *opt)
834{
835
836 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this,600, 900);
837
838 TString option(opt);
839 option.ToLower();
840
841 Int_t win = 1;
842
843 if (option.Contains("events"))
844 {
845 option.ReplaceAll("events","");
846 win += 1;
847 }
848 if (option.Contains("fourier"))
849 {
850 option.ReplaceAll("fourier","");
851 win += 2;
852 }
853
854 pad->SetTicks();
855 pad->SetBorderMode(0);
856 pad->Divide(1,win);
857 pad->cd(1);
858
859 if (!IsEmpty())
860 gPad->SetLogy();
861
862 fHGausHist.Draw(opt);
863
864 if (fFGausFit)
865 {
866 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
867 fFGausFit->Draw("same");
868 }
869 switch (win)
870 {
871 case 2:
872 pad->cd(2);
873 DrawEvents();
874 break;
875 case 3:
876 pad->cd(2);
877 DrawPowerSpectrum(*pad,3);
878 break;
879 case 4:
880 pad->cd(2);
881 DrawEvents();
882 pad->cd(3);
883 DrawPowerSpectrum(*pad,4);
884 break;
885 }
886}
887
888void MHGausEvents::DrawEvents()
889{
890
891 if (!fGraphEvents)
892 CreateGraphEvents();
893
894 fGraphEvents->SetBit(kCanDelete);
895 fGraphEvents->SetTitle("Events with time");
896 fGraphEvents->Draw("AL");
897
898}
899
900
901void MHGausEvents::DrawPowerSpectrum(TVirtualPad &pad, Int_t i)
902{
903
904 if (fPowerSpectrum)
905 {
906 if (!fGraphPowerSpectrum)
907 CreateGraphPowerSpectrum();
908
909 fGraphPowerSpectrum->Draw("AL");
910 fGraphPowerSpectrum->SetBit(kCanDelete);
911 }
912
913 pad.cd(i);
914
915 if (fHPowerProbability && fHPowerProbability->GetEntries() > 0)
916 {
917 gPad->SetLogy();
918 fHPowerProbability->Draw();
919 if (fFExpFit)
920 {
921 fFExpFit->SetLineColor(IsExpFitOK() ? kGreen : kRed);
922 fFExpFit->Draw("same");
923 }
924 }
925}
926
927const Double_t MHGausEvents::GetPickup() const
928{
929
930 if ((fMean == 0.) && (fSigma == 0.))
931 return -1.;
932
933 const Int_t first = fHGausHist.GetXaxis()->FindBin(fMean+fPickupLimit*fSigma);
934 const Int_t last = fHGausHist.GetXaxis()->GetLast();
935
936 if (first >= last)
937 return 0.;
938
939 return fHGausHist.Integral(first, last, "width");
940}
941
942const Double_t MHGausEvents::GetBlackout() const
943{
944
945 if ((fMean == 0.) && (fSigma == 0.))
946 return -1.;
947
948 const Int_t first = fHGausHist.GetXaxis()->GetFirst();
949 const Int_t last = fHGausHist.GetXaxis()->FindBin(fMean-fBlackoutLimit*fSigma);
950
951 if (first >= last)
952 return 0.;
953
954 return fHGausHist.Integral(first, last, "width");
955}
956
957
Note: See TracBrowser for help on using the repository browser.