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

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