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 follow a Gaussian distribution
|
---|
30 | // with time, (e.g. calibration events, observables containing white noise, etc.)
|
---|
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. with the function TH1F::SetBins(..) )
|
---|
43 | // - The histogram is filled with the functions FillHist or FillHistAndArray.
|
---|
44 | // - The histogram can be fitted with the function FitGaus(). This involves stripping
|
---|
45 | // of all zeros at the lower and upper end of the histogram and re-binning to
|
---|
46 | // a new number of bins, specified in the variable fBinsAfterStripping.
|
---|
47 | // - The fit result's probability is compared to a reference probability fProbLimit
|
---|
48 | // The NDF is compared to fNDFLimit and a check is made whether results are NaNs.
|
---|
49 | // Accordingly the flag GausFitOK is set.
|
---|
50 | //
|
---|
51 | // 2) The TArrayF fEvents:
|
---|
52 | // ==========================
|
---|
53 | //
|
---|
54 | // It is created with 0 entries and not expanded unless FillArray or FillHistAndArray is called.
|
---|
55 | // - A first call to FillArray or FillHistAndArray initializes fEvents by default to 512 entries.
|
---|
56 | // - Any later call to FillArray or FillHistAndArray fills up the array. Reaching the limit,
|
---|
57 | // the array is expanded by a factor 2.
|
---|
58 | // - The array can be fourier-transformed into the array fPowerSpectrum. Note that any FFT accepts
|
---|
59 | // only number of events which are a power of 2. Thus, fEvents is cut to the next power of 2 smaller
|
---|
60 | // than its actual number of entries. You might lose information at the end of your analysis.
|
---|
61 | // - Calling the function CreateFourierTransform creates the array fPowerSpectrum and its projection
|
---|
62 | // fHPowerProbability which in turn is fit to an exponential.
|
---|
63 | // - The fit result's probability is compared to a referenc probability fProbLimit
|
---|
64 | // and accordingly the flag ExpFitOK is set.
|
---|
65 | // - The flag FourierSpectrumOK is set accordingly to ExpFitOK. Later, a closer check will be
|
---|
66 | // installed.
|
---|
67 | // - You can display all arrays by calls to: CreateGraphEvents() and CreateGraphPowerSpectrum() and
|
---|
68 | // successive calls to the corresponding Getters.
|
---|
69 | //
|
---|
70 | //////////////////////////////////////////////////////////////////////////////
|
---|
71 | #include "MHGausEvents.h"
|
---|
72 |
|
---|
73 | #include <TH1.h>
|
---|
74 | #include <TF1.h>
|
---|
75 | #include <TGraph.h>
|
---|
76 | #include <TPad.h>
|
---|
77 | #include <TVirtualPad.h>
|
---|
78 | #include <TCanvas.h>
|
---|
79 | #include <TStyle.h>
|
---|
80 |
|
---|
81 | #include "MFFT.h"
|
---|
82 | #include "MArray.h"
|
---|
83 |
|
---|
84 | #include "MH.h"
|
---|
85 |
|
---|
86 | #include "MLog.h"
|
---|
87 | #include "MLogManip.h"
|
---|
88 |
|
---|
89 | ClassImp(MHGausEvents);
|
---|
90 |
|
---|
91 | using namespace std;
|
---|
92 |
|
---|
93 | const Float_t MHGausEvents::fgProbLimit = 0.005;
|
---|
94 | const Int_t MHGausEvents::fgNDFLimit = 2;
|
---|
95 | const Int_t MHGausEvents::fgPowerProbabilityBins = 20;
|
---|
96 | const Int_t MHGausEvents::fgBinsAfterStripping = 40;
|
---|
97 | // --------------------------------------------------------------------------
|
---|
98 | //
|
---|
99 | // Default Constructor.
|
---|
100 | //
|
---|
101 | MHGausEvents::MHGausEvents(const char *name, const char *title)
|
---|
102 | : fHPowerProbability(NULL),
|
---|
103 | fPowerSpectrum(NULL),
|
---|
104 | fGraphEvents(NULL), fGraphPowerSpectrum(NULL),
|
---|
105 | fHGausHist(), fEvents(0),
|
---|
106 | fFGausFit(NULL), fFExpFit(NULL)
|
---|
107 | {
|
---|
108 |
|
---|
109 | fName = name ? name : "MHGausEvents";
|
---|
110 | fTitle = title ? title : "Events with expected Gaussian distributions";
|
---|
111 |
|
---|
112 | Clear();
|
---|
113 |
|
---|
114 | SetPowerProbabilityBins();
|
---|
115 | SetEventFrequency();
|
---|
116 | SetProbLimit();
|
---|
117 | SetNDFLimit();
|
---|
118 | SetBinsAfterStripping();
|
---|
119 |
|
---|
120 | fHGausHist.SetName("HGausHist");
|
---|
121 | fHGausHist.SetTitle("Histogram of Events with Gaussian Distribution");
|
---|
122 | // important, other ROOT will not draw the axes:
|
---|
123 | fHGausHist.UseCurrentStyle();
|
---|
124 | fHGausHist.SetDirectory(NULL);
|
---|
125 | }
|
---|
126 |
|
---|
127 |
|
---|
128 |
|
---|
129 |
|
---|
130 | MHGausEvents::~MHGausEvents()
|
---|
131 | {
|
---|
132 |
|
---|
133 | // delete histograms
|
---|
134 | if (fHPowerProbability)
|
---|
135 | delete fHPowerProbability;
|
---|
136 |
|
---|
137 | // delete fits
|
---|
138 | if (fFGausFit)
|
---|
139 | delete fFGausFit;
|
---|
140 | if (fFExpFit)
|
---|
141 | delete fFExpFit;
|
---|
142 |
|
---|
143 | // delete arrays
|
---|
144 | if (fPowerSpectrum)
|
---|
145 | delete fPowerSpectrum;
|
---|
146 |
|
---|
147 | // delete graphs
|
---|
148 | if (fGraphEvents)
|
---|
149 | delete fGraphEvents;
|
---|
150 | if (fGraphPowerSpectrum)
|
---|
151 | delete fGraphPowerSpectrum;
|
---|
152 | }
|
---|
153 |
|
---|
154 |
|
---|
155 | void MHGausEvents::Clear(Option_t *o)
|
---|
156 | {
|
---|
157 |
|
---|
158 | SetGausFitOK ( kFALSE );
|
---|
159 | SetExpFitOK ( kFALSE );
|
---|
160 | SetFourierSpectrumOK( kFALSE );
|
---|
161 |
|
---|
162 | fMean = 0.;
|
---|
163 | fSigma = 0.;
|
---|
164 | fMeanErr = 0.;
|
---|
165 | fSigmaErr = 0.;
|
---|
166 | fProb = 0.;
|
---|
167 |
|
---|
168 | fCurrentSize = 0;
|
---|
169 |
|
---|
170 | if (fHPowerProbability)
|
---|
171 | delete fHPowerProbability;
|
---|
172 |
|
---|
173 | // delete fits
|
---|
174 | if (fFGausFit)
|
---|
175 | delete fFGausFit;
|
---|
176 | if (fFExpFit)
|
---|
177 | delete fFExpFit;
|
---|
178 |
|
---|
179 | // delete arrays
|
---|
180 | if (fPowerSpectrum)
|
---|
181 | delete fPowerSpectrum;
|
---|
182 |
|
---|
183 | // delete graphs
|
---|
184 | if (fGraphEvents)
|
---|
185 | delete fGraphEvents;
|
---|
186 | if (fGraphPowerSpectrum)
|
---|
187 | delete fGraphPowerSpectrum;
|
---|
188 | }
|
---|
189 |
|
---|
190 |
|
---|
191 | void MHGausEvents::Reset()
|
---|
192 | {
|
---|
193 |
|
---|
194 | Clear();
|
---|
195 | fHGausHist.Reset();
|
---|
196 | fEvents.Set(0);
|
---|
197 |
|
---|
198 | }
|
---|
199 |
|
---|
200 |
|
---|
201 | Bool_t MHGausEvents::FillHistAndArray(const Float_t f)
|
---|
202 | {
|
---|
203 |
|
---|
204 | FillArray(f);
|
---|
205 | return FillHist(f);
|
---|
206 | }
|
---|
207 |
|
---|
208 | Bool_t MHGausEvents::FillHist(const Float_t f)
|
---|
209 | {
|
---|
210 |
|
---|
211 | if (fHGausHist.Fill(f) == -1)
|
---|
212 | return kFALSE;
|
---|
213 |
|
---|
214 | return kTRUE;
|
---|
215 | }
|
---|
216 |
|
---|
217 | void MHGausEvents::FillArray(const Float_t f)
|
---|
218 | {
|
---|
219 | if (fEvents.GetSize() == 0)
|
---|
220 | fEvents.Set(512);
|
---|
221 |
|
---|
222 | if (fCurrentSize >= fEvents.GetSize())
|
---|
223 | fEvents.Set(fEvents.GetSize()*2);
|
---|
224 |
|
---|
225 | fEvents.AddAt(f,fCurrentSize++);
|
---|
226 | }
|
---|
227 |
|
---|
228 |
|
---|
229 | const Double_t MHGausEvents::GetChiSquare() const
|
---|
230 | {
|
---|
231 | return ( fFGausFit ? fFGausFit->GetChisquare() : 0.);
|
---|
232 | }
|
---|
233 |
|
---|
234 | const Int_t MHGausEvents::GetNdf() const
|
---|
235 | {
|
---|
236 | return ( fFGausFit ? fFGausFit->GetNDF() : 0);
|
---|
237 | }
|
---|
238 |
|
---|
239 |
|
---|
240 | const Double_t MHGausEvents::GetExpChiSquare() const
|
---|
241 | {
|
---|
242 | return ( fFExpFit ? fFExpFit->GetChisquare() : 0.);
|
---|
243 | }
|
---|
244 |
|
---|
245 |
|
---|
246 | const Int_t MHGausEvents::GetExpNdf() const
|
---|
247 | {
|
---|
248 | return ( fFExpFit ? fFExpFit->GetNDF() : 0);
|
---|
249 | }
|
---|
250 |
|
---|
251 |
|
---|
252 | const Double_t MHGausEvents::GetExpProb() const
|
---|
253 | {
|
---|
254 | return ( fFExpFit ? fFExpFit->GetProb() : 0.);
|
---|
255 | }
|
---|
256 |
|
---|
257 |
|
---|
258 | const Double_t MHGausEvents::GetOffset() const
|
---|
259 | {
|
---|
260 | return ( fFExpFit ? fFExpFit->GetParameter(0) : 0.);
|
---|
261 | }
|
---|
262 |
|
---|
263 |
|
---|
264 | const Double_t MHGausEvents::GetSlope() const
|
---|
265 | {
|
---|
266 | return ( fFExpFit ? fFExpFit->GetParameter(1) : 0.);
|
---|
267 | }
|
---|
268 |
|
---|
269 |
|
---|
270 | const Bool_t MHGausEvents::IsEmpty() const
|
---|
271 | {
|
---|
272 | return !(fHGausHist.GetEntries());
|
---|
273 | }
|
---|
274 |
|
---|
275 |
|
---|
276 | const Bool_t MHGausEvents::IsFourierSpectrumOK() const
|
---|
277 | {
|
---|
278 | return TESTBIT(fFlags,kFourierSpectrumOK);
|
---|
279 | }
|
---|
280 |
|
---|
281 |
|
---|
282 | const Bool_t MHGausEvents::IsGausFitOK() const
|
---|
283 | {
|
---|
284 | return TESTBIT(fFlags,kGausFitOK);
|
---|
285 | }
|
---|
286 |
|
---|
287 |
|
---|
288 | const Bool_t MHGausEvents::IsExpFitOK() const
|
---|
289 | {
|
---|
290 | return TESTBIT(fFlags,kExpFitOK);
|
---|
291 | }
|
---|
292 |
|
---|
293 |
|
---|
294 | // -------------------------------------------------------------------
|
---|
295 | //
|
---|
296 | // The flag setters are to be used ONLY for Monte-Carlo!!
|
---|
297 | //
|
---|
298 | void MHGausEvents::SetGausFitOK(const Bool_t b)
|
---|
299 | {
|
---|
300 | b ? SETBIT(fFlags,kGausFitOK) : CLRBIT(fFlags,kGausFitOK);
|
---|
301 |
|
---|
302 | }
|
---|
303 | // -------------------------------------------------------------------
|
---|
304 | //
|
---|
305 | // The flag setters are to be used ONLY for Monte-Carlo!!
|
---|
306 | //
|
---|
307 | void MHGausEvents::SetExpFitOK(const Bool_t b)
|
---|
308 | {
|
---|
309 |
|
---|
310 | b ? SETBIT(fFlags,kExpFitOK) : CLRBIT(fFlags,kExpFitOK);
|
---|
311 | }
|
---|
312 |
|
---|
313 | // -------------------------------------------------------------------
|
---|
314 | //
|
---|
315 | // The flag setters are to be used ONLY for Monte-Carlo!!
|
---|
316 | //
|
---|
317 | void MHGausEvents::SetFourierSpectrumOK(const Bool_t b)
|
---|
318 | {
|
---|
319 |
|
---|
320 | b ? SETBIT(fFlags,kFourierSpectrumOK) : CLRBIT(fFlags,kFourierSpectrumOK);
|
---|
321 | }
|
---|
322 |
|
---|
323 |
|
---|
324 | // -------------------------------------------------------------------
|
---|
325 | //
|
---|
326 | // Create the fourier spectrum using the class MFFT.
|
---|
327 | // The result is projected into a histogram and fitted by an exponential
|
---|
328 | //
|
---|
329 | void MHGausEvents::CreateFourierSpectrum()
|
---|
330 | {
|
---|
331 |
|
---|
332 | if (fFExpFit)
|
---|
333 | return;
|
---|
334 |
|
---|
335 | //
|
---|
336 | // The number of entries HAS to be a potence of 2,
|
---|
337 | // so we can only cut out from the last potence of 2 to the rest.
|
---|
338 | // Another possibility would be to fill everything with
|
---|
339 | // zeros, but that gives a low frequency peak, which we would
|
---|
340 | // have to cut out later again.
|
---|
341 | //
|
---|
342 | // So, we have to live with the possibility that at the end
|
---|
343 | // of the calibration run, something has happened without noticing
|
---|
344 | // it...
|
---|
345 | //
|
---|
346 |
|
---|
347 | // This cuts only the non-used zero's, but MFFT will later cut the rest
|
---|
348 | MArray::StripZeros(fEvents);
|
---|
349 |
|
---|
350 | MFFT fourier;
|
---|
351 |
|
---|
352 | fPowerSpectrum = fourier.PowerSpectrumDensity(&fEvents);
|
---|
353 | fHPowerProbability = ProjectArray(*fPowerSpectrum, fPowerProbabilityBins,
|
---|
354 | "PowerProbability",
|
---|
355 | "Probability of Power occurrance");
|
---|
356 | fHPowerProbability->SetXTitle("P(f)");
|
---|
357 | fHPowerProbability->SetDirectory(NULL);
|
---|
358 | //
|
---|
359 | // First guesses for the fit (should be as close to reality as possible,
|
---|
360 | //
|
---|
361 | const Double_t xmax = fHPowerProbability->GetXaxis()->GetXmax();
|
---|
362 |
|
---|
363 | fFExpFit = new TF1("FExpFit","exp([0]-[1]*x)",0.,xmax);
|
---|
364 |
|
---|
365 | const Double_t slope_guess = (TMath::Log(fHPowerProbability->GetEntries())+1.)/xmax;
|
---|
366 | const Double_t offset_guess = slope_guess*xmax;
|
---|
367 |
|
---|
368 | fFExpFit->SetParameters(offset_guess, slope_guess);
|
---|
369 | fFExpFit->SetParNames("Offset","Slope");
|
---|
370 | fFExpFit->SetParLimits(0,offset_guess/2.,2.*offset_guess);
|
---|
371 | fFExpFit->SetParLimits(1,slope_guess/1.5,1.5*slope_guess);
|
---|
372 | fFExpFit->SetRange(0.,xmax);
|
---|
373 |
|
---|
374 | fHPowerProbability->Fit(fFExpFit,"RQL0");
|
---|
375 |
|
---|
376 | if (GetExpProb() > fProbLimit)
|
---|
377 | SetExpFitOK(kTRUE);
|
---|
378 |
|
---|
379 | //
|
---|
380 | // For the moment, this is the only check, later we can add more...
|
---|
381 | //
|
---|
382 | SetFourierSpectrumOK(IsExpFitOK());
|
---|
383 |
|
---|
384 | return;
|
---|
385 | }
|
---|
386 |
|
---|
387 | // -------------------------------------------------------------------
|
---|
388 | //
|
---|
389 | // Fit fGausHist with a Gaussian after stripping zeros from both ends
|
---|
390 | // and rebinned to the number of bins specified in fBinsAfterStripping
|
---|
391 | //
|
---|
392 | // The fit results are retrieved and stored in class-own variables.
|
---|
393 | //
|
---|
394 | // A flag IsGausFitOK() is set according to whether the fit probability
|
---|
395 | // is smaller or bigger than fProbLimit, whether the NDF is bigger than
|
---|
396 | // fNDFLimit and whether results are NaNs.
|
---|
397 | //
|
---|
398 | Bool_t MHGausEvents::FitGaus(Option_t *option)
|
---|
399 | {
|
---|
400 |
|
---|
401 | if (IsGausFitOK())
|
---|
402 | return kTRUE;
|
---|
403 |
|
---|
404 | //
|
---|
405 | // First, strip the zeros from the edges which contain only zeros and rebin
|
---|
406 | // to about fBinsAfterStripping bins.
|
---|
407 | //
|
---|
408 | // (ATTENTION: The Chisquare method is more sensitive,
|
---|
409 | // the _less_ bins, you have!)
|
---|
410 | //
|
---|
411 | StripZeros(&fHGausHist,fBinsAfterStripping);
|
---|
412 |
|
---|
413 | //
|
---|
414 | // Get the fitting ranges
|
---|
415 | //
|
---|
416 | Axis_t rmin = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst());
|
---|
417 | Axis_t rmax = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast());
|
---|
418 |
|
---|
419 | //
|
---|
420 | // First guesses for the fit (should be as close to reality as possible,
|
---|
421 | //
|
---|
422 | const Stat_t entries = fHGausHist.Integral("width");
|
---|
423 | const Double_t mu_guess = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
|
---|
424 | const Double_t sigma_guess = fHGausHist.GetRMS();
|
---|
425 | const Double_t area_guess = entries/TMath::Sqrt(TMath::TwoPi())/sigma_guess;
|
---|
426 |
|
---|
427 | fFGausFit = new TF1("GausFit","gaus",rmin,rmax);
|
---|
428 |
|
---|
429 | if (!fFGausFit)
|
---|
430 | {
|
---|
431 | *fLog << warn << dbginf << "WARNING: Could not create fit function for Gauss fit" << endl;
|
---|
432 | return kFALSE;
|
---|
433 | }
|
---|
434 |
|
---|
435 | fFGausFit->SetParameters(area_guess,mu_guess,sigma_guess);
|
---|
436 | fFGausFit->SetParNames("Area","#mu","#sigma");
|
---|
437 | fFGausFit->SetParLimits(0,0.,area_guess*1.5);
|
---|
438 | fFGausFit->SetParLimits(1,rmin,rmax);
|
---|
439 | fFGausFit->SetParLimits(2,0.,rmax-rmin);
|
---|
440 | fFGausFit->SetRange(rmin,rmax);
|
---|
441 |
|
---|
442 | fHGausHist.Fit(fFGausFit,option);
|
---|
443 |
|
---|
444 |
|
---|
445 | fMean = fFGausFit->GetParameter(1);
|
---|
446 | fSigma = fFGausFit->GetParameter(2);
|
---|
447 | fMeanErr = fFGausFit->GetParError(1);
|
---|
448 | fSigmaErr = fFGausFit->GetParError(2);
|
---|
449 | fProb = fFGausFit->GetProb();
|
---|
450 | //
|
---|
451 | // The fit result is accepted under condition:
|
---|
452 | // 1) The results are not nan's
|
---|
453 | // 2) The NDF is not smaller than fNDFLimit (5)
|
---|
454 | // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
|
---|
455 | //
|
---|
456 | if ( TMath::IsNaN(fMean)
|
---|
457 | || TMath::IsNaN(fMeanErr)
|
---|
458 | || TMath::IsNaN(fProb)
|
---|
459 | || TMath::IsNaN(fSigma)
|
---|
460 | || TMath::IsNaN(fSigmaErr)
|
---|
461 | || fFGausFit->GetNDF() < fNDFLimit
|
---|
462 | || fProb < fProbLimit )
|
---|
463 | return kFALSE;
|
---|
464 |
|
---|
465 | SetGausFitOK(kTRUE);
|
---|
466 | return kTRUE;
|
---|
467 | }
|
---|
468 |
|
---|
469 | // -----------------------------------------------------------------------------------
|
---|
470 | //
|
---|
471 | // A default print
|
---|
472 | //
|
---|
473 | void MHGausEvents::Print(const Option_t *o) const
|
---|
474 | {
|
---|
475 |
|
---|
476 | *fLog << all << endl;
|
---|
477 | *fLog << all << "Results of the Gauss Fit: " << endl;
|
---|
478 | *fLog << all << "Mean: " << GetMean() << endl;
|
---|
479 | *fLog << all << "Sigma: " << GetSigma() << endl;
|
---|
480 | *fLog << all << "Chisquare: " << GetChiSquare() << endl;
|
---|
481 | *fLog << all << "DoF: " << GetNdf() << endl;
|
---|
482 | *fLog << all << "Probability: " << GetProb() << endl;
|
---|
483 | *fLog << all << endl;
|
---|
484 |
|
---|
485 | }
|
---|
486 |
|
---|
487 | // ----------------------------------------------------------------------------------
|
---|
488 | //
|
---|
489 | // Create a graph to display the array fEvents
|
---|
490 | // If the variable fEventFrequency is set, the x-axis is transformed into real time.
|
---|
491 | //
|
---|
492 | void MHGausEvents::CreateGraphEvents()
|
---|
493 | {
|
---|
494 |
|
---|
495 | MArray::StripZeros(fEvents);
|
---|
496 |
|
---|
497 | const Int_t n = fEvents.GetSize();
|
---|
498 |
|
---|
499 | fGraphEvents = new TGraph(n,CreateXaxis(n),fEvents.GetArray());
|
---|
500 | fGraphEvents->SetTitle("Evolution of Events with time");
|
---|
501 | fGraphEvents->GetXaxis()->SetTitle((fEventFrequency) ? "Time [s]" : "Event Nr.");
|
---|
502 | }
|
---|
503 |
|
---|
504 | // ----------------------------------------------------------------------------------
|
---|
505 | //
|
---|
506 | // Create a graph to display the array fPowerSpectrum
|
---|
507 | // If the variable fEventFrequency is set, the x-axis is transformed into real frequency.
|
---|
508 | //
|
---|
509 | void MHGausEvents::CreateGraphPowerSpectrum()
|
---|
510 | {
|
---|
511 |
|
---|
512 | MArray::StripZeros(*fPowerSpectrum);
|
---|
513 |
|
---|
514 | const Int_t n = fPowerSpectrum->GetSize();
|
---|
515 |
|
---|
516 | fGraphPowerSpectrum = new TGraph(n,CreateXaxis(n),fPowerSpectrum->GetArray());
|
---|
517 | fGraphPowerSpectrum->SetTitle("Power Spectrum Density");
|
---|
518 | fGraphPowerSpectrum->GetXaxis()->SetTitle((fEventFrequency) ? "Frequency [Hz]" : "Frequency");
|
---|
519 | fGraphPowerSpectrum->GetYaxis()->SetTitle("P(f)");
|
---|
520 | }
|
---|
521 |
|
---|
522 | // -----------------------------------------------------------------------------
|
---|
523 | //
|
---|
524 | // Create the x-axis for the graph
|
---|
525 | //
|
---|
526 | Float_t *MHGausEvents::CreateXaxis(Int_t n)
|
---|
527 | {
|
---|
528 |
|
---|
529 | Float_t *xaxis = new Float_t[n];
|
---|
530 |
|
---|
531 | if (fEventFrequency)
|
---|
532 | for (Int_t i=0;i<n;i++)
|
---|
533 | xaxis[i] = (Float_t)i/fEventFrequency;
|
---|
534 | else
|
---|
535 | for (Int_t i=0;i<n;i++)
|
---|
536 | xaxis[i] = (Float_t)i;
|
---|
537 |
|
---|
538 | return xaxis;
|
---|
539 |
|
---|
540 | }
|
---|
541 |
|
---|
542 | // -----------------------------------------------------------------------------
|
---|
543 | //
|
---|
544 | // Default draw:
|
---|
545 | //
|
---|
546 | // The following options can be chosen:
|
---|
547 | //
|
---|
548 | // "EVENTS": displays a TGraph of the array fEvents
|
---|
549 | // "FOURIER": display a TGraph of the fourier transform of fEvents
|
---|
550 | // displays the projection of the fourier transform with the fit
|
---|
551 | //
|
---|
552 | void MHGausEvents::Draw(const Option_t *opt)
|
---|
553 | {
|
---|
554 |
|
---|
555 | TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this,600, 900);
|
---|
556 |
|
---|
557 | TString option(opt);
|
---|
558 | option.ToLower();
|
---|
559 |
|
---|
560 | Int_t win = 1;
|
---|
561 |
|
---|
562 | if (option.Contains("events"))
|
---|
563 | {
|
---|
564 | option.ReplaceAll("events","");
|
---|
565 | win += 1;
|
---|
566 | }
|
---|
567 | if (option.Contains("fourier"))
|
---|
568 | {
|
---|
569 | option.ReplaceAll("fourier","");
|
---|
570 | win += 2;
|
---|
571 | }
|
---|
572 |
|
---|
573 | pad->SetTicks();
|
---|
574 | pad->SetBorderMode(0);
|
---|
575 | pad->Divide(1,win);
|
---|
576 | pad->cd(1);
|
---|
577 |
|
---|
578 | if (!IsEmpty())
|
---|
579 | gPad->SetLogy();
|
---|
580 |
|
---|
581 | fHGausHist.Draw(opt);
|
---|
582 |
|
---|
583 | if (fFGausFit)
|
---|
584 | {
|
---|
585 | fFGausFit->SetLineColor(IsGausFitOK() ? kGreen : kRed);
|
---|
586 | fFGausFit->Draw("same");
|
---|
587 | }
|
---|
588 | switch (win)
|
---|
589 | {
|
---|
590 | case 2:
|
---|
591 | pad->cd(2);
|
---|
592 | DrawEvents();
|
---|
593 | break;
|
---|
594 | case 3:
|
---|
595 | pad->cd(2);
|
---|
596 | DrawPowerSpectrum(*pad,3);
|
---|
597 | break;
|
---|
598 | case 4:
|
---|
599 | pad->cd(2);
|
---|
600 | DrawEvents();
|
---|
601 | pad->cd(3);
|
---|
602 | DrawPowerSpectrum(*pad,4);
|
---|
603 | break;
|
---|
604 | }
|
---|
605 | }
|
---|
606 |
|
---|
607 | void MHGausEvents::DrawEvents()
|
---|
608 | {
|
---|
609 |
|
---|
610 | if (!fGraphEvents)
|
---|
611 | CreateGraphEvents();
|
---|
612 |
|
---|
613 | fGraphEvents->SetBit(kCanDelete);
|
---|
614 | fGraphEvents->SetTitle("Events with time");
|
---|
615 | fGraphEvents->Draw("AL");
|
---|
616 |
|
---|
617 | }
|
---|
618 |
|
---|
619 |
|
---|
620 | void MHGausEvents::DrawPowerSpectrum(TVirtualPad &pad, Int_t i)
|
---|
621 | {
|
---|
622 |
|
---|
623 | if (fPowerSpectrum)
|
---|
624 | {
|
---|
625 | if (!fGraphPowerSpectrum)
|
---|
626 | CreateGraphPowerSpectrum();
|
---|
627 |
|
---|
628 | fGraphPowerSpectrum->Draw("AL");
|
---|
629 | fGraphPowerSpectrum->SetBit(kCanDelete);
|
---|
630 | }
|
---|
631 |
|
---|
632 | pad.cd(i);
|
---|
633 |
|
---|
634 | if (fHPowerProbability && fHPowerProbability->GetEntries() > 0)
|
---|
635 | {
|
---|
636 | gPad->SetLogy();
|
---|
637 | fHPowerProbability->Draw();
|
---|
638 | if (fFExpFit)
|
---|
639 | {
|
---|
640 | fFExpFit->SetLineColor(IsExpFitOK() ? kGreen : kRed);
|
---|
641 | fFExpFit->Draw("same");
|
---|
642 | }
|
---|
643 | }
|
---|
644 | }
|
---|
645 |
|
---|
646 |
|
---|