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

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