source: trunk/MagicSoft/Mars/mhflux/MHCollectionArea.cc@ 8049

Last change on this file since 8049 was 8049, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 12.7 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): Thomas Bretz 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Harald Kornmayer 1/2001
20!
21! Copyright: MAGIC Software Development, 2000-2002
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27//
28// MHCollectionArea
29//
30// Class Version 2:
31// ----------------
32// + added //! to fMcEvt which was missing before
33// - removed obsolete data member fEnergy
34//
35//////////////////////////////////////////////////////////////////////////////
36#include "MHCollectionArea.h"
37
38#include <TLatex.h>
39#include <TCanvas.h>
40#include <TPaveStats.h>
41
42#include "MLog.h"
43#include "MLogManip.h"
44
45#include "MBinning.h"
46
47#include "MMcEvt.hxx"
48#include "MMcRunHeader.hxx"
49#include "MMcCorsikaRunHeader.h"
50
51#include "MParList.h"
52#include "MParameters.h"
53
54ClassImp(MHCollectionArea);
55
56using namespace std;
57
58// --------------------------------------------------------------------------
59//
60// Creates the three necessary histograms:
61// - selected showers (input)
62// - all showers (input)
63// - collection area (result)
64//
65MHCollectionArea::MHCollectionArea(const char *name, const char *title)
66 : fMcEvt(0), fMcAreaRadius(-1), fIsExtern(kFALSE)
67{
68 // initialize the histogram for the distribution r vs E
69 //
70 // we set the energy range from 2 Gev to 20000 GeV (in log 4 orders
71 // of magnitude) and for each order we take 25 subdivision --> 100 xbins
72 //
73 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
74 //
75 fName = name ? name : "MHCollectionArea";
76 fTitle = title ? title : "Collection Area vs. Energy/Theta";
77
78 fHistSel.SetName("SelEvts");
79 fHistSel.SetTitle("Number of Events after cuts");
80 fHistSel.SetXTitle("\\Theta [deg]");
81 fHistSel.SetYTitle("E [GeV]");
82 fHistSel.SetDirectory(NULL);
83 fHistSel.UseCurrentStyle();
84 fHistSel.SetLineColor(kBlue);
85
86 fHistAll.SetName("AllEvts");
87 fHistAll.SetTitle("Number of events produced");
88 fHistAll.SetXTitle("\\Theta [deg]");
89 fHistAll.SetYTitle("E_{mc} [GeV]");
90 fHistAll.SetDirectory(NULL);
91 fHistAll.UseCurrentStyle();
92
93 fHEnergy.SetName("CollEnergy");
94 fHEnergy.SetTitle("Collection Area");
95 fHEnergy.SetXTitle("E [GeV]");
96 fHEnergy.SetYTitle("A_{eff} [m^{2}]");
97 fHEnergy.SetDirectory(NULL);
98 fHEnergy.UseCurrentStyle();
99
100 MBinning binsa, binse, binst;
101 binse.SetEdgesLog(15, 10, 1000000);
102 binst.SetEdgesASin(67, -0.005, 0.665);
103
104 binse.Apply(fHEnergy);
105
106 MH::SetBinning(&fHistSel, &binst, &binse);
107 MH::SetBinning(&fHistAll, &binst, &binse);
108
109 // For some unknown reasons this must be called after
110 // the binning has been initialized at least once
111 fHistSel.Sumw2();
112 fHistAll.Sumw2();
113 fHEnergy.Sumw2();
114}
115
116// --------------------------------------------------------------------------
117//
118// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
119// flag
120//
121void MHCollectionArea::CalcEfficiency()
122{
123 TH1D *hsel = fHistSel.ProjectionY("Spy", -1, 9999, "E");;
124 TH1D *hall = fHistAll.ProjectionY("Apy", -1, 9999, "E");
125
126 //
127 // Impact parameter range.
128 //
129 const Float_t totalarea = GetCollectionAreaAbs();//TMath::Pi() * (r2*r2 - r1*r1);
130
131 // "b" option: calculate binomial errors
132 fHEnergy.Divide(hsel, hall, totalarea, 1, "b");
133#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
134 MH::SetBinomialErrors(fHEnergy, *hsel, *hall, totalarea, 1);
135#endif
136
137 delete hsel;
138 delete hall;
139}
140
141
142Bool_t MHCollectionArea::SetupFill(const MParList *pl)
143{
144 fHistSel.Reset();
145 if (!fIsExtern)
146 fHistAll.Reset();
147
148 fMcEvt = (MMcEvt*)pl->FindObject("MMcEvt");
149 if (!fMcEvt)
150 {
151 *fLog << err << "MMcEvt not found... abort." << endl;
152 return kFALSE;
153 }
154 /*
155 fEnergy = (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
156 if (!fEnergy)
157 {
158 *fLog << err << "MEnergyEst [MParameterD] not found... abort." << endl;
159 return kFALSE;
160 }
161 */
162 MBinning binst, binse;
163 binst.SetEdges(fHistAll, 'x');
164 binse.SetEdges(fHistAll, 'y');
165
166 //if (!fIsExtern)
167 {
168 MBinning *bins = (MBinning*)pl->FindObject("BinningTheta", "MBinning");
169 if (bins)
170 binst.SetEdges(*bins);
171
172 bins = (MBinning*)pl->FindObject("BinningEnergyEst", "MBinning");
173 if (bins)
174 binse.SetEdges(*bins);
175 }
176
177 binse.Apply(fHEnergy);
178
179 MH::SetBinning(&fHistSel, &binst, &binse);
180 MH::SetBinning(&fHistAll, &binst, &binse);
181
182 fMcAreaRadius = -1;
183
184 return kTRUE;
185}
186
187Bool_t MHCollectionArea::ReInit(MParList *plist)
188{
189 MMcRunHeader *runheader = (MMcRunHeader*)plist->FindObject("MMcRunHeader");
190 if (!runheader)
191 {
192 *fLog << err << "MMcRunHeader not found... abort." << endl;
193 return kFALSE;
194 }
195
196 // FIXME: Does this need some weighting with the number of produced events?
197 if (runheader->GetImpactMax()>fMcAreaRadius*100)
198 {
199 fMcAreaRadius = 0.01*runheader->GetImpactMax(); // cm->m
200 *fLog << inf << "Maximum simulated impact: " << fMcAreaRadius << "m" << endl;
201 }
202
203 if (fCorsikaVersion!=0 && fCorsikaVersion!=runheader->GetCorsikaVersion())
204 *fLog << warn << "Warning - Read files have different Corsika versions..." << endl;
205 fCorsikaVersion = runheader->GetCorsikaVersion();
206
207 if (fIsExtern)
208 return kTRUE;
209
210 fTotalNumSimulatedShowers += runheader->GetNumSimulatedShowers();
211 *fLog << inf << "Total Number of Simulated showers: " << fTotalNumSimulatedShowers << endl;
212
213 fAllEvtsTriggered |= runheader->GetAllEvtsTriggered();
214 *fLog << inf << "Only triggered events avail: " << (fAllEvtsTriggered?"yes":"no") << endl;
215
216 MMcCorsikaRunHeader *crh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
217 if (!crh)
218 {
219 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
220 return kFALSE;
221 }
222
223 //
224 // Calculate approximately the original number of events in each
225 // energy bin:
226 //
227 const Float_t emin = crh->GetELowLim();
228 const Float_t emax = crh->GetEUppLim();
229 const Float_t expo = 1 + crh->GetSlopeSpec();
230 const Float_t k = runheader->GetNumSimulatedShowers() /
231 (pow(emax,expo) - pow(emin,expo));
232
233 const Int_t nbiny = fHistAll.GetNbinsY();
234
235 TAxis &axe = *fHistAll.GetYaxis();
236 for (Int_t i = 1; i <= nbiny; i++)
237 {
238 const Float_t e1 = axe.GetBinLowEdge(i);
239 const Float_t e2 = axe.GetBinLowEdge(i+1);
240
241 if (e1 < emin || e2 > emax)
242 continue;
243
244 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
245 //
246 // We fill the i-th energy bin, with the total number of events
247 // Second argument of Fill would be impact parameter of each
248 // event, but we don't really need it for the collection area,
249 // so we just put a dummy value (1.)
250 //
251
252 const Float_t energy = (e1+e2)/2.;
253 fHistAll.Fill(20, energy, events);
254 // you have MMcRunHeader.fShowerThetaMin and MMcRunHeader.fShowerThetaMax
255 }
256
257 return kTRUE;
258}
259
260void MHCollectionArea::Paint(Option_t *option)
261{
262 if (TString(option)=="paint3")
263 {
264 /*
265 TH1 *h = dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
266 if (h)
267 {
268 const TString txt = Form("N/N_{0}=%.2f",
269 GetCollectionAreaEff(),
270 GetCollectionAreaAbs(), fMcAreaRadius);
271
272 TLatex text(0.31, 0.95, txt);
273 text.SetBit(TLatex::kTextNDC);
274 text.SetTextSize(0.04);
275 text.Paint();*/
276 return;
277 }
278 if (TString(option)=="paint4")
279 {
280 //const TString txt = Form("A_{eff}=%.0fm^{2} A_{abs}=%.0fm^{2} r=%.0fm",
281 // GetCollectionAreaEff(),
282 // GetCollectionAreaAbs(), fMcAreaRadius);
283 const TString txt = Form("A_{abs}=%.0fm^{2} r=%.0fm",
284 GetCollectionAreaAbs(), fMcAreaRadius);
285
286 TLatex text(0.31, 0.95, txt);
287 text.SetBit(TLatex::kTextNDC);
288 text.SetTextSize(0.04);
289 text.Paint();
290 return;
291 }
292
293 TVirtualPad *pad;
294
295 TPaveStats *st=0;
296 for (int x=0; x<4; x++)
297 {
298 pad=gPad->GetPad(x+1);
299 if (!pad || !(st = (TPaveStats*)pad->GetPrimitive("stats")))
300 continue;
301
302 if (st->GetOptStat()==11)
303 continue;
304
305 const Double_t y1 = st->GetY1NDC();
306 const Double_t y2 = st->GetY2NDC();
307 const Double_t x1 = st->GetX1NDC();
308 const Double_t x2 = st->GetX2NDC();
309
310 st->SetY1NDC((y2-y1)/3+y1);
311 st->SetX1NDC((x2-x1)/3+x1);
312 st->SetOptStat(11);
313 }
314
315 pad = gPad;
316
317 TH1 *h1=0, *h2=0;
318
319 pad->cd(1);
320 if (gPad->FindObject("ProjSelX"))
321 fHistSel.ProjectionX("ProjSelX", -1, 9999, "E");
322
323 pad->cd(2);
324 if (gPad->FindObject("ProjAllY"))
325 h1=fHistAll.ProjectionY("ProjAllY", -1, 9999, "E");
326 if (gPad->FindObject("ProjSelY"))
327 h2=fHistSel.ProjectionY("ProjSelY", -1, 9999, "E");
328
329 if (h1 && h1->GetMaximum()>0)
330 {
331 gPad->SetLogx();
332 gPad->SetLogy();
333 }
334
335 pad->cd(3);
336 TH1 *h=dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
337 if (h1 && h2 && h)
338 {
339 h->Divide(h2, h1, 1, 1, "b");
340#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
341 MH::SetBinomialErrors(*h, *h2, *h1);
342#endif
343 h->SetMinimum(0);
344 }
345
346 pad->cd(4);
347 CalcEfficiency();
348 if (fHEnergy.GetMaximum()>0)
349 {
350 gPad->SetLogx();
351 gPad->SetLogy();
352 }
353}
354
355void MHCollectionArea::Draw(Option_t *option)
356{
357 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
358
359 // Do the projection before painting the histograms into
360 // the individual pads
361 AppendPad();
362
363 pad->SetBorderMode(0);
364 pad->Divide(2,2);
365
366 TH1 *h=0, *h1=0, *h2=0;
367
368 if (fHistSel.GetNbinsX()>1)
369 {
370 pad->cd(1);
371 gPad->SetBorderMode(0);
372 gPad->SetGridx();
373 gPad->SetGridy();
374 /*
375 h = fHistAll.ProjectionX("ProjAllX", -1, 9999, "E");
376 h->SetXTitle("\\Theta [\\circ]");
377 h->SetDirectory(NULL);
378 h->SetLineColor(kGreen);
379 h->SetBit(kCanDelete);
380 h->Draw();
381 */
382 h = fHistSel.ProjectionX("ProjSelX", -1, 9999, "E");
383 h->SetXTitle("\\Theta [\\circ]");
384 h->SetDirectory(NULL);
385 h->SetLineColor(kRed);
386 h->SetBit(kCanDelete);
387 h->Draw("hist"/*"same"*/);
388 }
389 else
390 delete pad->GetPad(1);
391
392 if (fHistSel.GetNbinsY()>1)
393 {
394 pad->cd(2);
395 gPad->SetBorderMode(0);
396 gPad->SetGridx();
397 gPad->SetGridy();
398
399 h1 = fHistAll.ProjectionY("ProjAllY", -1, 9999, "E");
400 h1->SetDirectory(NULL);
401 h1->SetLineColor(kGreen);
402 h1->SetXTitle("E [GeV]");
403 h1->SetBit(kCanDelete);
404 h1->Draw();
405
406 h2 = fHistSel.ProjectionY("ProjSelY", -1, 9999, "E");
407 h2->SetDirectory(NULL);
408 h2->SetLineColor(kRed);
409 h2->SetBit(kCanDelete);
410 h2->Draw("same");
411 }
412 else
413 delete pad->GetPad(2);
414
415 if (h1 && h2)
416 {
417 pad->cd(3);
418 gPad->SetBorderMode(0);
419 gPad->SetGridx();
420 gPad->SetGridy();
421 gPad->SetLogx();
422 h = h2->DrawCopy();
423 h->Divide(h2, h1, 1, 1, "b");
424#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
425 MH::SetBinomialErrors(*h, *h2, *h1);
426#endif
427 h->SetNameTitle("Efficiency", "Combined cut and trigger efficiency");
428 h->SetDirectory(NULL);
429 AppendPad("paint3");
430 }
431 else
432 delete pad->GetPad(4);
433
434 if (fHEnergy.GetNbinsX()>1)
435 {
436 pad->cd(4);
437 gPad->SetBorderMode(0);
438 gPad->SetGridx();
439 gPad->SetGridy();
440 fHEnergy.Draw();
441 AppendPad("paint4");
442 }
443 else
444 delete pad->GetPad(4);
445}
446
447Bool_t MHCollectionArea::Fill(const MParContainer *par, const Stat_t weight)
448{
449 const Double_t energy = fMcEvt->GetEnergy();
450 const Double_t theta = fMcEvt->GetTelescopeTheta()*TMath::RadToDeg();
451
452 fHistSel.Fill(theta, energy, weight);
453
454 return kTRUE;
455}
Note: See TracBrowser for help on using the repository browser.