source: trunk/MagicSoft/Mars/mhist/MHMcCT1CollectionArea.cc@ 1889

Last change on this file since 1889 was 1854, checked in by moralejo, 22 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 10.4 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): 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 "MH.h"
37#include "MBinning.h"
38
39ClassImp(MHMcCT1CollectionArea);
40
41// --------------------------------------------------------------------------
42//
43// Creates the three necessary histograms:
44// - selected showers (input)
45// - all showers (input)
46// - collection area (result)
47//
48MHMcCT1CollectionArea::MHMcCT1CollectionArea(const char *name, const char *title, Int_t nbins, Axis_t minEnergy, Axis_t maxEnergy)
49{
50 //
51 // nbins, minEnergy, maxEnergy defaults:
52 // we set the energy range from 100 Gev to 30000 GeV (in log, 3.5 orders
53 // of magnitude) and for each order we take 10 subdivisions --> 35 xbins
54 // we set the theta range from 12.5 to 48 deg, with 6 bins (the latter
55 // choice has been done to make the bin centers as close as possible to
56 // the actual zenith angles in the CT1 MC sample).
57 //
58
59 fName = name ? name : "MHMcCT1CollectionArea";
60 fTitle = title ? title : "Collection Area vs. log10 Energy";
61
62 fHistAll = new TH2D;
63 fHistSel = new TH2D;
64 fHistCol = new TH2D;
65
66 SetBins(nbins, minEnergy, maxEnergy);
67
68 fHistCol->SetName(fName);
69 fHistAll->SetName("AllEvents");
70 fHistSel->SetName("SelectedEvents");
71
72 fHistCol->SetTitle(fTitle);
73 fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
74 fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
75
76 fHistAll->SetDirectory(NULL);
77 fHistSel->SetDirectory(NULL);
78 fHistCol->SetDirectory(NULL);
79
80 fHistAll->SetXTitle("log10 E [GeV]");
81 fHistAll->SetYTitle("theta [deg]");
82 fHistAll->SetZTitle("N");
83
84 fHistSel->SetXTitle("log10 E [GeV]");
85 fHistSel->SetYTitle("theta [deg]");
86 fHistSel->SetZTitle("N");
87
88 fHistCol->SetXTitle("log10 E [GeV]");
89 fHistCol->SetYTitle("theta [deg]");
90 fHistCol->SetZTitle("A [m^{2}]");
91}
92
93// --------------------------------------------------------------------------
94//
95// Set binning of histograms. Binning of energy axis can be changed, that
96// of theta axis is always the same (to match the CT1 MC sample theta dist):
97//
98void MHMcCT1CollectionArea::SetBins(Int_t nbins, Axis_t minEnergy, Axis_t maxEnergy)
99{
100 MBinning binsx;
101 binsx.SetEdges(nbins, minEnergy, maxEnergy);
102
103 MBinning binsy;
104 const Double_t yedge[7] = {12.5, 17.5, 23.5, 29.5, 35.5, 42., 48.};
105 const TArrayD yed(7,yedge);
106 binsy.SetEdges(yed);
107
108 MH::SetBinning(fHistAll, &binsx, &binsy);
109 MH::SetBinning(fHistSel, &binsx, &binsy);
110 MH::SetBinning(fHistCol, &binsx, &binsy);
111}
112
113
114// --------------------------------------------------------------------------
115//
116// Delete the three histograms
117//
118MHMcCT1CollectionArea::~MHMcCT1CollectionArea()
119{
120 delete fHistAll;
121 delete fHistSel;
122 delete fHistCol;
123}
124
125// --------------------------------------------------------------------------
126//
127// Fill data into the histogram which contains the selected showers
128//
129void MHMcCT1CollectionArea::FillSel(Double_t energy, Double_t theta)
130{
131 fHistSel->Fill(log10(energy), theta);
132}
133
134// --------------------------------------------------------------------------
135//
136// Draw the histogram with all showers
137//
138void MHMcCT1CollectionArea::DrawAll(Option_t* option)
139{
140 if (!gPad)
141 MH::MakeDefCanvas(fHistAll);
142
143 fHistAll->Draw(option);
144
145 gPad->Modified();
146 gPad->Update();
147}
148
149// --------------------------------------------------------------------------
150//
151// Draw the histogram with the selected showers only.
152//
153void MHMcCT1CollectionArea::DrawSel(Option_t* option)
154{
155 if (!gPad)
156 MH::MakeDefCanvas(fHistSel);
157
158 fHistSel->Draw(option);
159
160 gPad->Modified();
161 gPad->Update();
162}
163
164// --------------------------------------------------------------------------
165//
166// Creates a new canvas and draws the histogram into it.
167// Be careful: The histogram belongs to this object and won't get deleted
168// together with the canvas.
169//
170TObject *MHMcCT1CollectionArea::DrawClone(Option_t* option) const
171{
172 TCanvas *c = MH::MakeDefCanvas(fHistCol);
173
174 //
175 // This is necessary to get the expected behaviour of DrawClone
176 //
177 gROOT->SetSelectedPad(NULL);
178
179 fHistCol->DrawCopy(option);
180
181 c->Modified();
182 c->Update();
183
184 return c;
185}
186
187void MHMcCT1CollectionArea::Draw(Option_t* option)
188{
189 if (!gPad)
190 MH::MakeDefCanvas(fHistCol);
191
192 fHistCol->Draw(option);
193
194 gPad->Modified();
195 gPad->Update();
196}
197
198//
199// Calculate the Efficiency (collection area) for the CT1 MC sample
200// and set the 'ReadyToSave' flag
201//
202void MHMcCT1CollectionArea::CalcEfficiency()
203{
204 //
205 // Here we estimate the total number of showers in each energy bin
206 // from the known the energy range and spectral index of the generated
207 // showers. This procedure is intended for the CT1 MC files. The total
208 // number of generated events, collection area, spectral index etc will be
209 // set here by hand, so make sure that the MC sample you are using is the
210 // right one (check all these quantities in your files and compare with
211 // is written below. In some theta bins, there are two different
212 // productions, with different energy limits but with the same spectral
213 // slope. We account for this when calculating the original number of
214 // events in each energy bin.
215 //
216
217 for (Int_t thetabin = 1; thetabin <= fHistAll->GetNbinsY(); thetabin++)
218 {
219 // This theta is not exactly the one of the MC events, just about
220 // the same:
221 Float_t theta = fHistAll->GetYaxis()->GetBinCenter(thetabin);
222
223 Float_t emin1, emax1, emin2, emax2;
224 Float_t index, expo, k1, k2;
225 Float_t numevts1, numevts2;
226 Float_t r1, r2; // Impact parameter range (on ground).
227
228 emin1 = 0; emax1 = 0; emin2 = 0; emax2 = 0;
229 expo = 0.; k1 = 0.; k2 = 0.; r1 = 0.; r2 = 0.;
230 numevts1 = 0; numevts2 = 0;
231
232 if (theta > 14 && theta < 16) // 15 deg
233 {
234 r1 = 0.;
235 r2 = 250.; //meters
236 emin1 = 300.;
237 emax1 = 400.; // Energies in GeV.
238 emin2 = 400.;
239 emax2 = 30000.;
240 numevts1 = 4000.;
241 numevts2 = 25740.;
242 }
243 else if (theta > 20 && theta < 21) // 20.5 deg
244 {
245 r1 = 0.;
246 r2 = 263.; //meters
247 emin1 = 300.;
248 emax1 = 400.; // Energies in GeV.
249 emin2 = 400.;
250 emax2 = 30000.;
251 numevts1 = 6611.;
252 numevts2 = 24448.;
253 }
254 else if (theta > 26 && theta < 27) // 26.5 degrees
255 {
256 r1 = 0.;
257 r2 = 290.; //meters
258 emin1 = 300.;
259 emax1 = 400.; // Energies in GeV.
260 emax2 = emax1; emin2 = 400.;
261 emax2 = 30000.;
262 numevts1 = 4000.;
263 numevts2 = 26316.;
264 }
265 else if (theta > 32 && theta < 33) // 32.5 degrees
266 {
267 r1 = 0.;
268 r2 = 350.; //meters
269 emin1 = 300.;
270 emax1 = 30000.; // Energies in GeV.
271 emax2 = emax1;
272 numevts1 = 33646.;
273 }
274 else if (theta > 38 && theta < 39) // 38.75 degrees
275 {
276 r1 = 0.;
277 r2 = 380.; //meters
278 emin1 = 300.;
279 emax1 = 30000.; // Energies in GeV.
280 emax2 = emax1;
281 numevts1 = 38415.;
282 }
283 else if (theta > 44 && theta < 46) // 45 degrees
284 {
285 r1 = 0.;
286 r2 = 565.; //meters
287 emin1 = 300.;
288 emax1 = 50000.; // Energies in GeV.
289 emax2 = emax1;
290 numevts1 = 30197.;
291 }
292
293 index = 1.5; // Differential spectral Index.
294 expo = 1.-index;
295 k1 = numevts1 / (pow(emax1,expo) - pow(emin1,expo));
296 k2 = numevts2 / (pow(emax2,expo) - pow(emin2,expo));
297
298 for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
299 {
300 const Float_t e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
301 const Float_t e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
302
303 if (e1 < emin1 || e2 > emax2)
304 continue;
305
306 Float_t events;
307
308 if (e2 <= emax1)
309 events = k1 * (pow(e2, expo) - pow(e1, expo));
310 else if (e1 >= emin2)
311 events = k2 * (pow(e2, expo) - pow(e1, expo));
312 else
313 events =
314 k1 * (pow(emax1, expo) - pow(e1, expo))+
315 k2 * (pow(e2, expo) - pow(emin2, expo));
316
317 fHistAll->SetBinContent(i, thetabin, events);
318 fHistAll->SetBinError(i, thetabin, sqrt(events));
319 }
320
321 // -----------------------------------------------------------
322
323 const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
324
325 for (Int_t ix = 1; ix <= fHistAll->GetNbinsX(); ix++)
326 {
327 const Float_t Na = fHistAll->GetBinContent(ix,thetabin);
328
329 if (Na <= 0)
330 {
331 //
332 // If energy is large, this case means that no or very few events
333 // were generated at this energy bin. In this case we assign it
334 // the effective area of the bin below it in energy. If energy is
335 // below 1E4, it means that no events triggered -> eff area = 0
336 //
337
338 if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
339 {
340 fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
341 fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
342 }
343 continue;
344 }
345
346 const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
347
348 // Since Na is an estimate of the total number of showers generated
349 // in the energy bin, it may happen that Ns (triggered showers) is
350 // larger than Na. In that case, the bin is skipped:
351
352 if (Na < Ns)
353 continue;
354
355 const Double_t eff = Ns/Na;
356 const Double_t err = sqrt((1.-eff)*Ns)/Na;
357
358
359 const Float_t area = dr * cos(theta*TMath::Pi()/180.);
360
361 fHistCol->SetBinContent(ix, thetabin, eff*area);
362 fHistCol->SetBinError(ix, thetabin, err*area);
363
364 }
365 }
366
367 SetReadyToSave();
368}
Note: See TracBrowser for help on using the repository browser.