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

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