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

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