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

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