source: trunk/MagicSoft/Mars/mhist/MHMcCollectionArea.cc@ 971

Last change on this file since 971 was 970, checked in by tbretz, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 6.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): Harald Kornmayer 1/2001 (harald@mppmu.mpg.de)
19! Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
20!
21! Copyright: MAGIC Software Development, 2000-2001
22!
23!
24\* ======================================================================== */
25
26#include "MHMcCollectionArea.h"
27
28#include <TH2.h>
29#include <TCanvas.h>
30
31ClassImp(MHMcCollectionArea);
32
33// --------------------------------------------------------------------------
34//
35// Creates the three necessary histograms:
36// - selected showers (input)
37// - all showers (input)
38// - collection area (result)
39//
40MHMcCollectionArea::MHMcCollectionArea(const char *name, const char *title)
41{
42 // initialize the histogram for the distribution r vs E
43 //
44 // we set the energy range from 1 Gev to 10000 GeV (in log 5 orders
45 // of magnitude) and for each order we take 10 subdivision --> 50 xbins
46 //
47 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
48
49
50 *fName = name ? name : "MHMcCollectionArea";
51 *fTitle = title ? title : "Data to Calculate Coll-Area";
52
53
54 fHistAll = new TH2D("AllEvents", "All showers - Radius vs log(E) distribution",
55 50, 0., 5.,
56 50, 0., 500.);
57
58 fHistAll->SetXTitle("log(E/GeV)");
59 fHistAll->SetYTitle("r/m");
60 fHistAll->SetZTitle("N");
61
62 fHistSel = new TH2D("SelectedEvents", "Selected showers - Radius vs log(E) distribution",
63 50, 0., 5.,
64 50, 0., 500.);
65
66 fHistSel->SetXTitle("log(E/GeV)");
67 fHistSel->SetYTitle("r/m");
68 fHistSel->SetYTitle("N");
69
70 fHistCol = new TH1D("CollectionArea", "Collection Area vs. log(E)",
71 50, 0., 5.);
72
73 fHistCol->SetXTitle("log(E/GeV)");
74 fHistCol->SetYTitle("A/m^{2}");
75}
76
77// --------------------------------------------------------------------------
78//
79// Delete the three histograms
80//
81MHMcCollectionArea::~MHMcCollectionArea()
82{
83 delete fHistAll;
84 delete fHistSel;
85 delete fHistCol;
86}
87
88// --------------------------------------------------------------------------
89//
90// Fill data into the histogram which contains all showers
91//
92void MHMcCollectionArea::FillAll(Float_t log10E, Float_t radius)
93{
94 fHistAll->Fill(log10E, radius);
95}
96
97// --------------------------------------------------------------------------
98//
99// Fill data into the histogram which contains the selected showers
100//
101void MHMcCollectionArea::FillSel(Float_t log10E, Float_t radius)
102{
103 fHistSel->Fill(log10E, radius);
104}
105
106// --------------------------------------------------------------------------
107//
108// Draw the histogram with all showers
109//
110void MHMcCollectionArea::DrawAll(Option_t* option)
111{
112 if (!gPad)
113 {
114 if (!gROOT->GetMakeDefCanvas())
115 return;
116 (gROOT->GetMakeDefCanvas())();
117 }
118
119 //TCanvas *c=new TCanvas(fHistAll->GetName(), fHistAll->GetTitle());
120
121 fHistAll->Draw(option);
122
123 gPad->Modified();
124 gPad->Update();
125}
126
127// --------------------------------------------------------------------------
128//
129// Draw the histogram with the selected showers only.
130//
131void MHMcCollectionArea::DrawSel(Option_t* option)
132{
133 if (!gPad)
134 {
135 if (!gROOT->GetMakeDefCanvas())
136 return;
137 (gROOT->GetMakeDefCanvas())();
138 }
139
140 //TCanvas *c=new TCanvas(fHistSel->GetName(), fHistSel->GetTitle());
141
142 fHistSel->Draw(option);
143
144 gPad->Modified();
145 gPad->Update();
146}
147
148// --------------------------------------------------------------------------
149//
150// Creates a new canvas and draws the histogram into it.
151// Be careful: The histogram belongs to this object and won't get deleted
152// together with the canvas.
153//
154TObject *MHMcCollectionArea::DrawClone(Option_t* option) const
155{
156 TCanvas *c=new TCanvas(fHistCol->GetName(), fHistCol->GetTitle());
157
158 //
159 // This is necessary to get the expected bahviour of DrawClone
160 //
161 gROOT->SetSelectedPad(NULL);
162
163 fHistCol->DrawClone(option);
164
165 c->Modified();
166 c->Update();
167
168 return c;
169}
170
171void MHMcCollectionArea::Draw(Option_t* option)
172{
173 if (!gPad)
174 {
175 if (!gROOT->GetMakeDefCanvas())
176 return;
177 (gROOT->GetMakeDefCanvas())();
178 }
179 // TCanvas *c=new TCanvas(fHistCol->GetName(), fHistCol->GetTitle());
180
181 fHistCol->Draw(option);
182
183 gPad->Modified();
184 gPad->Update();
185}
186
187// --------------------------------------------------------------------------
188//
189// Calculate the Efficiency (collection area)
190//
191void MHMcCollectionArea::CalcEfficiency()
192{
193 // Description!
194
195 //
196 // first of all calculate the efficency
197 // do it here by hand to get the right error of efficency
198 //
199 const Int_t nbinx = fHistSel->GetXaxis()->GetNbins();
200 const Int_t nbiny = fHistSel->GetYaxis()->GetNbins();
201
202 for (Int_t ix=1; ix<=nbinx; ix++)
203 {
204 for (Int_t iy=1; iy<=nbiny; iy++)
205 {
206 const Float_t Na = fHistAll->GetCellContent(ix, iy);
207
208 if (Na <= 0)
209 continue;
210
211 const Float_t Ns = fHistSel->GetCellContent(ix, iy);
212
213 const Double_t eff = Ns/Na;
214 const Double_t err = sqrt((1.-eff)*Ns)/Na;
215
216 // old calculation from Harald:
217 // const Double_t eff = Ns/Na;
218 // const Double_t err = sqrt(Na + Na*Ns - Ns*Ns - Ns) / (Na*Na);
219
220 fHistSel->SetCellContent(ix, iy, eff);
221 fHistSel->SetCellError(ix, iy, err);
222 }
223 }
224
225 //
226 // now calculate the Collection area for different
227 // energies
228 //
229 for (Int_t ix=1; ix<=nbinx; ix++)
230 {
231 Double_t errA = 0;
232 Double_t colA = 0;
233
234 for (Int_t iy=1; iy<=nbiny; iy++)
235 {
236 TAxis *yaxis = fHistSel->GetYaxis();
237
238 const Double_t r1 = yaxis->GetBinLowEdge(iy);
239 const Double_t r2 = yaxis->GetBinLowEdge(iy+1);
240
241 const Double_t A = TMath::Pi() * (r2*r2 - r1*r1);
242
243 const Double_t eff = fHistSel->GetCellContent(ix, iy);
244 const Double_t err = fHistSel->GetCellError(ix, iy);
245
246 colA += eff*A;
247 errA += A*A * err*err;
248 }
249
250 fHistCol->SetBinContent(ix, colA);
251 fHistCol->SetBinError(ix, sqrt(errA));
252 }
253}
Note: See TracBrowser for help on using the repository browser.