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

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