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

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