source: trunk/MagicSoft/Mars/mmontecarlo/MMcCollectionAreaCalc.cc@ 3975

Last change on this file since 3975 was 2252, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 6.3 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// MHMcCollectionAreaCalc
29//
30// Remark: The initialization is mainly done in the ReInit function.
31// Please make sure, that you don't use MReadTree when processing
32// a file. Use a 'ReInit'-calling task like MReadMarsFile
33//
34//////////////////////////////////////////////////////////////////////////////
35
36#include "MMcCollectionAreaCalc.h"
37
38#include "MParList.h"
39
40#include "MLog.h"
41#include "MLogManip.h"
42
43#include "MMcEvt.hxx"
44#include "MMcTrig.hxx"
45#include "MMcRunHeader.hxx"
46#include "MMcCorsikaRunHeader.h"
47
48#include "MHMcCollectionArea.h"
49
50ClassImp(MMcCollectionAreaCalc);
51
52using namespace std;
53
54MMcCollectionAreaCalc::MMcCollectionAreaCalc(const char *input,
55 const char *name, const char *title)
56{
57 fName = name ? name : "MMcCollectionAreaCalc";
58 fTitle = title ? title : "Task to calculate the collection area";
59
60 fObjName = input ? input : "MMcTrig";
61 AddToBranchList(Form("%s.fNumFirstLevel", input));
62
63 AddToBranchList("MMcEvt.fEnergy");
64 AddToBranchList("MMcEvt.fImpact");
65
66 fCorsikaVersion = 0;
67 fAllEvtsTriggered = kFALSE;
68 fTotalNumSimulatedShowers = 0;
69 fTheta = -1.;
70
71}
72
73Int_t MMcCollectionAreaCalc::PreProcess (MParList *pList)
74{
75 // connect the raw data with this task
76
77 fMcEvt = (MMcEvt*)pList->FindObject("MMcEvt");
78 if (!fMcEvt)
79 {
80 *fLog << err << dbginf << "MMcEvt not found... exit." << endl;
81 return kFALSE;
82 }
83
84 return kTRUE;
85}
86
87Bool_t MMcCollectionAreaCalc::ReInit(MParList *plist)
88{
89 MMcRunHeader *runheader = (MMcRunHeader*)plist->FindObject("MMcRunHeader");
90 if (!runheader)
91 {
92 *fLog << err << dbginf << "Error - MMcRunHeader not found... exit." << endl;
93 return kFALSE;
94 }
95
96 MMcCorsikaRunHeader *corrunheader = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
97 if (!corrunheader)
98 {
99 *fLog << err << dbginf << "Error - MMcCorsikaRunHeader not found... exit." << endl;
100 return kFALSE;
101 }
102
103 fTotalNumSimulatedShowers += runheader->GetNumSimulatedShowers();
104 *fLog << inf << "Total Number of Simulated showers: " << fTotalNumSimulatedShowers << endl;
105
106
107 if (fTheta>=0 && fTheta!=runheader->GetTelesTheta())
108 {
109 *fLog << warn << dbginf << "Warning - Read files have different TelesTheta (";
110 *fLog << fTheta << ", " << runheader->GetTelesTheta() << ")..." << endl;
111 }
112 fTheta = runheader->GetTelesTheta();
113
114
115 if (fCorsikaVersion!=0 && fCorsikaVersion!=runheader->GetCorsikaVersion())
116 *fLog << warn << dbginf << "Warning - Read files have different Corsika versions..." << endl;
117 fCorsikaVersion = runheader->GetCorsikaVersion();
118
119
120 fAllEvtsTriggered |= runheader->GetAllEvtsTriggered();
121 *fLog << inf << "Only triggered events avail: " << (fAllEvtsTriggered?"yes":"no") << endl;
122
123
124 fCollArea = (MHMcCollectionArea*)plist->FindCreateObj("MHMcCollectionArea");
125 if (!fCollArea)
126 return kFALSE;
127
128 if (fAllEvtsTriggered)
129 {
130 //
131 // Calculate approximately the original number of events in each
132 // energy bin:
133 //
134
135 const Float_t emin = corrunheader->GetELowLim();
136 const Float_t emax = corrunheader->GetEUppLim();
137 const Float_t expo = 1 + corrunheader->GetSlopeSpec();
138 const Float_t k = runheader->GetNumSimulatedShowers() /
139 (pow(emax,expo) - pow(emin,expo)) ;
140
141 TH2 &hall = *fCollArea->GetHistAll();
142
143 const Int_t nbinx = hall.GetNbinsX();
144
145 TAxis &axe = *hall.GetXaxis();
146 for (Int_t i = 1; i <= nbinx; i++)
147 {
148 const Float_t e1 = axe.GetBinLowEdge(i);
149 const Float_t e2 = axe.GetBinLowEdge(i+1);
150
151 if (e1 < emin || e2 > emax)
152 continue;
153
154 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
155 //
156 // We fill the i-th energy bin, with the total number of events
157 // Second argument of Fill would be impact parameter of each
158 // event, but we don't really need it for the collection area,
159 // so we just put a dummy value (1.)
160 //
161
162 const Float_t energy = (e1+e2)/2.;
163 hall.Fill(energy, 1., events);
164 }
165
166 return kTRUE;
167 }
168
169 fMcTrig = (MMcTrig*)plist->FindObject(fObjName);
170 if (!fMcTrig)
171 {
172 *fLog << err << dbginf << fObjName << " [MMcTrig] not found... exit." << endl;
173 return kFALSE;
174 }
175
176 return kTRUE;
177}
178
179Int_t MMcCollectionAreaCalc::Process()
180{
181 const Double_t energy = fMcEvt->GetEnergy();
182 const Double_t impact = fMcEvt->GetImpact()/100.;
183
184 //
185 // This happens for camera files created with Camera 0.5
186 //
187 if (TMath::IsNaN(impact))
188 {
189 *fLog << err << dbginf << "ERROR - Impact=NaN (Not a number)... abort." << endl;
190 return kERROR;
191 }
192
193 if (!fAllEvtsTriggered)
194 {
195 fCollArea->FillAll(energy, impact);
196
197 if (fMcTrig->GetFirstLevel() <= 0)
198 return kTRUE;
199 }
200
201 fCollArea->FillSel(energy, impact);
202
203 return kTRUE;
204}
205
206Int_t MMcCollectionAreaCalc::PostProcess()
207{
208 //
209 // do the calculation of the effectiv area
210 //
211 *fLog << inf << "Calculation Collection Area..." << endl;
212
213 if (!fAllEvtsTriggered)
214 fCollArea->CalcEfficiency();
215 else
216 {
217 *fLog << inf << "Total number of showers: " << fTotalNumSimulatedShowers << endl;
218 fCollArea->CalcEfficiency2();
219 }
220
221 return kTRUE;
222}
223
Note: See TracBrowser for help on using the repository browser.