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

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