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

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