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

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