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

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