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

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