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

Last change on this file since 3616 was 3616, checked in by gaug, 20 years ago
*** empty log message ***
File size: 20.0 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. with the function TH1F::SetBins(..) )
43// - The histogram is filled with the functions FillHist() or FillHistAndArray().
44// - The histogram can be fitted with the function FitGaus(). This involves stripping
45// of all zeros at the lower and upper end of the histogram and re-binning to
46// a new number of bins, specified in the variable fBinsAfterStripping.
47// - The fit result's probability is compared to a reference probability fProbLimit
48// The NDF is compared to fNDFLimit and a check is made whether results are NaNs.
49// Accordingly, a flag IsGausFitOK() is set.
50//
51// 2) The TArrayF fEvents:
52// ==========================
53//
54// It is created with 0 entries and not expanded unless FillArray() or FillHistAndArray()
55// are called.
56// - A first call to FillArray() or FillHistAndArray() initializes fEvents by default
57// to 512 entries.
58// - Any later call to FillArray() or FillHistAndArray() fills up the array.
59// Reaching the limit, the array is expanded by a factor 2.
60// - The array can be fourier-transformed into the array fPowerSpectrum.
61// Note that any FFT accepts only number of events which are a power of 2.
62// Thus, fEvents is cut to the next power of 2 smaller than its actual number of entries.
63// Be aware that you might lose information at the end of your analysis.
64// - Calling the function CreateFourierSpectrum() creates the array fPowerSpectrum
65// and its projection fHPowerProbability which in turn is fit to an exponential.
66// - The fit result's probability is compared to a referenc probability fProbLimit
67// and accordingly the flag IsExpFitOK() is set.
68// - The flag IsFourierSpectrumOK() is set accordingly to IsExpFitOK().
69// Later, a closer check will be installed.
70// - You can display all arrays by calls to: CreateGraphEvents() and
71// CreateGraphPowerSpectrum() and successive calls to the corresponding Getters.
72//
73//////////////////////////////////////////////////////////////////////////////
74#include "MHGausEvents.h"
75
76#include <TH1.h>
77#include <TF1.h>
78#include <TGraph.h>
79#include <TPad.h>
80#include <TVirtualPad.h>
81#include <TCanvas.h>
82#include <TStyle.h>
83
84#include "MFFT.h"
85#include "MArray.h"
86
87#include "MH.h"
88
89#include "MLog.h"
90#include "MLogManip.h"
91
92ClassImp(MHGausEvents);
93
94using namespace std;
95
96const Float_t MHGausEvents::fgProbLimit = 0.001;
97const Int_t MHGausEvents::fgNDFLimit = 2;
98const Int_t MHGausEvents::fgPowerProbabilityBins = 20;
99const Int_t MHGausEvents::fgBinsAfterStripping = 40;
100// --------------------------------------------------------------------------
101//
102// Default Constructor.
103// Sets:
104// - the default Probability Bins for fPowerProbabilityBins (fgPowerProbabilityBins)
105// - the default Probability Limit for fProbLimit (fgProbLimit)
106// - the default NDF Limit for fNDFLimit (fgNDFLimit)
107// - the default number of bins after stripping for fBinsAfterStipping (fgBinsAfterStipping)
108// - the default name of the fHGausHist ("HGausHist")
109// - the default title of the fHGausHist ("Histogram of Events with Gaussian Distribution")
110// - the default directory of the fHGausHist (NULL)
111//
112// Initializes:
113// - fEvents to 0 entries
114// - fHGausHist()
115// - all other pointers to NULL
116// - all variables to 0.
117// - all flags to kFALSE
118//
119MHGausEvents::MHGausEvents(const char *name, const char *title)
120 : fEventFrequency(0), fHPowerProbability(NULL),
121 fPowerSpectrum(NULL),
122 fGraphEvents(NULL), fGraphPowerSpectrum(NULL),
123 fHGausHist(), fEvents(0),
124 fFGausFit(NULL), fFExpFit(NULL)
125{
126
127 fName = name ? name : "MHGausEvents";
128 fTitle = title ? title : "Events with expected Gaussian distributions";
129
130 Clear();
131
132 SetPowerProbabilityBins();
133 SetProbLimit();
134 SetNDFLimit();
135 SetBinsAfterStripping();
136
137 fHGausHist.SetName("HGausHist");
138 fHGausHist.SetTitle("Histogram of Events with Gaussian Distribution");
139 // important, other ROOT will not draw the axes:
140 fHGausHist.UseCurrentStyle();
141 fHGausHist.SetDirectory(NULL);
142}
143
144
145
146// --------------------------------------------------------------------------
147//
148// Default Destructor.
149//
150// Deletes (if Pointer is not NULL):
151//
152// - fHPowerProbability
153// - fFGausFit
154// - fFExpFit
155// - fPowerSpectrum
156// - fGraphEvents
157// - fGraphPowerSpectrum
158//
159MHGausEvents::~MHGausEvents()
160{
161
162 // delete histograms
163 if (fHPowerProbability)
164 delete fHPowerProbability;
165
166 // delete fits
167 if (fFGausFit)
168 delete fFGausFit;
169 if (fFExpFit)
170 delete fFExpFit;
171
172 // delete arrays
173 if (fPowerSpectrum)
174 delete fPowerSpectrum;
175
176 // delete graphs
177 if (fGraphEvents)
178 delete fGraphEvents;
179 if (fGraphPowerSpectrum)
180 delete fGraphPowerSpectrum;
181}
182
183// --------------------------------------------------------------------------
184//
185// Sets:
186// - all other pointers to NULL
187// - all variables to 0., except fEventFrequency
188// - all flags to kFALSE
189//
190// Deletes (if not NULL):
191// - all pointers
192//
193void MHGausEvents::Clear(Option_t *o)
194{
195
196 SetGausFitOK ( kFALSE );
197 SetExpFitOK ( kFALSE );
198 SetFourierSpectrumOK( kFALSE );
199
200 fMean = 0.;
201 fSigma = 0.;
202 fMeanErr = 0.;
203 fSigmaErr = 0.;
204 fProb = 0.;
205
206 fCurrentSize = 0;
207
208 if (fHPowerProbability)
209 {
210 delete fHPowerProbability;
211 fHPowerProbability = NULL;
212 }
213
214 // delete fits
215 if (fFGausFit)
216 {
217 delete fFGausFit;
218 fFGausFit = NULL;
219 }
220
221 if (fFExpFit)
222 {
223 delete fFExpFit;
224 fFExpFit = NULL;
225 }
226
227 // delete arrays
228 if (fPowerSpectrum)
229 {
230 delete fPowerSpectrum;
231 fPowerSpectrum = NULL;
232 }
233
234 // delete graphs
235 if (fGraphEvents)
236 {
237 delete fGraphEvents;
238 fGraphEvents = NULL;
239 }
240
241 if (fGraphPowerSpectrum)
242 {
243 delete fGraphPowerSpectrum;
244 fGraphPowerSpectrum = NULL;
245 }
246}
247
248// --------------------------------------------------------------------------
249//
250// Executes:
251// - Clear()
252// - fHGausHist.Reset()
253// - fEvents.Set(0)
254//
255void MHGausEvents::Reset()
256{
257
258 Clear();
259 fHGausHist.Reset();
260 fEvents.Set(0);
261
262}
263
264// --------------------------------------------------------------------------
265//
266// Executes:
267// - FillArray()
268// - FillHist()
269//
270Bool_t MHGausEvents::FillHistAndArray(const Float_t f)
271{
272
273 FillArray(f);
274 return FillHist(f);
275}
276
277// --------------------------------------------------------------------------
278//
279// Fills fHGausHist with f
280// Returns kFALSE, if overflow or underflow occurred, else kTRUE
281//
282Bool_t MHGausEvents::FillHist(const Float_t f)
283{
284
285 if (fHGausHist.Fill(f) == -1)
286 return kFALSE;
287
288 return kTRUE;
289}
290
291// --------------------------------------------------------------------------
292//
293// Fill fEvents with f
294// If size of fEvents is 0, initializes it to 512
295// If size of fEvents is smaller than fCurrentSize, double the size
296// Increase fCurrentSize by 1
297//
298void MHGausEvents::FillArray(const Float_t f)
299{
300 if (fEvents.GetSize() == 0)
301 fEvents.Set(512);
302
303 if (fCurrentSize >= fEvents.GetSize())
304 fEvents.Set(fEvents.GetSize()*2);
305
306 fEvents.AddAt(f,fCurrentSize++);
307}
308
309
310const Double_t MHGausEvents::GetChiSquare() const
311{
312 return ( fFGausFit ? fFGausFit->GetChisquare() : 0.);
313}
314
315const Int_t MHGausEvents::GetNdf() const
316{
317 return ( fFGausFit ? fFGausFit->GetNDF() : 0);
318}
319
320
321const Double_t MHGausEvents::GetExpChiSquare() const
322{
323 return ( fFExpFit ? fFExpFit->GetChisquare() : 0.);
324}
325
326
327const Int_t MHGausEvents::GetExpNdf() const
328{
329 return ( fFExpFit ? fFExpFit->GetNDF() : 0);
330}
331
332
333const Double_t MHGausEvents::GetExpProb() const
334{
335 return ( fFExpFit ? fFExpFit->GetProb() : 0.);
336}
337
338
339const Double_t MHGausEvents::GetOffset() const
340{
341 return ( fFExpFit ? fFExpFit->GetParameter(0) : 0.);
342}
343
344
345const Double_t MHGausEvents::GetSlope() const
346{
347 return ( fFExpFit ? fFExpFit->GetParameter(1) : 0.);
348}
349
350
351const Bool_t MHGausEvents::IsEmpty() const
352{
353 return !(fHGausHist.GetEntries());
354}
355
356
357const Bool_t MHGausEvents::IsFourierSpectrumOK() const
358{
359 return TESTBIT(fFlags,kFourierSpectrumOK);
360}
361
362
363const Bool_t MHGausEvents::IsGausFitOK() const
364{
365 return TESTBIT(fFlags,kGausFitOK);
366}
367
368
369const Bool_t MHGausEvents::IsExpFitOK() const
370{
371 return TESTBIT(fFlags,kExpFitOK);
372}
373
374
375// -------------------------------------------------------------------
376//
377// The flag setters are to be used ONLY for Monte-Carlo!!
378//
379void MHGausEvents::SetGausFitOK(const Bool_t b)
380{
381 b ? SETBIT(fFlags,kGausFitOK) : CLRBIT(fFlags,kGausFitOK);
382
383}
384// -------------------------------------------------------------------
385//
386// The flag setters are to be used ONLY for Monte-Carlo!!
387//
388void MHGausEvents::SetExpFitOK(const Bool_t b)
389{
390
391 b ? SETBIT(fFlags,kExpFitOK) : CLRBIT(fFlags,kExpFitOK);
392}
393
394// -------------------------------------------------------------------
395//
396// The flag setters are to be used ONLY for Monte-Carlo!!
397//
398void MHGausEvents::SetFourierSpectrumOK(const Bool_t b)
399{
400
401 b ? SETBIT(fFlags,kFourierSpectrumOK) : CLRBIT(fFlags,kFourierSpectrumOK);
402}
403
404
405// -------------------------------------------------------------------
406//
407// Create the fourier spectrum using the class MFFT.
408// The result is projected into a histogram and fitted by an exponential
409//
410void MHGausEvents::CreateFourierSpectrum()
411{
412
413 if (fFExpFit)
414 return;
415
416 //
417 // The number of entries HAS to be a potence of 2,
418 // so we can only cut out from the last potence of 2 to the rest.
419 // Another possibility would be to fill everything with
420 // zeros, but that gives a low frequency peak, which we would
421 // have to cut out later again.
422 //
423 // So, we have to live with the possibility that at the end
424 // of the calibration run, something has happened without noticing
425 // it...
426 //
427
428 // This cuts only the non-used zero's, but MFFT will later cut the rest
429 MArray::StripZeros(fEvents);
430
431 MFFT fourier;
432
433 fPowerSpectrum = fourier.PowerSpectrumDensity(&fEvents);
434 fHPowerProbability = ProjectArray(*fPowerSpectrum, fPowerProbabilityBins,
435 "PowerProbability",
436 "Probability of Power occurrance");
437 fHPowerProbability->SetXTitle("P(f)");
438 fHPowerProbability->SetDirectory(NULL);
439 //
440 // First guesses for the fit (should be as close to reality as possible,
441 //
442 const Double_t xmax = fHPowerProbability->GetXaxis()->GetXmax();
443
444 fFExpFit = new TF1("FExpFit","exp([0]-[1]*x)",0.,xmax);
445
446 const Double_t slope_guess = (TMath::Log(fHPowerProbability->GetEntries())+1.)/xmax;
447 const Double_t offset_guess = slope_guess*xmax;
448
449 fFExpFit->SetParameters(offset_guess, slope_guess);
450 fFExpFit->SetParNames("Offset","Slope");
451 fFExpFit->SetParLimits(0,offset_guess/2.,2.*offset_guess);
452 fFExpFit->SetParLimits(1,slope_guess/1.5,1.5*slope_guess);
453 fFExpFit->SetRange(0.,xmax);
454
455 fHPowerProbability->Fit(fFExpFit,"RQL0");
456
457 if (GetExpProb() > fProbLimit)
458 SetExpFitOK(kTRUE);
459
460 //
461 // For the moment, this is the only check, later we can add more...
462 //
463 SetFourierSpectrumOK(IsExpFitOK());
464
465 return;
466}
467
468// -------------------------------------------------------------------
469//
470// Fit fGausHist with a Gaussian after stripping zeros from both ends
471// and rebinned to the number of bins specified in fBinsAfterStripping
472//
473// The fit results are retrieved and stored in class-own variables.
474//
475// A flag IsGausFitOK() is set according to whether the fit probability
476// is smaller or bigger than fProbLimit, whether the NDF is bigger than
477// fNDFLimit and whether results are NaNs.
478//
479Bool_t MHGausEvents::FitGaus(Option_t *option, const Double_t xmin, const Double_t xmax)
480{
481
482 if (IsGausFitOK())
483 return kTRUE;
484
485 //
486 // First, strip the zeros from the edges which contain only zeros and rebin
487 // to about fBinsAfterStripping bins.
488 //
489 // (ATTENTION: The Chisquare method is more sensitive,
490 // the _less_ bins, you have!)
491 //
492 StripZeros(&fHGausHist,fBinsAfterStripping);
493
494 //
495 // Get the fitting ranges
496 //
497 Axis_t rmin = (xmin==0.) && (xmax==0.) ? fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst()) : xmin;
498 Axis_t rmax = (xmin==0.) && (xmax==0.) ? fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast()) : xmax;
499
500 //
501 // First guesses for the fit (should be as close to reality as possible,
502 //
503 const Stat_t entries = fHGausHist.Integral("width");
504 const Double_t mu_guess = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
505 const Double_t sigma_guess = fHGausHist.GetRMS();
506 const Double_t area_guess = entries/TMath::Sqrt(TMath::TwoPi())/sigma_guess;
507
508 fFGausFit = new TF1("GausFit","gaus",rmin,rmax);
509
510 if (!fFGausFit)
511 {
512 *fLog << warn << dbginf << "WARNING: Could not create fit function for Gauss fit" << endl;
513 return kFALSE;
514 }
515
516 fFGausFit->SetParameters(area_guess,mu_guess,sigma_guess);
517 fFGausFit->SetParNames("Area","#mu","#sigma");
518 fFGausFit->SetParLimits(0,0.,area_guess*1.5);
519 fFGausFit->SetParLimits(1,rmin,rmax);
520 fFGausFit->SetParLimits(2,0.,rmax-rmin);
521 fFGausFit->SetRange(rmin,rmax);
522
523 fHGausHist.Fit(fFGausFit,option);
524
525
526 fMean = fFGausFit->GetParameter(1);
527 fSigma = fFGausFit->GetParameter(2);
528 fMeanErr = fFGausFit->GetParError(1);
529 fSigmaErr = fFGausFit->GetParError(2);
530 fProb = fFGausFit->GetProb();
531 //
532 // The fit result is accepted under condition:
533 // 1) The results are not nan's
534 // 2) The NDF is not smaller than fNDFLimit (5)
535 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
536 //
537 if ( TMath::IsNaN(fMean)
538 || TMath::IsNaN(fMeanErr)
539 || TMath::IsNaN(fProb)
540 || TMath::IsNaN(fSigma)
541 || TMath::IsNaN(fSigmaErr)
542 || fFGausFit->GetNDF() < fNDFLimit
543 || fProb < fProbLimit )
544 return kFALSE;
545
546 SetGausFitOK(kTRUE);
547 return kTRUE;
548}
549
550// -----------------------------------------------------------------------------------
551//
552// A default print
553//
554void MHGausEvents::Print(const Option_t *o) const
555{
556
557 *fLog << all << endl;
558 *fLog << all << "Results of the Gauss Fit: " << endl;
559 *fLog << all << "Mean: " << GetMean() << endl;
560 *fLog << all << "Sigma: " << GetSigma() << endl;
561 *fLog << all << "Chisquare: " << GetChiSquare() << endl;
562 *fLog << all << "DoF: " << GetNdf() << endl;
563 *fLog << all << "Probability: " << GetProb() << endl;
564 *fLog << all << endl;
565
566}
567
568// ----------------------------------------------------------------------------------
569//
570// Create a graph to display the array fEvents
571// If the variable fEventFrequency is set, the x-axis is transformed into real time.
572//
573void MHGausEvents::CreateGraphEvents()
574{
575
576 MArray::StripZeros(fEvents);
577
578 const Int_t n = fEvents.GetSize();
579
580 fGraphEvents = new TGraph(n,CreateEventXaxis(n),fEvents.GetArray());
581 fGraphEvents->SetTitle("Evolution of Events with time");
582 fGraphEvents->GetXaxis()->SetTitle((fEventFrequency) ? "Time [s]" : "Event Nr.");
583}
584
585// ----------------------------------------------------------------------------------
586//
587// Create a graph to display the array fPowerSpectrum
588// If the variable fEventFrequency is set, the x-axis is transformed into real frequency.
589//
590void MHGausEvents::CreateGraphPowerSpectrum()
591{
592
593 MArray::StripZeros(*fPowerSpectrum);
594
595 const Int_t n = fPowerSpectrum->GetSize();
596
597 fGraphPowerSpectrum = new TGraph(n,CreatePSDXaxis(n),fPowerSpectrum->GetArray());
598 fGraphPowerSpectrum->SetTitle("Power Spectrum Density");
599 fGraphPowerSpectrum->GetXaxis()->SetTitle((fEventFrequency) ? "Frequency [Hz]" : "Frequency");
600 fGraphPowerSpectrum->GetYaxis()->SetTitle("P(f)");
601}
602
603// -----------------------------------------------------------------------------
604//
605// Create the x-axis for the event graph
606//
607Float_t *MHGausEvents::CreateEventXaxis(Int_t n)
608{
609
610 Float_t *xaxis = new Float_t[n];
611
612 if (fEventFrequency)
613 for (Int_t i=0;i<n;i++)
614 xaxis[i] = (Float_t)i/fEventFrequency;
615 else
616 for (Int_t i=0;i<n;i++)
617 xaxis[i] = (Float_t)i;
618
619 return xaxis;
620
621}
622
623// -----------------------------------------------------------------------------
624//
625// Create the x-axis for the event graph
626//
627Float_t *MHGausEvents::CreatePSDXaxis(Int_t n)
628{
629
630 Float_t *xaxis = new Float_t[n];
631
632 if (fEventFrequency)
633 for (Int_t i=0;i<n;i++)
634 xaxis[i] = 0.5*(Float_t)i*fEventFrequency/n;
635 else
636 for (Int_t i=0;i<n;i++)
637 xaxis[i] = 0.5*(Float_t)i/n;
638
639 return xaxis;
640
641}
642
643// -----------------------------------------------------------------------------
644//
645// Default draw:
646//
647// The following options can be chosen:
648//
649// "EVENTS": displays a TGraph of the array fEvents
650// "FOURIER": display a TGraph of the fourier transform of fEvents
651// displays the projection of the fourier transform with the fit
652//
653void MHGausEvents::Draw(const Option_t *opt)
654{
655
656 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this,600, 900);
657
658 TString option(opt);
659 option.ToLower();
660
661 Int_t win = 1;
662
663 if (option.Contains("events"))
664 {
665 option.ReplaceAll("events","");
666 win += 1;
667 }
668 if (option.Contains("fourier"))
669 {
670 option.ReplaceAll("fourier","");
671 win += 2;
672 }
673
674 pad->SetTicks();
675 pad->SetBorderMode(0);
676 pad->Divide(1,win);
677 pad->cd(1);
678
679 if (!IsEmpty())
680 gPad->SetLogy();
681
682 fHGausHist.Draw(opt);
683
684 if (fFGausFit)
685 {
686 fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
687 fFGausFit->Draw("same");
688 }
689 switch (win)
690 {
691 case 2:
692 pad->cd(2);
693 DrawEvents();
694 break;
695 case 3:
696 pad->cd(2);
697 DrawPowerSpectrum(*pad,3);
698 break;
699 case 4:
700 pad->cd(2);
701 DrawEvents();
702 pad->cd(3);
703 DrawPowerSpectrum(*pad,4);
704 break;
705 }
706}
707
708void MHGausEvents::DrawEvents()
709{
710
711 if (!fGraphEvents)
712 CreateGraphEvents();
713
714 fGraphEvents->SetBit(kCanDelete);
715 fGraphEvents->SetTitle("Events with time");
716 fGraphEvents->Draw("AL");
717
718}
719
720
721void MHGausEvents::DrawPowerSpectrum(TVirtualPad &pad, Int_t i)
722{
723
724 if (fPowerSpectrum)
725 {
726 if (!fGraphPowerSpectrum)
727 CreateGraphPowerSpectrum();
728
729 fGraphPowerSpectrum->Draw("AL");
730 fGraphPowerSpectrum->SetBit(kCanDelete);
731 }
732
733 pad.cd(i);
734
735 if (fHPowerProbability && fHPowerProbability->GetEntries() > 0)
736 {
737 gPad->SetLogy();
738 fHPowerProbability->Draw();
739 if (fFExpFit)
740 {
741 fFExpFit->SetLineColor(IsExpFitOK() ? kGreen : kRed);
742 fFExpFit->Draw("same");
743 }
744 }
745}
746
747
Note: See TracBrowser for help on using the repository browser.