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