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

Last change on this file since 8590 was 8265, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 12.8 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-2006
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 // Do not use totalarea inside the binomial error calculation:
133 // it is not a weight.
134 fHEnergy.Divide(hsel, hall, 1, 1, "b");
135#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
136 MH::SetBinomialErrors(fHEnergy, *hsel, *hall);
137#endif
138 fHEnergy.Scale(totalarea);
139
140 delete hsel;
141 delete hall;
142}
143
144
145Bool_t MHCollectionArea::SetupFill(const MParList *pl)
146{
147 fHistSel.Reset();
148 if (!fIsExtern)
149 fHistAll.Reset();
150
151 fMcEvt = (MMcEvt*)pl->FindObject("MMcEvt");
152 if (!fMcEvt)
153 {
154 *fLog << err << "MMcEvt not found... abort." << endl;
155 return kFALSE;
156 }
157 /*
158 fEnergy = (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
159 if (!fEnergy)
160 {
161 *fLog << err << "MEnergyEst [MParameterD] not found... abort." << endl;
162 return kFALSE;
163 }
164 */
165 MBinning binst, binse;
166 binst.SetEdges(fHistAll, 'x');
167 binse.SetEdges(fHistAll, 'y');
168
169 //if (!fIsExtern)
170 {
171 MBinning *bins = (MBinning*)pl->FindObject("BinningTheta", "MBinning");
172 if (bins)
173 binst.SetEdges(*bins);
174
175 bins = (MBinning*)pl->FindObject("BinningEnergyEst", "MBinning");
176 if (bins)
177 binse.SetEdges(*bins);
178 }
179
180 binse.Apply(fHEnergy);
181
182 MH::SetBinning(&fHistSel, &binst, &binse);
183 MH::SetBinning(&fHistAll, &binst, &binse);
184
185 fMcAreaRadius = -1;
186
187 return kTRUE;
188}
189
190Bool_t MHCollectionArea::ReInit(MParList *plist)
191{
192 MMcRunHeader *runheader = (MMcRunHeader*)plist->FindObject("MMcRunHeader");
193 if (!runheader)
194 {
195 *fLog << err << "MMcRunHeader not found... abort." << endl;
196 return kFALSE;
197 }
198
199 // FIXME: Does this need some weighting with the number of produced events?
200 if (runheader->GetImpactMax()>fMcAreaRadius*100)
201 {
202 fMcAreaRadius = 0.01*runheader->GetImpactMax(); // cm->m
203 *fLog << inf << "Maximum simulated impact: " << fMcAreaRadius << "m" << endl;
204 }
205
206 if (fCorsikaVersion!=0 && fCorsikaVersion!=runheader->GetCorsikaVersion())
207 *fLog << warn << "Warning - Read files have different Corsika versions..." << endl;
208 fCorsikaVersion = runheader->GetCorsikaVersion();
209
210 if (fIsExtern)
211 return kTRUE;
212
213 fTotalNumSimulatedShowers += runheader->GetNumSimulatedShowers();
214 *fLog << inf << "Total Number of Simulated showers: " << fTotalNumSimulatedShowers << endl;
215
216 fAllEvtsTriggered |= runheader->GetAllEvtsTriggered();
217 *fLog << inf << "Only triggered events avail: " << (fAllEvtsTriggered?"yes":"no") << endl;
218
219 MMcCorsikaRunHeader *crh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
220 if (!crh)
221 {
222 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
223 return kFALSE;
224 }
225
226 //
227 // Calculate approximately the original number of events in each
228 // energy bin:
229 //
230 const Float_t emin = crh->GetELowLim();
231 const Float_t emax = crh->GetEUppLim();
232 const Float_t expo = 1 + crh->GetSlopeSpec();
233 const Float_t k = runheader->GetNumSimulatedShowers() /
234 (pow(emax,expo) - pow(emin,expo));
235
236 const Int_t nbiny = fHistAll.GetNbinsY();
237
238 TAxis &axe = *fHistAll.GetYaxis();
239 for (Int_t i = 1; i <= nbiny; i++)
240 {
241 const Float_t e1 = axe.GetBinLowEdge(i);
242 const Float_t e2 = axe.GetBinLowEdge(i+1);
243
244 if (e1 < emin || e2 > emax)
245 continue;
246
247 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
248 //
249 // We fill the i-th energy bin, with the total number of events
250 // Second argument of Fill would be impact parameter of each
251 // event, but we don't really need it for the collection area,
252 // so we just put a dummy value (1.)
253 //
254
255 const Float_t energy = (e1+e2)/2.;
256 fHistAll.Fill(20, energy, events);
257 // you have MMcRunHeader.fShowerThetaMin and MMcRunHeader.fShowerThetaMax
258 }
259
260 return kTRUE;
261}
262
263void MHCollectionArea::Paint(Option_t *option)
264{
265 if (TString(option)=="paint3")
266 {
267 /*
268 TH1 *h = dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
269 if (h)
270 {
271 const TString txt = Form("N/N_{0}=%.2f",
272 GetCollectionAreaEff(),
273 GetCollectionAreaAbs(), fMcAreaRadius);
274
275 TLatex text(0.31, 0.95, txt);
276 text.SetBit(TLatex::kTextNDC);
277 text.SetTextSize(0.04);
278 text.Paint();*/
279 return;
280 }
281 if (TString(option)=="paint4")
282 {
283 //const TString txt = Form("A_{eff}=%.0fm^{2} A_{abs}=%.0fm^{2} r=%.0fm",
284 // GetCollectionAreaEff(),
285 // GetCollectionAreaAbs(), fMcAreaRadius);
286 const TString txt = Form("A_{abs}=%.0fm^{2} r=%.0fm",
287 GetCollectionAreaAbs(), fMcAreaRadius);
288
289 TLatex text(0.31, 0.95, txt);
290 text.SetBit(TLatex::kTextNDC);
291 text.SetTextSize(0.04);
292 text.Paint();
293 return;
294 }
295
296 TVirtualPad *pad;
297
298 TPaveStats *st=0;
299 for (int x=0; x<4; x++)
300 {
301 pad=gPad->GetPad(x+1);
302 if (!pad || !(st = (TPaveStats*)pad->GetPrimitive("stats")))
303 continue;
304
305 if (st->GetOptStat()==11)
306 continue;
307
308 const Double_t y1 = st->GetY1NDC();
309 const Double_t y2 = st->GetY2NDC();
310 const Double_t x1 = st->GetX1NDC();
311 const Double_t x2 = st->GetX2NDC();
312
313 st->SetY1NDC((y2-y1)/3+y1);
314 st->SetX1NDC((x2-x1)/3+x1);
315 st->SetOptStat(11);
316 }
317
318 pad = gPad;
319
320 TH1 *h1=0, *h2=0;
321
322 pad->cd(1);
323 if (gPad->FindObject("ProjSelX"))
324 fHistSel.ProjectionX("ProjSelX", -1, 9999, "E");
325
326 pad->cd(2);
327 if (gPad->FindObject("ProjAllY"))
328 h1=fHistAll.ProjectionY("ProjAllY", -1, 9999, "E");
329 if (gPad->FindObject("ProjSelY"))
330 h2=fHistSel.ProjectionY("ProjSelY", -1, 9999, "E");
331
332 if (h1 && h1->GetMaximum()>0)
333 {
334 gPad->SetLogx();
335 gPad->SetLogy();
336 }
337
338 pad->cd(3);
339 TH1 *h=dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
340 if (h1 && h2 && h)
341 {
342 h->Divide(h2, h1, 1, 1, "b");
343#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
344 MH::SetBinomialErrors(*h, *h2, *h1);
345#endif
346 h->SetMinimum(0);
347 }
348
349 pad->cd(4);
350 CalcEfficiency();
351 if (fHEnergy.GetMaximum()>0)
352 {
353 gPad->SetLogx();
354 gPad->SetLogy();
355 }
356}
357
358void MHCollectionArea::Draw(Option_t *option)
359{
360 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
361
362 // Do the projection before painting the histograms into
363 // the individual pads
364 AppendPad();
365
366 pad->SetBorderMode(0);
367 pad->Divide(2,2);
368
369 TH1 *h=0, *h1=0, *h2=0;
370
371 if (fHistSel.GetNbinsX()>1)
372 {
373 pad->cd(1);
374 gPad->SetBorderMode(0);
375 gPad->SetGridx();
376 gPad->SetGridy();
377 /*
378 h = fHistAll.ProjectionX("ProjAllX", -1, 9999, "E");
379 h->SetXTitle("\\Theta [\\circ]");
380 h->SetDirectory(NULL);
381 h->SetLineColor(kGreen);
382 h->SetBit(kCanDelete);
383 h->Draw();
384 */
385 h = fHistSel.ProjectionX("ProjSelX", -1, 9999, "E");
386 h->SetXTitle("\\Theta [\\circ]");
387 h->SetDirectory(NULL);
388 h->SetLineColor(kRed);
389 h->SetBit(kCanDelete);
390 h->Draw("hist"/*"same"*/);
391 }
392 else
393 delete pad->GetPad(1);
394
395 if (fHistSel.GetNbinsY()>1)
396 {
397 pad->cd(2);
398 gPad->SetBorderMode(0);
399 gPad->SetGridx();
400 gPad->SetGridy();
401
402 h1 = fHistAll.ProjectionY("ProjAllY", -1, 9999, "E");
403 h1->SetDirectory(NULL);
404 h1->SetLineColor(kGreen);
405 h1->SetXTitle("E [GeV]");
406 h1->SetBit(kCanDelete);
407 h1->Draw();
408
409 h2 = fHistSel.ProjectionY("ProjSelY", -1, 9999, "E");
410 h2->SetDirectory(NULL);
411 h2->SetLineColor(kRed);
412 h2->SetBit(kCanDelete);
413 h2->Draw("same");
414 }
415 else
416 delete pad->GetPad(2);
417
418 if (h1 && h2)
419 {
420 pad->cd(3);
421 gPad->SetBorderMode(0);
422 gPad->SetGridx();
423 gPad->SetGridy();
424 gPad->SetLogx();
425 h = h2->DrawCopy();
426 h->Divide(h2, h1, 1, 1, "b");
427#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
428 MH::SetBinomialErrors(*h, *h2, *h1);
429#endif
430 h->SetNameTitle("Efficiency", "Combined cut and trigger efficiency");
431 h->SetDirectory(NULL);
432 AppendPad("paint3");
433 }
434 else
435 delete pad->GetPad(4);
436
437 if (fHEnergy.GetNbinsX()>1)
438 {
439 pad->cd(4);
440 gPad->SetBorderMode(0);
441 gPad->SetGridx();
442 gPad->SetGridy();
443 fHEnergy.Draw();
444 AppendPad("paint4");
445 }
446 else
447 delete pad->GetPad(4);
448}
449
450Bool_t MHCollectionArea::Fill(const MParContainer *par, const Stat_t weight)
451{
452 const Double_t energy = fMcEvt->GetEnergy();
453 const Double_t theta = fMcEvt->GetTelescopeTheta()*TMath::RadToDeg();
454
455 fHistSel.Fill(theta, energy, weight);
456
457 return kTRUE;
458}
Note: See TracBrowser for help on using the repository browser.