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

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