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

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