source: trunk/MagicSoft/Mars/mmontecarlo/MCollArea.cc@ 769

Last change on this file since 769 was 749, checked in by tbretz, 24 years ago
*** empty log message ***
File size: 4.2 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! Author(s): Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
20!
21! Copyright: MAGIC Software Development, 2000-2001
22!
23!
24\* ======================================================================== */
25
26#include "MCollArea.h"
27
28#include <MLog.h>
29#include <TH2.h>
30
31ClassImp(MCollArea)
32
33
34
35MCollArea::MCollArea(const char *name, const char *title)
36{
37 //
38 // default constructor
39 //
40
41 // initialize the histogram for the distribution r vs E
42 //
43 // we set the energy range from 1 Gev to 10000 GeV (in log 5 orders
44 // of magnitude) and for each order we take 10 subdivision --> 50 xbins
45 //
46 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
47
48
49 *fName = name ? name : "MCollArea";
50 *fTitle = title ? title : "Data to Calculate Coll-Area";
51
52
53 fHistAll = new TH2D("collAreaAll", "all showers - Radius vs log(E) distribution",
54 50, 0., 5.,
55 50, 0., 500. ) ;
56
57 fHistSel = new TH2D("collAreaSel", "selected showers - Radius vs log(E) distribution",
58 50, 0., 5.,
59 50, 0., 500. ) ;
60
61 fHistCol = new TH1D("collArea", "Collection Area",
62 50, 0., 5.) ;
63
64}
65
66MCollArea::~MCollArea()
67{
68 delete fHistAll ;
69 delete fHistSel ;
70 delete fHistCol ;
71}
72
73void MCollArea::FillAll(Float_t log10E, Float_t radius)
74{
75 fHistAll->Fill(log10E, radius ) ;
76}
77
78void MCollArea::FillSel(Float_t log10E, Float_t radius)
79{
80 fHistSel->Fill(log10E, radius ) ;
81}
82
83void MCollArea::DrawAll()
84{
85 fHistAll->Draw() ;
86}
87
88void MCollArea::DrawSel()
89{
90 fHistSel->Draw() ;
91}
92
93void MCollArea::Draw(Option_t* option)
94{
95 fHistCol->Draw(option) ;
96}
97
98void MCollArea::CalcEfficiency()
99{
100 // Description!
101
102 //
103 // first of all calculate the efficency
104 // do it here by hand to get the right error of efficency
105 //
106 const Int_t iBinx = ((TAxis *)fHistSel->GetXaxis())->GetNbins();
107 const Int_t iBiny = ((TAxis *)fHistSel->GetYaxis())->GetNbins();
108
109 for (Int_t ix=1; ix<=iBiny; ix++ )
110 {
111 for (Int_t iy=1; iy<=iBiny; iy++ )
112 {
113 const Float_t N = fHistSel->GetCellContent(ix, iy);
114 const Float_t Nall = fHistAll->GetCellContent(ix, iy);
115
116 if ( Nall <= 0 ) {
117 // cout << ix << " " << iy << endl ;
118 continue;
119 }
120
121 const Double_t eff = N / Nall ;
122 const Double_t err = sqrt(Nall + Nall*N - N*N - N) / (Nall*Nall);
123 /*
124 cout << ix << " " << iy
125 << " N " << N
126 << " Nall " << Nall
127 << " effi " << eff
128 << " error " << err
129 << endl ;
130 */
131 fHistSel->SetCellContent(ix, iy, eff);
132 fHistSel->SetCellError(ix, iy, err);
133 }
134 }
135
136 //
137 // now calculate the Collection area for different
138 // energies
139 //
140 for (Int_t ix=1; ix<=iBiny; ix++ )
141 {
142 Double_t errA = 0;
143 Double_t colA = 0;
144
145 for (Int_t iy=1; iy<=iBiny; iy++ )
146 {
147 const Double_t r1 = ((TAxis *)fHistSel->GetYaxis())->GetBinLowEdge(iy);
148 const Double_t r2 = ((TAxis *)fHistSel->GetYaxis())->GetBinLowEdge(iy+1);
149
150 const Double_t A = kPI * (r2*r2 - r1*r1);
151
152 const Double_t eff = fHistSel->GetCellContent(ix, iy);
153 const Double_t err = fHistSel->GetCellError(ix, iy);
154
155 colA += eff*A;
156 errA += (A*A) * err * err;
157 }
158
159 fHistCol->SetBinContent(ix, colA);
160 fHistCol->SetBinError(ix, sqrt(errA));
161 }
162}
Note: See TracBrowser for help on using the repository browser.