source: tags/Mars-V1.1/mhcalib/MHGausEvents.cc

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