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

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