source: trunk/MagicSoft/Mars/mhcalib/MHGausEvents.cc@ 4959

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