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

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