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): A. Moralejo 3/2003 <mailto:moralejo@pd.infn.it>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2003
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | // //
|
---|
27 | // MHMcCT1CollectionArea //
|
---|
28 | // //
|
---|
29 | //////////////////////////////////////////////////////////////////////////////
|
---|
30 |
|
---|
31 | #include "MHMcCT1CollectionArea.h"
|
---|
32 |
|
---|
33 | #include <TH2.h>
|
---|
34 | #include <TCanvas.h>
|
---|
35 |
|
---|
36 | #include "MMcEvt.hxx"
|
---|
37 | #include "MH.h"
|
---|
38 | #include "MBinning.h"
|
---|
39 | #include "MParList.h"
|
---|
40 | #include "MLog.h"
|
---|
41 | #include "MLogManip.h"
|
---|
42 |
|
---|
43 |
|
---|
44 | ClassImp(MHMcCT1CollectionArea);
|
---|
45 |
|
---|
46 | // --------------------------------------------------------------------------
|
---|
47 | //
|
---|
48 | // Creates the three necessary histograms:
|
---|
49 | // - selected showers (input)
|
---|
50 | // - all showers (input)
|
---|
51 | // - collection area (result)
|
---|
52 | //
|
---|
53 | MHMcCT1CollectionArea::MHMcCT1CollectionArea(const char *name, const char *title)
|
---|
54 | {
|
---|
55 | //
|
---|
56 | // nbins, minEnergy, maxEnergy defaults:
|
---|
57 | // we set the energy range from 100 Gev to 30000 GeV (in log, 3.5 orders
|
---|
58 | // of magnitude) and for each order we take 10 subdivisions --> 35 xbins
|
---|
59 | // we set the theta range from 12.5 to 48 deg, with 6 bins (the latter
|
---|
60 | // choice has been done to make the bin centers as close as possible to
|
---|
61 | // the actual zenith angles in the CT1 MC sample).
|
---|
62 | //
|
---|
63 |
|
---|
64 | fName = name ? name : "MHMcCT1CollectionArea";
|
---|
65 | fTitle = title ? title : "Collection Area vs. log10 Energy";
|
---|
66 |
|
---|
67 | fHistAll = new TH2D;
|
---|
68 | fHistSel = new TH2D;
|
---|
69 | fHistCol = new TH2D;
|
---|
70 |
|
---|
71 | fHistCol->SetName(fName);
|
---|
72 | fHistAll->SetName("AllEvents");
|
---|
73 | fHistSel->SetName("SelectedEvents");
|
---|
74 |
|
---|
75 | fHistCol->SetTitle(fTitle);
|
---|
76 | fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
|
---|
77 | fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
|
---|
78 |
|
---|
79 | fHistAll->SetDirectory(NULL);
|
---|
80 | fHistSel->SetDirectory(NULL);
|
---|
81 | fHistCol->SetDirectory(NULL);
|
---|
82 |
|
---|
83 | fHistAll->SetXTitle("log10 E [GeV]");
|
---|
84 | fHistAll->SetYTitle("theta [deg]");
|
---|
85 | fHistAll->SetZTitle("N");
|
---|
86 |
|
---|
87 | fHistSel->SetXTitle("log10 E [GeV]");
|
---|
88 | fHistSel->SetYTitle("theta [deg]");
|
---|
89 | fHistSel->SetZTitle("N");
|
---|
90 |
|
---|
91 | fHistCol->SetXTitle("log10 E [GeV]");
|
---|
92 | fHistCol->SetYTitle("theta [deg]");
|
---|
93 | fHistCol->SetZTitle("A [m^{2}]");
|
---|
94 | }
|
---|
95 |
|
---|
96 | // --------------------------------------------------------------------------
|
---|
97 | //
|
---|
98 | // Delete the three histograms
|
---|
99 | //
|
---|
100 | MHMcCT1CollectionArea::~MHMcCT1CollectionArea()
|
---|
101 | {
|
---|
102 | delete fHistAll;
|
---|
103 | delete fHistSel;
|
---|
104 | delete fHistCol;
|
---|
105 | }
|
---|
106 |
|
---|
107 | // --------------------------------------------------------------------------
|
---|
108 | //
|
---|
109 | // Set the binnings and prepare the filling of the histograms
|
---|
110 | //
|
---|
111 | Bool_t MHMcCT1CollectionArea::SetupFill(const MParList *plist)
|
---|
112 | {
|
---|
113 | const MBinning* binsenergy = (MBinning*)plist->FindObject("BinningE");
|
---|
114 | const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
|
---|
115 |
|
---|
116 | if (!binsenergy || !binstheta)
|
---|
117 | {
|
---|
118 | *fLog << err << dbginf << "At least one MBinning not found... aborting.";
|
---|
119 | *fLog << endl;
|
---|
120 | return kFALSE;
|
---|
121 | }
|
---|
122 |
|
---|
123 | SetBinning(fHistAll, binsenergy, binstheta);
|
---|
124 | SetBinning(fHistSel, binsenergy, binstheta);
|
---|
125 | SetBinning(fHistCol, binsenergy, binstheta);
|
---|
126 |
|
---|
127 | fHistAll->Sumw2();
|
---|
128 | fHistSel->Sumw2();
|
---|
129 | fHistCol->Sumw2();
|
---|
130 |
|
---|
131 | return kTRUE;
|
---|
132 | }
|
---|
133 |
|
---|
134 |
|
---|
135 | // --------------------------------------------------------------------------
|
---|
136 | //
|
---|
137 | // Fill data into the histogram which contains the selected showers
|
---|
138 | //
|
---|
139 | Bool_t MHMcCT1CollectionArea::Fill(const MParContainer *par, const Stat_t w)
|
---|
140 | {
|
---|
141 | MMcEvt &mcevt = *(MMcEvt*)par;
|
---|
142 |
|
---|
143 | fHistSel->Fill(log10(mcevt.GetEnergy()), kRad2Deg*mcevt.GetTelescopeTheta(), w);
|
---|
144 | return kTRUE;
|
---|
145 | }
|
---|
146 |
|
---|
147 | // --------------------------------------------------------------------------
|
---|
148 | //
|
---|
149 | // Draw the histogram with all showers
|
---|
150 | //
|
---|
151 | void MHMcCT1CollectionArea::DrawAll(Option_t* option)
|
---|
152 | {
|
---|
153 | if (!gPad)
|
---|
154 | MH::MakeDefCanvas(fHistAll);
|
---|
155 |
|
---|
156 | fHistAll->Draw(option);
|
---|
157 |
|
---|
158 | gPad->Modified();
|
---|
159 | gPad->Update();
|
---|
160 | }
|
---|
161 |
|
---|
162 | // --------------------------------------------------------------------------
|
---|
163 | //
|
---|
164 | // Draw the histogram with the selected showers only.
|
---|
165 | //
|
---|
166 | void MHMcCT1CollectionArea::DrawSel(Option_t* option)
|
---|
167 | {
|
---|
168 | if (!gPad)
|
---|
169 | MH::MakeDefCanvas(fHistSel);
|
---|
170 |
|
---|
171 | fHistSel->Draw(option);
|
---|
172 |
|
---|
173 | gPad->Modified();
|
---|
174 | gPad->Update();
|
---|
175 | }
|
---|
176 |
|
---|
177 | // --------------------------------------------------------------------------
|
---|
178 | //
|
---|
179 | // Creates a new canvas and draws the histogram into it.
|
---|
180 | // Be careful: The histogram belongs to this object and won't get deleted
|
---|
181 | // together with the canvas.
|
---|
182 | //
|
---|
183 | TObject *MHMcCT1CollectionArea::DrawClone(Option_t* option) const
|
---|
184 | {
|
---|
185 | TCanvas &c = *MakeDefCanvas("CollArea", "Collection area plots", 600, 600);
|
---|
186 | c.Divide(2,2);
|
---|
187 |
|
---|
188 | //
|
---|
189 | // This is necessary to get the expected behaviour of DrawClone
|
---|
190 | //
|
---|
191 | gROOT->SetSelectedPad(NULL);
|
---|
192 |
|
---|
193 | c.cd(1);
|
---|
194 | fHistCol->SetDirectory(NULL);
|
---|
195 | fHistCol->DrawCopy(option);
|
---|
196 |
|
---|
197 | c.cd(2);
|
---|
198 | fHistSel->SetDirectory(NULL);
|
---|
199 | fHistSel->DrawCopy(option);
|
---|
200 |
|
---|
201 | c.cd(3);
|
---|
202 | fHistAll->SetDirectory(NULL);
|
---|
203 | fHistAll->DrawCopy(option);
|
---|
204 |
|
---|
205 |
|
---|
206 | c.Modified();
|
---|
207 | c.Update();
|
---|
208 |
|
---|
209 | return &c;
|
---|
210 | }
|
---|
211 |
|
---|
212 | void MHMcCT1CollectionArea::Draw(Option_t* option)
|
---|
213 | {
|
---|
214 | if (!gPad)
|
---|
215 | MH::MakeDefCanvas(fHistCol);
|
---|
216 |
|
---|
217 | fHistCol->Draw(option);
|
---|
218 |
|
---|
219 | gPad->Modified();
|
---|
220 | gPad->Update();
|
---|
221 | }
|
---|
222 |
|
---|
223 | //
|
---|
224 | // Calculate the Efficiency (collection area) for the CT1 MC sample
|
---|
225 | // and set the 'ReadyToSave' flag
|
---|
226 | //
|
---|
227 | void MHMcCT1CollectionArea::CalcEfficiency()
|
---|
228 | {
|
---|
229 | //
|
---|
230 | // Here we estimate the total number of showers in each energy bin
|
---|
231 | // from the known the energy range and spectral index of the generated
|
---|
232 | // showers. This procedure is intended for the CT1 MC files. The total
|
---|
233 | // number of generated events, collection area, spectral index etc will be
|
---|
234 | // set here by hand, so make sure that the MC sample you are using is the
|
---|
235 | // right one (check all these quantities in your files and compare with
|
---|
236 | // what is written below. In some theta bins, there are two different
|
---|
237 | // productions, with different energy limits but with the same spectral
|
---|
238 | // slope. We account for this when calculating the original number of
|
---|
239 | // events in each energy bin.
|
---|
240 | //
|
---|
241 |
|
---|
242 | for (Int_t thetabin = 1; thetabin <= fHistAll->GetNbinsY(); thetabin++)
|
---|
243 | {
|
---|
244 | // This theta is not exactly the one of the MC events, just about
|
---|
245 | // the same:
|
---|
246 | Float_t theta = fHistAll->GetYaxis()->GetBinCenter(thetabin);
|
---|
247 |
|
---|
248 | Float_t emin1, emax1, emin2, emax2;
|
---|
249 | Float_t index, expo, k1, k2;
|
---|
250 | Float_t numevts1, numevts2;
|
---|
251 | Float_t r1, r2; // Impact parameter range (on ground).
|
---|
252 |
|
---|
253 | emin1 = 0; emax1 = 0; emin2 = 0; emax2 = 0;
|
---|
254 | expo = 0.; k1 = 0.; k2 = 0.; r1 = 0.; r2 = 0.;
|
---|
255 | numevts1 = 0; numevts2 = 0;
|
---|
256 |
|
---|
257 | if (theta > 14 && theta < 16) // 15 deg
|
---|
258 | {
|
---|
259 | r1 = 0.;
|
---|
260 | r2 = 250.; //meters
|
---|
261 | emin1 = 300.;
|
---|
262 | emax1 = 400.; // Energies in GeV.
|
---|
263 | emin2 = 400.;
|
---|
264 | emax2 = 30000.;
|
---|
265 | numevts1 = 4000.;
|
---|
266 | numevts2 = 25740.;
|
---|
267 | }
|
---|
268 | else if (theta > 20 && theta < 21) // 20.5 deg
|
---|
269 | {
|
---|
270 | r1 = 0.;
|
---|
271 | r2 = 263.; //meters
|
---|
272 | emin1 = 300.;
|
---|
273 | emax1 = 400.; // Energies in GeV.
|
---|
274 | emin2 = 400.;
|
---|
275 | emax2 = 30000.;
|
---|
276 | numevts1 = 6611.;
|
---|
277 | numevts2 = 24448.;
|
---|
278 | }
|
---|
279 | else if (theta > 26 && theta < 27) // 26.5 degrees
|
---|
280 | {
|
---|
281 | r1 = 0.;
|
---|
282 | r2 = 290.; //meters
|
---|
283 | emin1 = 300.;
|
---|
284 | emax1 = 400.; // Energies in GeV.
|
---|
285 | emax2 = emax1; emin2 = 400.;
|
---|
286 | emax2 = 30000.;
|
---|
287 | numevts1 = 4000.;
|
---|
288 | numevts2 = 26316.;
|
---|
289 | }
|
---|
290 | else if (theta > 32 && theta < 33) // 32.5 degrees
|
---|
291 | {
|
---|
292 | r1 = 0.;
|
---|
293 | r2 = 350.; //meters
|
---|
294 | emin1 = 300.;
|
---|
295 | emax1 = 30000.; // Energies in GeV.
|
---|
296 | emax2 = emax1;
|
---|
297 | numevts1 = 33646.;
|
---|
298 | }
|
---|
299 | else if (theta > 38 && theta < 39) // 38.75 degrees
|
---|
300 | {
|
---|
301 | r1 = 0.;
|
---|
302 | r2 = 380.; //meters
|
---|
303 | emin1 = 300.;
|
---|
304 | emax1 = 30000.; // Energies in GeV.
|
---|
305 | emax2 = emax1;
|
---|
306 | numevts1 = 38415.;
|
---|
307 | }
|
---|
308 | else if (theta > 44 && theta < 46) // 45 degrees
|
---|
309 | {
|
---|
310 | r1 = 0.;
|
---|
311 | r2 = 565.; //meters
|
---|
312 | emin1 = 300.;
|
---|
313 | emax1 = 50000.; // Energies in GeV.
|
---|
314 | emax2 = emax1;
|
---|
315 | numevts1 = 30197.;
|
---|
316 | }
|
---|
317 |
|
---|
318 | index = 1.5; // Differential spectral Index.
|
---|
319 | expo = 1.-index;
|
---|
320 | k1 = numevts1 / (pow(emax1,expo) - pow(emin1,expo));
|
---|
321 | k2 = numevts2 / (pow(emax2,expo) - pow(emin2,expo));
|
---|
322 |
|
---|
323 | for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
|
---|
324 | {
|
---|
325 | const Float_t e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
|
---|
326 | const Float_t e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
|
---|
327 |
|
---|
328 | if (e1 < emin1 || e2 > emax2)
|
---|
329 | continue;
|
---|
330 |
|
---|
331 | Float_t events;
|
---|
332 |
|
---|
333 | if (e2 <= emax1)
|
---|
334 | events = k1 * (pow(e2, expo) - pow(e1, expo));
|
---|
335 | else if (e1 >= emin2)
|
---|
336 | events = k2 * (pow(e2, expo) - pow(e1, expo));
|
---|
337 | else
|
---|
338 | events =
|
---|
339 | k1 * (pow(emax1, expo) - pow(e1, expo))+
|
---|
340 | k2 * (pow(e2, expo) - pow(emin2, expo));
|
---|
341 |
|
---|
342 | fHistAll->SetBinContent(i, thetabin, events);
|
---|
343 | fHistAll->SetBinError(i, thetabin, sqrt(events));
|
---|
344 | }
|
---|
345 |
|
---|
346 | // -----------------------------------------------------------
|
---|
347 |
|
---|
348 | const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
|
---|
349 |
|
---|
350 | for (Int_t ix = 1; ix <= fHistAll->GetNbinsX(); ix++)
|
---|
351 | {
|
---|
352 | const Float_t Na = fHistAll->GetBinContent(ix,thetabin);
|
---|
353 |
|
---|
354 | if (Na <= 0)
|
---|
355 | {
|
---|
356 | //
|
---|
357 | // If energy is large, this case means that no or very few events
|
---|
358 | // were generated at this energy bin. In this case we assign it
|
---|
359 | // the effective area of the bin below it in energy. If energy is
|
---|
360 | // below 1E4, it means that no events triggered -> eff area = 0
|
---|
361 | //
|
---|
362 |
|
---|
363 | if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
|
---|
364 | {
|
---|
365 | fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
|
---|
366 | fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
|
---|
367 | }
|
---|
368 | continue;
|
---|
369 | }
|
---|
370 |
|
---|
371 | const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
|
---|
372 |
|
---|
373 | // Since Na is an estimate of the total number of showers generated
|
---|
374 | // in the energy bin, it may happen that Ns (triggered showers) is
|
---|
375 | // larger than Na. In that case, the bin is skipped:
|
---|
376 |
|
---|
377 | if (Na < Ns)
|
---|
378 | continue;
|
---|
379 |
|
---|
380 | const Double_t eff = Ns/Na;
|
---|
381 | const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
|
---|
382 |
|
---|
383 |
|
---|
384 | const Float_t area = dr * cos(theta*TMath::Pi()/180.);
|
---|
385 |
|
---|
386 | fHistCol->SetBinContent(ix, thetabin, eff*area);
|
---|
387 | fHistCol->SetBinError(ix, thetabin, efferr*area);
|
---|
388 |
|
---|
389 | }
|
---|
390 | }
|
---|
391 |
|
---|
392 | SetReadyToSave();
|
---|
393 | }
|
---|