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): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2005
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MHAlpha
|
---|
28 | //
|
---|
29 | // Create a single Alpha-Plot. The alpha-plot is fitted online. You can
|
---|
30 | // check the result when it is filles in the MStatusDisplay
|
---|
31 | // For more information see MHFalseSource::FitSignificance
|
---|
32 | //
|
---|
33 | // For convinience (fit) the output significance is stored in a
|
---|
34 | // container in the parlisrt
|
---|
35 | //
|
---|
36 | // PRELIMINARY!
|
---|
37 | //
|
---|
38 | //
|
---|
39 | // ToDo:
|
---|
40 | // =====
|
---|
41 | //
|
---|
42 | // - Replace time-binning by histogram (which is enhanced bin-wise)
|
---|
43 | //
|
---|
44 | //
|
---|
45 | //////////////////////////////////////////////////////////////////////////////
|
---|
46 | #include "MHAlpha.h"
|
---|
47 |
|
---|
48 | #include <TStyle.h>
|
---|
49 | #include <TLatex.h>
|
---|
50 | #include <TCanvas.h>
|
---|
51 | #include <TPaveStats.h>
|
---|
52 |
|
---|
53 | #include "MSrcPosCam.h"
|
---|
54 | #include "MHillas.h"
|
---|
55 | #include "MHillasSrc.h"
|
---|
56 | #include "MTime.h"
|
---|
57 | #include "MObservatory.h"
|
---|
58 | #include "MPointingPos.h"
|
---|
59 | #include "MAstroSky2Local.h"
|
---|
60 | #include "MStatusDisplay.h"
|
---|
61 | #include "MParameters.h"
|
---|
62 | #include "MHMatrix.h"
|
---|
63 |
|
---|
64 | #include "MString.h"
|
---|
65 | #include "MBinning.h"
|
---|
66 | #include "MParList.h"
|
---|
67 |
|
---|
68 | #include "MLog.h"
|
---|
69 | #include "MLogManip.h"
|
---|
70 |
|
---|
71 | ClassImp(MHAlpha);
|
---|
72 |
|
---|
73 | using namespace std;
|
---|
74 |
|
---|
75 | // --------------------------------------------------------------------------
|
---|
76 | //
|
---|
77 | // Default Constructor
|
---|
78 | //
|
---|
79 | MHAlpha::MHAlpha(const char *name, const char *title)
|
---|
80 | : fNameParameter("MHillasSrc"), fParameter(0),
|
---|
81 | fOffData(0), fResult(0), fSigma(0), fEnergy(0), fBin(0),
|
---|
82 | fPointPos(0), fTimeEffOn(0), fTime(0), fNumTimeBins(10),
|
---|
83 | fHillas(0), fMatrix(0), fSkipHistTime(kFALSE), fSkipHistTheta(kFALSE),
|
---|
84 | fSkipHistEnergy(kFALSE), fForceUsingSize(kFALSE)
|
---|
85 | {
|
---|
86 | //
|
---|
87 | // set the name and title of this object
|
---|
88 | //
|
---|
89 | fName = name ? name : "MHAlpha";
|
---|
90 | fTitle = title ? title : "Alpha plot";
|
---|
91 |
|
---|
92 | fHist.SetName("Alpha");
|
---|
93 | fHist.SetTitle("Alpha");
|
---|
94 | fHist.SetXTitle("\\Theta [deg]");
|
---|
95 | //fHist.SetYTitle("E_{est} [GeV]");
|
---|
96 | fHist.SetZTitle("|\\alpha| [\\circ]");
|
---|
97 | fHist.SetDirectory(NULL);
|
---|
98 | fHist.UseCurrentStyle();
|
---|
99 |
|
---|
100 | // Main histogram
|
---|
101 | fHistTime.SetName("Alpha");
|
---|
102 | fHistTime.SetXTitle("|\\alpha| [\\circ]");
|
---|
103 | fHistTime.SetYTitle("Counts");
|
---|
104 | fHistTime.UseCurrentStyle();
|
---|
105 | fHistTime.SetDirectory(NULL);
|
---|
106 |
|
---|
107 | fHEnergy.SetName("Excess");
|
---|
108 | //fHEnergy.SetTitle(" N_{exc} vs. E_{est} ");
|
---|
109 | //fHEnergy.SetXTitle("E_{est} [GeV]");
|
---|
110 | fHEnergy.SetYTitle("N_{exc}");
|
---|
111 | fHEnergy.SetDirectory(NULL);
|
---|
112 | fHEnergy.UseCurrentStyle();
|
---|
113 |
|
---|
114 | fHTheta.SetName("ExcessTheta");
|
---|
115 | fHTheta.SetTitle(" N_{exc} vs. \\Theta ");
|
---|
116 | fHTheta.SetXTitle("\\Theta [\\circ]");
|
---|
117 | fHTheta.SetYTitle("N_{exc}");
|
---|
118 | fHTheta.SetDirectory(NULL);
|
---|
119 | fHTheta.UseCurrentStyle();
|
---|
120 | fHTheta.SetMinimum(0);
|
---|
121 |
|
---|
122 | // effective on time versus time
|
---|
123 | fHTime.SetName("ExcessTime");
|
---|
124 | fHTime.SetTitle(" N_{exc} vs. Time ");
|
---|
125 | fHTime.SetXTitle("Time");
|
---|
126 | fHTime.SetYTitle("N_{exc} [s]");
|
---|
127 | fHTime.UseCurrentStyle();
|
---|
128 | fHTime.SetDirectory(NULL);
|
---|
129 | fHTime.GetYaxis()->SetTitleOffset(1.2);
|
---|
130 | fHTime.GetXaxis()->SetLabelSize(0.033);
|
---|
131 | fHTime.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
|
---|
132 | fHTime.GetXaxis()->SetTimeDisplay(1);
|
---|
133 | fHTime.SetMinimum(0);
|
---|
134 | fHTime.Sumw2();
|
---|
135 |
|
---|
136 | MBinning binsa, binse, binst;
|
---|
137 | binsa.SetEdges(18, 0, 90);
|
---|
138 | binse.SetEdgesLog(15, 10, 100000);
|
---|
139 | binst.SetEdgesASin(67, -0.005, 0.665);
|
---|
140 | binse.Apply(fHEnergy);
|
---|
141 | binst.Apply(fHTheta);
|
---|
142 | binsa.Apply(fHistTime);
|
---|
143 |
|
---|
144 | MH::SetBinning(&fHist, &binst, &binse, &binsa);
|
---|
145 | }
|
---|
146 |
|
---|
147 | Float_t MHAlpha::FitEnergyBins(Bool_t paint)
|
---|
148 | {
|
---|
149 | // Do not store the 'final' result in fFit
|
---|
150 | MAlphaFitter fit(fFit);
|
---|
151 |
|
---|
152 | const Int_t n = fHist.GetNbinsY();
|
---|
153 |
|
---|
154 | fHEnergy.SetEntries(0);
|
---|
155 |
|
---|
156 | Float_t mean = 0;
|
---|
157 | Int_t num = 0;
|
---|
158 | for (int i=1; i<=n; i++)
|
---|
159 | {
|
---|
160 | if (!fit.FitEnergy(fHist, fOffData, i))
|
---|
161 | continue;
|
---|
162 |
|
---|
163 | // FIXME: Calculate UL!
|
---|
164 | if (fit.GetEventsExcess()<=0)
|
---|
165 | continue;
|
---|
166 |
|
---|
167 | fHEnergy.SetBinContent(i, fit.GetEventsExcess());
|
---|
168 | fHEnergy.SetBinError(i, fit.GetErrorExcess());
|
---|
169 |
|
---|
170 | mean += fit.GetEventsExcess()*fit.GetEventsExcess()/fit.GetErrorExcess()/fit.GetErrorExcess();
|
---|
171 | num++;
|
---|
172 | }
|
---|
173 | return TMath::Sqrt(mean)/num;
|
---|
174 | }
|
---|
175 |
|
---|
176 | void MHAlpha::FitThetaBins(Bool_t paint)
|
---|
177 | {
|
---|
178 | // Do not store the 'final' result in fFit
|
---|
179 | MAlphaFitter fit(fFit);
|
---|
180 |
|
---|
181 | const Int_t n = fHist.GetNbinsX();
|
---|
182 |
|
---|
183 | fHTheta.SetEntries(0);
|
---|
184 |
|
---|
185 | for (int i=1; i<=n; i++)
|
---|
186 | {
|
---|
187 | if (!fit.FitTheta(fHist, fOffData, i))
|
---|
188 | continue;
|
---|
189 |
|
---|
190 | // FIXME: Calculate UL!
|
---|
191 | if (fit.GetEventsExcess()<=0)
|
---|
192 | continue;
|
---|
193 |
|
---|
194 | fHTheta.SetBinContent(i, fit.GetEventsExcess());
|
---|
195 | fHTheta.SetBinError(i, fit.GetErrorExcess());
|
---|
196 | }
|
---|
197 | }
|
---|
198 |
|
---|
199 | // --------------------------------------------------------------------------
|
---|
200 | //
|
---|
201 | // Return the value from fParemeter which should be filled into the plots
|
---|
202 | //
|
---|
203 | Double_t MHAlpha::GetVal() const
|
---|
204 | {
|
---|
205 | return static_cast<const MHillasSrc*>(fParameter)->GetAlpha();
|
---|
206 | }
|
---|
207 |
|
---|
208 |
|
---|
209 | // --------------------------------------------------------------------------
|
---|
210 | //
|
---|
211 | // Store the pointer to the parameter container storing the plotted value
|
---|
212 | // (here MHillasSrc) in fParameter.
|
---|
213 | //
|
---|
214 | // return whether it was found or not.
|
---|
215 | //
|
---|
216 | Bool_t MHAlpha::GetParameter(const MParList &pl)
|
---|
217 | {
|
---|
218 | fParameter = (MParContainer*)pl.FindObject(fNameParameter, "MHillasSrc");
|
---|
219 | if (fParameter)
|
---|
220 | return kTRUE;
|
---|
221 |
|
---|
222 | *fLog << err << fNameParameter << " [MHillasSrc] not found... abort." << endl;
|
---|
223 | return kFALSE;
|
---|
224 | }
|
---|
225 |
|
---|
226 | Bool_t MHAlpha::SetupFill(const MParList *pl)
|
---|
227 | {
|
---|
228 | fHist.Reset();
|
---|
229 | fHEnergy.Reset();
|
---|
230 | fHTheta.Reset();
|
---|
231 | fHTime.Reset();
|
---|
232 |
|
---|
233 | const TString off(MString::Format("%sOff", fName.Data()));
|
---|
234 | if (fName!=off && fOffData==NULL)
|
---|
235 | {
|
---|
236 | const TString desc(MString::Format("%s [%s] found... using ", off.Data(), ClassName()));
|
---|
237 | MHAlpha *hoff = (MHAlpha*)pl->FindObject(off, ClassName());
|
---|
238 | if (!hoff)
|
---|
239 | *fLog << inf << "No " << desc << "current data only!" << endl;
|
---|
240 | else
|
---|
241 | {
|
---|
242 | *fLog << inf << desc << "on-off mode!" << endl;
|
---|
243 | SetOffData(*hoff);
|
---|
244 | }
|
---|
245 | }
|
---|
246 |
|
---|
247 | if (!GetParameter(*pl))
|
---|
248 | return kFALSE;
|
---|
249 |
|
---|
250 | fHillas = 0;
|
---|
251 | fEnergy = fForceUsingSize ? 0 : (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
|
---|
252 | if (!fEnergy)
|
---|
253 | {
|
---|
254 | *fLog << warn << "MEnergyEst [MParameterD] not found... " << flush;
|
---|
255 |
|
---|
256 | if (!fHillas)
|
---|
257 | fHillas = (MHillas*)pl->FindObject("MHillas");
|
---|
258 | if (fHillas)
|
---|
259 | *fLog << "using SIZE instead!" << endl;
|
---|
260 | else
|
---|
261 | *fLog << "ignored." << endl;
|
---|
262 |
|
---|
263 | fHEnergy.SetTitle(" N_{exc} vs. Size ");
|
---|
264 | fHEnergy.SetXTitle("Size [phe]");
|
---|
265 | fHist.SetYTitle("Size [phe]");
|
---|
266 | }
|
---|
267 | else
|
---|
268 | {
|
---|
269 | fHEnergy.SetTitle(" N_{exc} vs. E_{est} ");
|
---|
270 | fHEnergy.SetXTitle("E_{est} [GeV]");
|
---|
271 | fHist.SetYTitle("E_{est} [GeV]");
|
---|
272 | }
|
---|
273 |
|
---|
274 | fPointPos = (MPointingPos*)pl->FindObject("MPointingPos");
|
---|
275 | if (!fPointPos)
|
---|
276 | *fLog << warn << "MPointingPos not found... ignored." << endl;
|
---|
277 |
|
---|
278 | fTimeEffOn = (MTime*)pl->FindObject("MTimeEffectiveOnTime", "MTime");
|
---|
279 | if (!fTimeEffOn)
|
---|
280 | *fLog << warn << "MTimeEffectiveOnTime [MTime] not found... ignored." << endl;
|
---|
281 | else
|
---|
282 | *fTimeEffOn = MTime(); // FIXME: How to do it different?
|
---|
283 |
|
---|
284 | fTime = (MTime*)pl->FindObject("MTime");
|
---|
285 | if (!fTime)
|
---|
286 | *fLog << warn << "MTime not found... ignored." << endl;
|
---|
287 |
|
---|
288 | fResult = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MinimizationValue");
|
---|
289 | if (!fResult)
|
---|
290 | return kFALSE;
|
---|
291 | fSigma = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "GaussSigma");
|
---|
292 | if (!fSigma)
|
---|
293 | return kFALSE;
|
---|
294 | fBin = (MParameterI*)const_cast<MParList*>(pl)->FindCreateObj("MParameterI", "Bin");
|
---|
295 | if (!fBin)
|
---|
296 | return kFALSE;
|
---|
297 |
|
---|
298 | //fExcess = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MExcess");
|
---|
299 | //if (!fExcess)
|
---|
300 | // return kFALSE;
|
---|
301 |
|
---|
302 | fLastTime = MTime();
|
---|
303 | fNumRebin = fNumTimeBins;
|
---|
304 |
|
---|
305 | MBinning binst, binse, binsa;
|
---|
306 | binst.SetEdges(fOffData ? *fOffData : fHist, 'x');
|
---|
307 | binse.SetEdges(fOffData ? *fOffData : fHist, 'y');
|
---|
308 | binsa.SetEdges(fOffData ? *fOffData : fHist, 'z');
|
---|
309 | if (!fOffData)
|
---|
310 | {
|
---|
311 | if (!fPointPos)
|
---|
312 | binst.SetEdges(1, 0, 60);
|
---|
313 | else
|
---|
314 | binst.SetEdges(*pl, "BinningTheta");
|
---|
315 |
|
---|
316 | if (!fEnergy && !fHillas)
|
---|
317 | binse.SetEdges(1, 10, 100000);
|
---|
318 | else
|
---|
319 | if (fEnergy)
|
---|
320 | binse.SetEdges(*pl, "BinningEnergyEst");
|
---|
321 | else
|
---|
322 | binse.SetEdges(*pl, "BinningSize");
|
---|
323 |
|
---|
324 | binsa.SetEdges(*pl, MString::Format("Binning%s", ClassName()+2));
|
---|
325 | }
|
---|
326 | else
|
---|
327 | {
|
---|
328 | fHEnergy.SetTitle(fOffData->GetTitle());
|
---|
329 | fHEnergy.SetXTitle(fOffData->GetYaxis()->GetTitle());
|
---|
330 | fHist.SetYTitle(fOffData->GetYaxis()->GetTitle());
|
---|
331 | }
|
---|
332 |
|
---|
333 | binse.Apply(fHEnergy);
|
---|
334 | binst.Apply(fHTheta);
|
---|
335 | binsa.Apply(fHistTime);
|
---|
336 | MH::SetBinning(&fHist, &binst, &binse, &binsa);
|
---|
337 |
|
---|
338 | MAlphaFitter *fit = (MAlphaFitter*)pl->FindObject("MAlphaFitter");
|
---|
339 | if (!fit)
|
---|
340 | *fLog << warn << "MAlphaFitter not found... ignored." << endl;
|
---|
341 | else
|
---|
342 | fit->Copy(fFit);
|
---|
343 |
|
---|
344 | *fLog << inf;
|
---|
345 | fFit.Print();
|
---|
346 |
|
---|
347 | return kTRUE;
|
---|
348 | }
|
---|
349 |
|
---|
350 | void MHAlpha::InitAlphaTime(const MTime &t)
|
---|
351 | {
|
---|
352 | //
|
---|
353 | // If this is the first call we have to initialize the time-histogram
|
---|
354 | //
|
---|
355 | MBinning bins;
|
---|
356 | bins.SetEdges(1, t.GetAxisTime()-60, t.GetAxisTime());
|
---|
357 | bins.Apply(fHTime);
|
---|
358 |
|
---|
359 | fLastTime=t;
|
---|
360 | }
|
---|
361 |
|
---|
362 | void MHAlpha::UpdateAlphaTime(Bool_t final)
|
---|
363 | {
|
---|
364 | if (!fTimeEffOn)
|
---|
365 | return;
|
---|
366 |
|
---|
367 | if (!final)
|
---|
368 | {
|
---|
369 | if (fLastTime==MTime() && *fTimeEffOn!=MTime())
|
---|
370 | {
|
---|
371 | InitAlphaTime(*fTimeEffOn);
|
---|
372 | return;
|
---|
373 | }
|
---|
374 |
|
---|
375 | if (fLastTime!=*fTimeEffOn)
|
---|
376 | {
|
---|
377 | fLastTime=*fTimeEffOn;
|
---|
378 | fNumRebin--;
|
---|
379 | }
|
---|
380 |
|
---|
381 | if (fNumRebin>0)
|
---|
382 | return;
|
---|
383 | }
|
---|
384 |
|
---|
385 | // Check new 'last time'
|
---|
386 | MTime *time = final ? fTime : fTimeEffOn;
|
---|
387 |
|
---|
388 | if (time->GetAxisTime()<=fHTime.GetXaxis()->GetXmax())
|
---|
389 | {
|
---|
390 | *fLog << warn << "WARNING - New time-stamp " << *time << " lower" << endl;
|
---|
391 | *fLog << "than upper edge of histogram... skipped." << endl;
|
---|
392 | *fLog << "This should not happen. Maybe you started you eventloop with" << endl;
|
---|
393 | *fLog << "an already initialized time stamp MTimeEffectiveOnTime?" << endl;
|
---|
394 | fNumRebin++;
|
---|
395 | return;
|
---|
396 | }
|
---|
397 |
|
---|
398 | // Fit time histogram
|
---|
399 | MAlphaFitter fit(fFit);
|
---|
400 |
|
---|
401 | TH1D *h = fOffData ? fOffData->ProjectionZ("ProjTimeTemp", 0, fOffData->GetNbinsX()+1, 0, fOffData->GetNbinsY()+1, "E") : 0;
|
---|
402 | const Bool_t rc = fit.ScaleAndFit(fHistTime, h);
|
---|
403 |
|
---|
404 | if (h)
|
---|
405 | delete h;
|
---|
406 |
|
---|
407 | if (!rc)
|
---|
408 | return;
|
---|
409 |
|
---|
410 | // Reset Histogram
|
---|
411 | fHistTime.Reset();
|
---|
412 | fHistTime.SetEntries(0);
|
---|
413 |
|
---|
414 | //
|
---|
415 | // Prepare Histogram
|
---|
416 | //
|
---|
417 | if (final)
|
---|
418 | time->Plus1ns();
|
---|
419 |
|
---|
420 | // Enhance binning
|
---|
421 | MBinning bins;
|
---|
422 | bins.SetEdges(fHTime, 'x');
|
---|
423 | bins.AddEdge(time->GetAxisTime());
|
---|
424 | bins.Apply(fHTime);
|
---|
425 |
|
---|
426 | //
|
---|
427 | // Fill histogram
|
---|
428 | //
|
---|
429 | // Get number of bins
|
---|
430 | const Int_t n = fHTime.GetNbinsX();
|
---|
431 |
|
---|
432 | fHTime.SetBinContent(n, fit.GetEventsExcess());
|
---|
433 | fHTime.SetBinError(n, fit.GetErrorExcess());
|
---|
434 |
|
---|
435 | //*fLog << inf3 << *fTimeEffOn << " (" << n << "): " << fit.GetEventsExcess() << endl;
|
---|
436 |
|
---|
437 | fNumRebin = fNumTimeBins;
|
---|
438 | }
|
---|
439 |
|
---|
440 | void MHAlpha::SetBin(Int_t ibin)
|
---|
441 | {
|
---|
442 | // Is this necessary?
|
---|
443 | // Could be speed up up searching for it only once.
|
---|
444 | const Float_t max = fFit.GetSignalIntegralMax();
|
---|
445 | const Int_t bin0 = fHist.GetZaxis()->FindFixBin(max);
|
---|
446 |
|
---|
447 | const Int_t nbinsx = fHist.GetNbinsX();
|
---|
448 | const Int_t nbinsy = fHist.GetNbinsY();
|
---|
449 | const Int_t nxy = (nbinsx+2)*(nbinsy+2);
|
---|
450 |
|
---|
451 | const Int_t binz = ibin/nxy;
|
---|
452 |
|
---|
453 | const Bool_t issignal = binz>0 && binz<bin0;
|
---|
454 |
|
---|
455 | fBin->SetVal(issignal ? binz : -binz);
|
---|
456 | }
|
---|
457 |
|
---|
458 | // --------------------------------------------------------------------------
|
---|
459 | //
|
---|
460 | // Fill the histogram. For details see the code or the class description
|
---|
461 | //
|
---|
462 | Int_t MHAlpha::Fill(const MParContainer *par, const Stat_t w)
|
---|
463 | {
|
---|
464 | Double_t alpha, energy, theta;
|
---|
465 | Double_t size=-1;
|
---|
466 |
|
---|
467 | if (fMatrix)
|
---|
468 | {
|
---|
469 | alpha = fMap[0]<0 ? GetVal() : (*fMatrix)[fMap[0]];
|
---|
470 | energy = fMap[1]<0 ? -1 : (*fMatrix)[fMap[1]];
|
---|
471 | size = fMap[2]<0 ? -1 : (*fMatrix)[fMap[2]];
|
---|
472 | //<0 ? 1000 : (*fMatrix)[fMap[1]];
|
---|
473 | theta = 0;
|
---|
474 |
|
---|
475 | if (energy<0)
|
---|
476 | energy=size;
|
---|
477 | if (size<0)
|
---|
478 | size=energy;
|
---|
479 |
|
---|
480 | if (energy<0 && size<0)
|
---|
481 | energy = size = 1000;
|
---|
482 | }
|
---|
483 | else
|
---|
484 | {
|
---|
485 | alpha = GetVal();
|
---|
486 |
|
---|
487 | if (fHillas)
|
---|
488 | size = fHillas->GetSize();
|
---|
489 | energy = fEnergy ? fEnergy->GetVal() : (fHillas?fHillas->GetSize():1000);
|
---|
490 | theta = fPointPos ? fPointPos->GetZd() : 0;
|
---|
491 | }
|
---|
492 |
|
---|
493 | // enhance histogram if necessary
|
---|
494 | UpdateAlphaTime();
|
---|
495 |
|
---|
496 | // Fill histograms
|
---|
497 | const Int_t ibin = fHist.Fill(theta, energy, TMath::Abs(alpha), w);
|
---|
498 | SetBin(ibin);
|
---|
499 |
|
---|
500 | if (!fSkipHistTime)
|
---|
501 | fHistTime.Fill(TMath::Abs(alpha), w);
|
---|
502 |
|
---|
503 | return kTRUE;
|
---|
504 | }
|
---|
505 |
|
---|
506 | // --------------------------------------------------------------------------
|
---|
507 | //
|
---|
508 | // Paint the integral and the error on top of the histogram
|
---|
509 | //
|
---|
510 | void MHAlpha::PaintText(Double_t val, Double_t error) const
|
---|
511 | {
|
---|
512 | TLatex text(0.45, 0.94, MString::Format("N_{exc} = %.1f \\pm %.1f", val, error));
|
---|
513 | text.SetBit(TLatex::kTextNDC);
|
---|
514 | text.SetTextSize(0.04);
|
---|
515 | text.Paint();
|
---|
516 | }
|
---|
517 |
|
---|
518 | // --------------------------------------------------------------------------
|
---|
519 | //
|
---|
520 | // Paint the integral and the error on top of the histogram
|
---|
521 | //
|
---|
522 | void MHAlpha::PaintText(const TH1D &h) const
|
---|
523 | {
|
---|
524 | Double_t sumv = 0;
|
---|
525 | Double_t sume = 0;
|
---|
526 |
|
---|
527 | for (int i=0; i<h.GetNbinsX(); i++)
|
---|
528 | {
|
---|
529 | sumv += h.GetBinContent(i+1);
|
---|
530 | sume += h.GetBinError(i+1);
|
---|
531 | }
|
---|
532 | PaintText(sumv, sume);
|
---|
533 | }
|
---|
534 | // --------------------------------------------------------------------------
|
---|
535 | //
|
---|
536 | // Update the projections
|
---|
537 | //
|
---|
538 | void MHAlpha::Paint(Option_t *opt)
|
---|
539 | {
|
---|
540 | // Note: Any object cannot be added to a pad twice!
|
---|
541 | // The object is searched by gROOT->FindObject only in
|
---|
542 | // gPad itself!
|
---|
543 | TVirtualPad *padsave = gPad;
|
---|
544 |
|
---|
545 | TH1D *h0=0;
|
---|
546 |
|
---|
547 | TString o(opt);
|
---|
548 |
|
---|
549 | if (o==(TString)"proj")
|
---|
550 | {
|
---|
551 | TPaveStats *st=0;
|
---|
552 | for (int x=0; x<4; x++)
|
---|
553 | {
|
---|
554 | TVirtualPad *p=padsave->GetPad(x+1);
|
---|
555 | if (!p || !(st = (TPaveStats*)p->GetPrimitive("stats")))
|
---|
556 | continue;
|
---|
557 |
|
---|
558 | if (st->GetOptStat()==11)
|
---|
559 | continue;
|
---|
560 |
|
---|
561 | const Double_t y1 = st->GetY1NDC();
|
---|
562 | const Double_t y2 = st->GetY2NDC();
|
---|
563 | const Double_t x1 = st->GetX1NDC();
|
---|
564 | const Double_t x2 = st->GetX2NDC();
|
---|
565 |
|
---|
566 | st->SetY1NDC((y2-y1)/3+y1);
|
---|
567 | st->SetX1NDC((x2-x1)/3+x1);
|
---|
568 | st->SetOptStat(11);
|
---|
569 | }
|
---|
570 |
|
---|
571 | padsave->cd(1);
|
---|
572 |
|
---|
573 | TH1D *hon = (TH1D*)gPad->FindObject("Proj");
|
---|
574 | if (hon)
|
---|
575 | {
|
---|
576 | TH1D *dum = fHist.ProjectionZ("dumab", 0, fHist.GetNbinsX()+1, 0, fHist.GetNbinsY()+1);
|
---|
577 | dum->SetDirectory(0);
|
---|
578 | hon->Reset();
|
---|
579 | hon->Add(dum);
|
---|
580 | delete dum;
|
---|
581 |
|
---|
582 | if (fOffData)
|
---|
583 | {
|
---|
584 | TH1D *hoff = (TH1D*)gPad->FindObject("ProjOff");
|
---|
585 | if (hoff)
|
---|
586 | {
|
---|
587 | TH1D *dum = fOffData->ProjectionZ("dumxy", 0, fOffData->GetNbinsX()+1, 0, fOffData->GetNbinsY()+1);
|
---|
588 | dum->SetDirectory(0);
|
---|
589 | hoff->Reset();
|
---|
590 | hoff->Add(dum);
|
---|
591 | delete dum;
|
---|
592 |
|
---|
593 | const Double_t alpha = fFit.Scale(*hoff, *hon);
|
---|
594 |
|
---|
595 | hon->SetMaximum();
|
---|
596 | hon->SetMaximum(TMath::Max(hon->GetMaximum(), hoff->GetMaximum())*1.05);
|
---|
597 |
|
---|
598 | // BE CARFEULL: This is a really weird workaround!
|
---|
599 | hoff->SetMaximum(alpha);
|
---|
600 |
|
---|
601 | // For some reason the line-color is resetted
|
---|
602 | hoff->SetLineColor(kRed);
|
---|
603 |
|
---|
604 | if ((h0=(TH1D*)gPad->FindObject("ProjOnOff")))
|
---|
605 | {
|
---|
606 | h0->Reset();
|
---|
607 | h0->Add(hoff, hon, -1);
|
---|
608 | const Float_t min = h0->GetMinimum()*1.05;
|
---|
609 | hon->SetMinimum(min<0 ? min : 0);
|
---|
610 | }
|
---|
611 |
|
---|
612 | }
|
---|
613 | }
|
---|
614 | else
|
---|
615 | hon->SetMinimum(0);
|
---|
616 | }
|
---|
617 | FitEnergyBins();
|
---|
618 | FitThetaBins();
|
---|
619 | }
|
---|
620 |
|
---|
621 | if (o==(TString)"variable")
|
---|
622 | if ((h0 = (TH1D*)gPad->FindObject("Proj")))
|
---|
623 | {
|
---|
624 | // Check whether we are dealing with on-off analysis
|
---|
625 | TH1D *hoff = (TH1D*)gPad->FindObject("ProjOff");
|
---|
626 |
|
---|
627 | // BE CARFEULL: This is a really weird workaround!
|
---|
628 | const Double_t scale = !hoff || hoff->GetMaximum()<0 ? 1 : hoff->GetMaximum();
|
---|
629 |
|
---|
630 | // Do not store the 'final' result in fFit
|
---|
631 | MAlphaFitter fit(fFit);
|
---|
632 | fit.Fit(*h0, hoff, scale, kTRUE);
|
---|
633 | fit.PaintResult();
|
---|
634 | }
|
---|
635 |
|
---|
636 | if (o==(TString)"time")
|
---|
637 | PaintText(fHTime);
|
---|
638 |
|
---|
639 | if (o==(TString)"theta")
|
---|
640 | {
|
---|
641 | TH1 *h = (TH1*)gPad->FindObject(MString::Format("%s_x", fHist.GetName()));
|
---|
642 | if (h)
|
---|
643 | {
|
---|
644 | TH1D *h2 = (TH1D*)fHist.Project3D("dum_x");
|
---|
645 | h2->SetDirectory(0);
|
---|
646 | h2->Scale(fHTheta.Integral()/h2->Integral());
|
---|
647 | h->Reset();
|
---|
648 | h->Add(h2);
|
---|
649 | delete h2;
|
---|
650 | }
|
---|
651 | PaintText(fHTheta);
|
---|
652 | }
|
---|
653 |
|
---|
654 | if (o==(TString)"energy")
|
---|
655 | {
|
---|
656 | TH1 *h = (TH1*)gPad->FindObject(MString::Format("%s_y", fHist.GetName()));
|
---|
657 | if (h)
|
---|
658 | {
|
---|
659 | TH1D *h2 = (TH1D*)fHist.Project3D("dum_y");
|
---|
660 | h2->SetDirectory(0);
|
---|
661 | h2->Scale(fHEnergy.Integral()/h2->Integral());
|
---|
662 | h->Reset();
|
---|
663 | h->Add(h2);
|
---|
664 | delete h2;
|
---|
665 | }
|
---|
666 | PaintText(fHEnergy);
|
---|
667 |
|
---|
668 | if (fHEnergy.GetMaximum()>1)
|
---|
669 | {
|
---|
670 | gPad->SetLogx();
|
---|
671 | gPad->SetLogy();
|
---|
672 | }
|
---|
673 | }
|
---|
674 |
|
---|
675 | gPad = padsave;
|
---|
676 | }
|
---|
677 |
|
---|
678 | // --------------------------------------------------------------------------
|
---|
679 | //
|
---|
680 | // Draw the histogram
|
---|
681 | //
|
---|
682 | void MHAlpha::Draw(Option_t *opt)
|
---|
683 | {
|
---|
684 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
685 |
|
---|
686 | /*
|
---|
687 | if (TString(opt).Contains("sizebins", TString::kIgnoreCase))
|
---|
688 | {
|
---|
689 | AppendPad("sizebins");
|
---|
690 | return;
|
---|
691 | }
|
---|
692 | */
|
---|
693 |
|
---|
694 | // Do the projection before painting the histograms into
|
---|
695 | // the individual pads
|
---|
696 | AppendPad("proj");
|
---|
697 |
|
---|
698 | pad->SetBorderMode(0);
|
---|
699 | pad->Divide(2,2);
|
---|
700 |
|
---|
701 | TH1D *h=0;
|
---|
702 |
|
---|
703 | pad->cd(1);
|
---|
704 | gPad->SetBorderMode(0);
|
---|
705 |
|
---|
706 | h = fHist.ProjectionZ("Proj", 0, fHist.GetNbinsX()+1, 0, fHist.GetNbinsY()+1, "E");
|
---|
707 | h->SetBit(TH1::kNoTitle);
|
---|
708 | h->SetStats(kTRUE);
|
---|
709 | h->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
710 | h->SetYTitle("Counts");
|
---|
711 | h->SetDirectory(NULL);
|
---|
712 | h->SetMarkerStyle(0);
|
---|
713 | h->SetBit(kCanDelete);
|
---|
714 | h->Draw("");
|
---|
715 |
|
---|
716 | if (fOffData)
|
---|
717 | {
|
---|
718 | // To get green on-data
|
---|
719 | //h->SetMarkerColor(kGreen);
|
---|
720 | //h->SetLineColor(kGreen);
|
---|
721 |
|
---|
722 | h = fOffData->ProjectionZ("ProjOff", 0, fOffData->GetNbinsX()+1, 0, fOffData->GetNbinsY()+1, "E");
|
---|
723 | h->SetBit(TH1::kNoTitle);
|
---|
724 | h->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
725 | h->SetYTitle("Counts");
|
---|
726 | h->SetDirectory(NULL);
|
---|
727 | h->SetMarkerStyle(0);
|
---|
728 | h->SetBit(kCanDelete);
|
---|
729 | h->SetMarkerColor(kRed);
|
---|
730 | h->SetLineColor(kRed);
|
---|
731 | //h->SetFillColor(18);
|
---|
732 | h->Draw("same"/*"bar same"*/);
|
---|
733 |
|
---|
734 | // This is the only way to make it work...
|
---|
735 | // Clone and copy constructor give strange crashes :(
|
---|
736 | h = fOffData->ProjectionZ("ProjOnOff", 0, fOffData->GetNbinsX()+1, 0, fOffData->GetNbinsY()+1, "E");
|
---|
737 | h->SetBit(TH1::kNoTitle);
|
---|
738 | h->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
739 | h->SetYTitle("Counts");
|
---|
740 | h->SetDirectory(NULL);
|
---|
741 | h->SetMarkerStyle(0);
|
---|
742 | h->SetBit(kCanDelete);
|
---|
743 | h->SetMarkerColor(kBlue);
|
---|
744 | h->SetLineColor(kBlue);
|
---|
745 | h->Draw("same");
|
---|
746 |
|
---|
747 | TLine lin;
|
---|
748 | lin.SetLineStyle(kDashed);
|
---|
749 | lin.DrawLine(h->GetXaxis()->GetXmin(), 0, h->GetXaxis()->GetXmax(), 0);
|
---|
750 | }
|
---|
751 |
|
---|
752 | // After the Alpha-projection has been drawn. Fit the histogram
|
---|
753 | // and paint the result into this pad
|
---|
754 | AppendPad("variable");
|
---|
755 |
|
---|
756 | if (fHEnergy.GetNbinsX()>1 || fHEnergy.GetBinContent(1)>0)
|
---|
757 | {
|
---|
758 | pad->cd(2);
|
---|
759 | gPad->SetBorderMode(0);
|
---|
760 | gPad->SetGridx();
|
---|
761 | gPad->SetGridy();
|
---|
762 | fHEnergy.Draw();
|
---|
763 |
|
---|
764 | AppendPad("energy");
|
---|
765 |
|
---|
766 | h = (TH1D*)fHist.Project3D("y");
|
---|
767 | h->SetBit(TH1::kNoTitle|TH1::kNoStats);
|
---|
768 | h->SetXTitle("E [GeV]");
|
---|
769 | h->SetYTitle("Counts");
|
---|
770 | h->SetDirectory(NULL);
|
---|
771 | h->SetMarkerStyle(kFullDotSmall);
|
---|
772 | h->SetBit(kCanDelete);
|
---|
773 | h->SetMarkerColor(kCyan);
|
---|
774 | h->SetLineColor(kCyan);
|
---|
775 | h->Draw("Psame");
|
---|
776 | }
|
---|
777 | else
|
---|
778 | delete pad->GetPad(2);
|
---|
779 |
|
---|
780 | if ((fTimeEffOn && fTime) || fHTime.GetNbinsX()>1 || fHTime.GetBinError(1)>0)
|
---|
781 | {
|
---|
782 | pad->cd(3);
|
---|
783 | gPad->SetBorderMode(0);
|
---|
784 | gPad->SetGridx();
|
---|
785 | gPad->SetGridy();
|
---|
786 | fHTime.Draw();
|
---|
787 | AppendPad("time");
|
---|
788 | }
|
---|
789 | else
|
---|
790 | delete pad->GetPad(3);
|
---|
791 |
|
---|
792 | if (fHTheta.GetNbinsX()>1 || fHTheta.GetBinContent(1)>0)
|
---|
793 | {
|
---|
794 | pad->cd(4);
|
---|
795 | gPad->SetGridx();
|
---|
796 | gPad->SetGridy();
|
---|
797 | gPad->SetBorderMode(0);
|
---|
798 | fHTheta.Draw();
|
---|
799 |
|
---|
800 | AppendPad("theta");
|
---|
801 |
|
---|
802 | h = (TH1D*)fHist.Project3D("x");
|
---|
803 | h->SetBit(TH1::kNoTitle|TH1::kNoStats);
|
---|
804 | h->SetXTitle("\\theta [\\circ]");
|
---|
805 | h->SetYTitle("Counts");
|
---|
806 | h->SetDirectory(NULL);
|
---|
807 | h->SetMarkerStyle(kFullDotSmall);
|
---|
808 | h->SetBit(kCanDelete);
|
---|
809 | h->SetMarkerColor(kCyan);
|
---|
810 | h->SetLineColor(kCyan);
|
---|
811 | h->Draw("Psame");
|
---|
812 | }
|
---|
813 | else
|
---|
814 | delete pad->GetPad(4);
|
---|
815 | }
|
---|
816 |
|
---|
817 | void MHAlpha::DrawAll(Bool_t newc)
|
---|
818 | {
|
---|
819 | if (!newc && !fDisplay)
|
---|
820 | return;
|
---|
821 |
|
---|
822 | // FIXME: Do in Paint if special option given!
|
---|
823 | TCanvas &c = newc ? *new TCanvas : fDisplay->AddTab("SizeBins");
|
---|
824 | Int_t n = fHist.GetNbinsY();
|
---|
825 | Int_t nc = (Int_t)(TMath::Sqrt((Float_t)n-1)+1);
|
---|
826 | c.Divide(nc, nc, 1e-10, 1e-10);
|
---|
827 | gPad->SetBorderMode(0);
|
---|
828 | gPad->SetFrameBorderMode(0);
|
---|
829 |
|
---|
830 | // Do not store the 'final' result in fFit
|
---|
831 | MAlphaFitter fit(fFit);
|
---|
832 |
|
---|
833 | for (int i=1; i<=fHist.GetNbinsY(); i++)
|
---|
834 | {
|
---|
835 | c.cd(i);
|
---|
836 | gPad->SetBorderMode(0);
|
---|
837 | gPad->SetFrameBorderMode(0);
|
---|
838 |
|
---|
839 | TH1D *hon = fHist.ProjectionZ("Proj", 0, fHist.GetNbinsX()+1, i, i, "E");
|
---|
840 | hon->SetBit(TH1::kNoTitle);
|
---|
841 | hon->SetStats(kTRUE);
|
---|
842 | hon->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
843 | hon->SetYTitle("Counts");
|
---|
844 | hon->SetDirectory(NULL);
|
---|
845 | hon->SetMarkerStyle(0);
|
---|
846 | hon->SetBit(kCanDelete);
|
---|
847 | hon->Draw("");
|
---|
848 |
|
---|
849 | TH1D *hof = 0;
|
---|
850 | Double_t alpha = 1;
|
---|
851 |
|
---|
852 | if (fOffData)
|
---|
853 | {
|
---|
854 | hon->SetMarkerColor(kGreen);
|
---|
855 |
|
---|
856 | hof = fOffData->ProjectionZ("ProjOff", 0, fOffData->GetNbinsX()+1, i, i, "E");
|
---|
857 | hof->SetBit(TH1::kNoTitle|TH1::kNoStats);
|
---|
858 | hof->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
859 | hof->SetYTitle("Counts");
|
---|
860 | hof->SetDirectory(NULL);
|
---|
861 | hof->SetMarkerStyle(0);
|
---|
862 | hof->SetBit(kCanDelete);
|
---|
863 | hof->SetMarkerColor(kRed);
|
---|
864 | hof->SetLineColor(kRed);
|
---|
865 | hof->Draw("same");
|
---|
866 |
|
---|
867 | alpha = fit.Scale(*hof, *hon);
|
---|
868 |
|
---|
869 | hon->SetMaximum();
|
---|
870 | hon->SetMaximum(TMath::Max(hon->GetMaximum(), hof->GetMaximum())*1.05);
|
---|
871 |
|
---|
872 | TH1D *diff = new TH1D(*hon);
|
---|
873 | diff->Add(hof, -1);
|
---|
874 | diff->SetBit(TH1::kNoTitle);
|
---|
875 | diff->SetXTitle(fHist.GetZaxis()->GetTitle());
|
---|
876 | diff->SetYTitle("Counts");
|
---|
877 | diff->SetDirectory(NULL);
|
---|
878 | diff->SetMarkerStyle(0);
|
---|
879 | diff->SetBit(kCanDelete);
|
---|
880 | diff->SetMarkerColor(kBlue);
|
---|
881 | diff->Draw("same");
|
---|
882 |
|
---|
883 | TLine lin;
|
---|
884 | lin.SetLineStyle(kDashed);
|
---|
885 | lin.DrawLine(diff->GetXaxis()->GetXmin(), 0, diff->GetXaxis()->GetXmax(), 0);
|
---|
886 |
|
---|
887 | const Float_t min = diff->GetMinimum()*1.05;
|
---|
888 | hon->SetMinimum(min<0 ? min : 0);
|
---|
889 | }
|
---|
890 |
|
---|
891 | if (hof ? fit.Fit(*hon, *hof, alpha) : fit.Fit(*hon))
|
---|
892 | {
|
---|
893 | *fLog << dbg << "Bin " << i << ": sigmaexc=" << fit.GetEventsExcess()/fit.GetErrorExcess() << " omega=" << fit.GetGausSigma() << " events=" << fit.GetEventsExcess() << " scale=" << fit.GetScaleFactor() << endl;
|
---|
894 | fit.DrawResult();
|
---|
895 | }
|
---|
896 | /*
|
---|
897 | if (fit.FitEnergy(fHist, fOffData, i, kTRUE))
|
---|
898 | {
|
---|
899 | fHEnergy.SetBinContent(i, fit.GetEventsExcess());
|
---|
900 | fHEnergy.SetBinError(i, fit.GetEventsExcess()*0.2);
|
---|
901 | }*/
|
---|
902 | }
|
---|
903 |
|
---|
904 | }
|
---|
905 |
|
---|
906 | void MHAlpha::DrawNicePlot(Bool_t newc, const char *title, const char *watermark, Int_t binlo, Int_t binhi)
|
---|
907 | {
|
---|
908 | if (!newc && !fDisplay)
|
---|
909 | return;
|
---|
910 |
|
---|
911 | if (!fOffData)
|
---|
912 | return;
|
---|
913 |
|
---|
914 | // Open and setup canvas/pad
|
---|
915 | TCanvas &c = newc ? *new TCanvas : fDisplay->AddTab("ThetsSq");
|
---|
916 |
|
---|
917 | c.SetBorderMode(0);
|
---|
918 | c.SetFrameBorderMode(0);
|
---|
919 | c.SetFillColor(kWhite);
|
---|
920 |
|
---|
921 | c.SetLeftMargin(0.12);
|
---|
922 | c.SetRightMargin(0.01);
|
---|
923 | c.SetBottomMargin(0.16);
|
---|
924 | c.SetTopMargin(0.18);
|
---|
925 |
|
---|
926 | c.SetGridy();
|
---|
927 |
|
---|
928 | gStyle->SetOptStat(0);
|
---|
929 |
|
---|
930 | // Get on-data
|
---|
931 | TH1D *hon = (TH1D*)fHist.ProjectionZ("Proj", 0, fHist.GetNbinsX()+1, binlo, binhi, "E");
|
---|
932 | hon->SetDirectory(NULL);
|
---|
933 | hon->SetBit(kCanDelete);
|
---|
934 | hon->SetMarkerSize(0);
|
---|
935 | hon->SetLineWidth(2);
|
---|
936 | hon->SetLineColor(kBlack);
|
---|
937 | hon->SetMarkerColor(kBlack);
|
---|
938 |
|
---|
939 | // Get off-data
|
---|
940 | TH1D *hoff = 0;
|
---|
941 | if (fOffData)
|
---|
942 | {
|
---|
943 | hoff = (TH1D*)fOffData->ProjectionZ("ProjOff", 0, fHist.GetNbinsX()+1, binlo, binhi, "E");
|
---|
944 | hoff->SetDirectory(NULL);
|
---|
945 | hoff->SetBit(kCanDelete);
|
---|
946 | hoff->SetFillColor(17);
|
---|
947 | hoff->SetMarkerSize(0);
|
---|
948 | hoff->SetLineColor(kBlack);
|
---|
949 | hoff->SetMarkerColor(kBlack);
|
---|
950 | }
|
---|
951 |
|
---|
952 | // Setup plot which is drawn first
|
---|
953 | TH1D *h = hoff ? hoff : hon;
|
---|
954 | h->GetXaxis()->SetLabelSize(0.06);
|
---|
955 | h->GetXaxis()->SetTitleSize(0.06);
|
---|
956 | h->GetXaxis()->SetTitleOffset(0.95);
|
---|
957 | h->GetYaxis()->SetLabelSize(0.06);
|
---|
958 | h->GetYaxis()->SetTitleSize(0.06);
|
---|
959 | h->GetYaxis()->CenterTitle();
|
---|
960 | h->SetYTitle("Counts");
|
---|
961 | h->SetTitleSize(0.07);
|
---|
962 | h->SetTitle("");
|
---|
963 |
|
---|
964 | const Double_t imax = fFit.GetSignalIntegralMax();
|
---|
965 | if (imax<1)
|
---|
966 | h->GetXaxis()->SetRangeUser(0, 0.6*0.6);
|
---|
967 |
|
---|
968 | // scale off-data
|
---|
969 |
|
---|
970 | MAlphaFitter fit(fFit);
|
---|
971 | fit.ScaleAndFit(*hon, hoff);
|
---|
972 |
|
---|
973 | hon->SetMinimum(0);
|
---|
974 | hoff->SetMinimum(0);
|
---|
975 |
|
---|
976 | // draw data
|
---|
977 | if (hoff)
|
---|
978 | {
|
---|
979 | hoff->SetMaximum(TMath::Max(hon->GetMaximum(),hoff->GetMaximum())*1.1);
|
---|
980 | hoff->Draw("bar");
|
---|
981 | hon->Draw("same");
|
---|
982 | }
|
---|
983 | else
|
---|
984 | {
|
---|
985 | hon->SetMaximum();
|
---|
986 | hon->Draw();
|
---|
987 | }
|
---|
988 |
|
---|
989 | // draw a preliminary tag
|
---|
990 | TLatex text;
|
---|
991 | text.SetTextColor(kWhite);
|
---|
992 | text.SetBit(TLatex::kTextNDC);
|
---|
993 | text.SetTextSize(0.07);
|
---|
994 | text.SetTextAngle(2.5);
|
---|
995 |
|
---|
996 | TString wm(watermark);
|
---|
997 | if (binlo>=1 || binhi<hon->GetNbinsX())
|
---|
998 | {
|
---|
999 | wm += wm.IsNull() ? "(" : " (";
|
---|
1000 | if (binlo>=1)
|
---|
1001 | wm += MString::Format("%.1fGeV", fHist.GetYaxis()->GetBinLowEdge(binlo));
|
---|
1002 | wm += "-";
|
---|
1003 | if (binhi<hon->GetNbinsX())
|
---|
1004 | wm += MString::Format("%.1fGeV", fHist.GetYaxis()->GetBinLowEdge(binhi+1));
|
---|
1005 | wm += ")";
|
---|
1006 | }
|
---|
1007 | if (!wm.IsNull())
|
---|
1008 | text.DrawLatex(0.45, 0.2, wm);
|
---|
1009 | //enum { kTextNDC = BIT(14) };
|
---|
1010 |
|
---|
1011 | // draw line showing cut
|
---|
1012 | TLine line;
|
---|
1013 | line.SetLineColor(14);
|
---|
1014 | line.SetLineStyle(7);
|
---|
1015 | line.DrawLine(imax, 0, imax, h->GetMaximum());
|
---|
1016 |
|
---|
1017 | // Add a title above the plot
|
---|
1018 | TPaveText *pave=new TPaveText(0.12, 0.83, 0.99, 0.975, "blNDC");
|
---|
1019 | pave->SetBorderSize(1);
|
---|
1020 | pave->SetLabel(title);
|
---|
1021 |
|
---|
1022 | TText *ptxt = pave->AddText(" ");
|
---|
1023 | ptxt->SetTextAlign(23);
|
---|
1024 |
|
---|
1025 | ptxt = pave->AddText(MString::Format("Significance %.1f\\sigma, off-scale %.2f",
|
---|
1026 | fit.GetSignificance(), fit.GetScaleFactor()));
|
---|
1027 | ptxt->SetTextAlign(23);
|
---|
1028 |
|
---|
1029 | ptxt = pave->AddText(MString::Format("%.1f excess events, %.1f background events",
|
---|
1030 | fit.GetEventsExcess(), fit.GetEventsBackground()));
|
---|
1031 | ptxt->SetTextAlign(23);
|
---|
1032 | pave->SetBit(kCanDelete);
|
---|
1033 | pave->Draw();
|
---|
1034 | }
|
---|
1035 |
|
---|
1036 | Bool_t MHAlpha::Finalize()
|
---|
1037 | {
|
---|
1038 | if (!FitAlpha())
|
---|
1039 | {
|
---|
1040 | *fLog << warn << "MAlphaFitter - Fit failed..." << endl;
|
---|
1041 | return kTRUE;
|
---|
1042 | }
|
---|
1043 |
|
---|
1044 | // Store the final result in fFit
|
---|
1045 | *fLog << all;
|
---|
1046 | fFit.Print("result");
|
---|
1047 |
|
---|
1048 | fResult->SetVal(fFit.GetMinimizationValue());
|
---|
1049 | fSigma->SetVal(fFit.GetGausSigma());
|
---|
1050 |
|
---|
1051 | if (!fSkipHistEnergy)
|
---|
1052 | {
|
---|
1053 | *fLog << inf3 << "Processing energy bins..." << endl;
|
---|
1054 | FitEnergyBins();
|
---|
1055 | }
|
---|
1056 | if (!fSkipHistTheta)
|
---|
1057 | {
|
---|
1058 | *fLog << inf3 << "Processing theta bins..." << endl;
|
---|
1059 | FitThetaBins();
|
---|
1060 | }
|
---|
1061 | if (!fSkipHistTime)
|
---|
1062 | {
|
---|
1063 | *fLog << inf3 << "Processing time bins..." << endl;
|
---|
1064 | UpdateAlphaTime(kTRUE);
|
---|
1065 | MH::RemoveFirstBin(fHTime);
|
---|
1066 | }
|
---|
1067 |
|
---|
1068 | if (fOffData)
|
---|
1069 | DrawAll(kFALSE);
|
---|
1070 |
|
---|
1071 | return kTRUE;
|
---|
1072 | }
|
---|
1073 |
|
---|
1074 | // --------------------------------------------------------------------------
|
---|
1075 | //
|
---|
1076 | // You can use this function if you want to use a MHMatrix instead of
|
---|
1077 | // MMcEvt. This function adds all necessary columns to the
|
---|
1078 | // given matrix. Afterward you should fill the matrix with the corresponding
|
---|
1079 | // data (eg from a file by using MHMatrix::Fill). If you now loop
|
---|
1080 | // through the matrix (eg using MMatrixLoop) MHHadronness::Fill
|
---|
1081 | // will take the values from the matrix instead of the containers.
|
---|
1082 | //
|
---|
1083 | // It takes fSkipHist* into account!
|
---|
1084 | //
|
---|
1085 | void MHAlpha::InitMapping(MHMatrix *mat, Int_t type)
|
---|
1086 | {
|
---|
1087 | if (fMatrix)
|
---|
1088 | return;
|
---|
1089 |
|
---|
1090 | fMatrix = mat;
|
---|
1091 |
|
---|
1092 | fMap[0] = fMatrix->AddColumn(GetParameterRule());
|
---|
1093 | fMap[1] = -1;
|
---|
1094 | fMap[2] = -1;
|
---|
1095 | fMap[3] = -1;
|
---|
1096 | fMap[4] = -1;
|
---|
1097 |
|
---|
1098 | if (!fSkipHistEnergy)
|
---|
1099 | {
|
---|
1100 | fMap[1] = type==0 ? fMatrix->AddColumn("MEnergyEst.fVal") : -1;
|
---|
1101 | fMap[2] = type==0 ? -1 : fMatrix->AddColumn("MHillas.fSize");
|
---|
1102 | }
|
---|
1103 |
|
---|
1104 | if (!fSkipHistTheta)
|
---|
1105 | fMap[3] = fMatrix->AddColumn("MPointingPos.fZd");
|
---|
1106 |
|
---|
1107 | // if (!fSkipHistTime)
|
---|
1108 | // fMap[4] = fMatrix->AddColumn("MTime.GetAxisTime");
|
---|
1109 | }
|
---|
1110 |
|
---|
1111 | void MHAlpha::StopMapping()
|
---|
1112 | {
|
---|
1113 | fMatrix = NULL;
|
---|
1114 | }
|
---|
1115 |
|
---|
1116 | void MHAlpha::ApplyScaling()
|
---|
1117 | {
|
---|
1118 | if (!fOffData)
|
---|
1119 | return;
|
---|
1120 |
|
---|
1121 | fFit.ApplyScaling(fHist, *const_cast<TH3D*>(fOffData));
|
---|
1122 | }
|
---|
1123 |
|
---|
1124 | Int_t MHAlpha::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
---|
1125 | {
|
---|
1126 | Bool_t rc = kFALSE;
|
---|
1127 | if (IsEnvDefined(env, prefix, "NumTimeBins", print))
|
---|
1128 | {
|
---|
1129 | SetNumTimeBins(GetEnvValue(env, prefix, "NumTimeBins", fNumTimeBins));
|
---|
1130 | rc = kTRUE;
|
---|
1131 | }
|
---|
1132 | if (IsEnvDefined(env, prefix, "ForceUsingSize", print))
|
---|
1133 | {
|
---|
1134 | fForceUsingSize = GetEnvValue(env, prefix, "ForceUsingSize", fForceUsingSize);
|
---|
1135 | rc = kTRUE;
|
---|
1136 | }
|
---|
1137 | return rc;
|
---|
1138 | }
|
---|
1139 |
|
---|
1140 | Int_t MHAlpha::DistancetoPrimitive(Int_t px, Int_t py)
|
---|
1141 | {
|
---|
1142 | // If pad has no subpad return (we are in one of the subpads)
|
---|
1143 | if (gPad->GetPad(1)==NULL)
|
---|
1144 | return 9999;
|
---|
1145 |
|
---|
1146 | // If pad has a subpad we are in the main pad. Check its value.
|
---|
1147 | return gPad->GetPad(1)->DistancetoPrimitive(px,py)==0 ? 0 : 9999;
|
---|
1148 | }
|
---|
1149 |
|
---|
1150 | void MHAlpha::RecursiveRemove(TObject *obj)
|
---|
1151 | {
|
---|
1152 | if (obj==fOffData)
|
---|
1153 | fOffData = 0;
|
---|
1154 | }
|
---|