source: tags/Mars-V0.9.3/mcalib/MMcCalibrationCalc.cc

Last change on this file was 7005, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 10.6 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): Abelardo Moralejo, 12/2003 <mailto:moralejo@pd.infn.it>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MMcCalibrationCalc
28//
29// Input Containers:
30// MMcConfigRunHeader
31// MRawRunHeader
32// MMcFadcHeader
33// MHillas
34// MNewImagePar
35// MMcEvt
36//
37// Output Containers:
38// (containers mus exist already, they are filled with new values).
39// MCalibrationChargeCam
40// MCalibrationQECam
41//
42/////////////////////////////////////////////////////////////////////////////
43#include "MMcCalibrationCalc.h"
44
45#include <TH1.h>
46
47#include "MLog.h"
48#include "MLogManip.h"
49
50#include "MParList.h"
51
52#include "MCalibrationChargePix.h"
53#include "MCalibrationChargeCam.h"
54
55#include "MCalibrationQEPix.h"
56#include "MCalibrationQECam.h"
57
58#include "MGeomCam.h"
59#include "MRawRunHeader.h"
60#include "MMcConfigRunHeader.h"
61
62#include "MHillas.h"
63#include "MImagePar.h"
64#include "MNewImagePar.h"
65
66#include "MMcEvt.hxx"
67#include "MMcFadcHeader.hxx"
68
69ClassImp(MMcCalibrationCalc);
70
71using namespace std;
72
73MMcCalibrationCalc::MMcCalibrationCalc(const char *name, const char *title)
74{
75 fName = name ? name : "MMcCalibrationCalc";
76 fTitle = title ? title : "Calculate and write conversion factors into MCalibrationChargeCam and MCalibrationQECam containers";
77
78 fHistADC2PhotEl = new TH1F(AddSerialNumber("ADC2PhotEl"), "log10(fPhotElfromShower/fSize)", 1500, -3., 3.);
79 fHistADC2PhotEl->SetXTitle("log_{10}(fPhotElfromShower / fSize) [photel/ADC count]");
80
81
82 fHistPhot2PhotEl = new TH1F(AddSerialNumber("Phot2PhotEl"), "Photon conversion efficiency", 1000, 0., 1.);
83 fHistPhot2PhotEl->SetXTitle("Overall photon conversion efficiency [photoelectron/photon]");
84
85}
86
87// --------------------------------------------------------------------------
88//
89// Check for the run type. Return kTRUE if it is a MC run or if there
90// is no MC run header (old camera files) kFALSE in case of a different
91// run type
92//
93Bool_t MMcCalibrationCalc::CheckRunType(MParList *pList) const
94{
95 const MRawRunHeader *run = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
96 if (!run)
97 {
98 *fLog << warn << "Warning - cannot check file type, MRawRunHeader not found." << endl;
99 return kTRUE;
100 }
101
102 return run->IsMonteCarloRun();
103}
104
105// --------------------------------------------------------------------------
106//
107// Make sure, that there is an MCalibrationCam Object in the Parameter List.
108//
109Int_t MMcCalibrationCalc::PreProcess(MParList *pList)
110{
111 fHistADC2PhotEl->Reset();
112 fHistPhot2PhotEl->Reset();
113
114 fADC2PhotEl = 0;
115 fPhot2PhotEl = 0;
116
117 fCalCam = (MCalibrationChargeCam*) pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
118 if (!fCalCam)
119 {
120 *fLog << err << AddSerialNumber("MCalibrationChargeCam") << "not found... aborting." << endl;
121 return kFALSE;
122 }
123
124 fQECam = (MCalibrationQECam*) pList->FindObject(AddSerialNumber("MCalibrationQECam"));
125 if (!fQECam)
126 {
127 *fLog << err << AddSerialNumber("MCalibrationQECam") << "not found... aborting." << endl;
128 return kFALSE;
129 }
130
131 fHillas = (MHillas*) pList->FindObject(AddSerialNumber("MHillas"));
132 if ( !fHillas)
133 {
134 *fLog << err << AddSerialNumber("MHillas") << "not found... aborting." << endl;
135 return kFALSE;
136 }
137
138 fNew = (MNewImagePar*)pList->FindObject(AddSerialNumber("MNewImagePar"));
139 if (!fNew)
140 {
141 *fLog << err << AddSerialNumber("MNewImagePar") << "not found... aborting." << endl;
142 return kFALSE;
143 }
144
145 fPar = (MImagePar*)pList->FindObject(AddSerialNumber("MImagePar"));
146 if (!fPar)
147 {
148 *fLog << err << AddSerialNumber("MImagePar") << "not found... aborting." << endl;
149 return kFALSE;
150 }
151
152 fMcEvt = (MMcEvt*) pList->FindObject(AddSerialNumber("MMcEvt"));
153 if (!fMcEvt)
154 {
155 *fLog << err << AddSerialNumber("MMcEvt") << "not found... aborting." << endl;
156 return kFALSE;
157 }
158
159 return kTRUE;
160}
161
162// --------------------------------------------------------------------------
163//
164// Check for the runtype.
165// Search for MGeomCam and MMcFadcHeader.
166//
167Bool_t MMcCalibrationCalc::ReInit(MParList *pList)
168{
169 //
170 // If it is no MC file display error and exit
171 //
172 if (!CheckRunType(pList))
173 {
174 *fLog << err << "MMcCalibrationCalc can only be used with MC files... aborting." << endl;
175 return kFALSE;
176 }
177
178 //
179 // Now check the existence of all necessary containers.
180 //
181 fGeom = (MGeomCam*) pList->FindObject(AddSerialNumber("MGeomCam"));
182 if (!fGeom)
183 {
184 *fLog << err << AddSerialNumber("MGeomCam") << " not found... aborting." << endl;
185 return kFALSE;
186 }
187
188 fHeaderFadc = (MMcFadcHeader*)pList->FindObject(AddSerialNumber("MMcFadcHeader"));
189 if (!fHeaderFadc)
190 {
191 *fLog << err << AddSerialNumber("MMcFadcHeader") << " not found... aborting." << endl;
192 return kFALSE;
193 }
194
195 for (UInt_t ipix = 0; ipix < fGeom->GetNumPixels(); ipix++)
196 {
197 if (fHeaderFadc->GetPedestalRmsHigh(ipix) > 0 ||
198 fHeaderFadc->GetPedestalRmsLow(ipix) > 0 )
199 {
200 *fLog << err << "Trying to calibrate the data using a Camera file produced with added noise." << endl;
201 *fLog << "Please use a noiseless file for calibration... aborting." << endl << endl;
202 return kFALSE;
203 }
204 }
205
206 // Now check the light collection for inner and outer pixels to
207 // calculate the ratio between the two. FIXME! Light collection
208 // depends on the incidence angle of the light w.r.t. the camera
209 // plane. For the moment we take the ratio for light impinging
210 // perpendicular to the camera plane.
211 //
212 // FIXME! We should look for AddSerialNumber("MMcConfigRunHeader") but
213 // for the moment the stereo version of camera does not write one such
214 // header per telescope (it should!)
215 //
216 MMcConfigRunHeader* mcconfig = (MMcConfigRunHeader*) pList->FindObject("MMcConfigRunHeader");
217 if (!mcconfig)
218 {
219 *fLog << err << "MMcConfigRunHeader" << " not found... aborting." << endl;
220 return kFALSE;
221 }
222 TArrayF innerlightcoll = mcconfig->GetLightCollectionFactor();
223 TArrayF outerlightcoll = mcconfig->GetLightCollectionFactorOuter();
224
225 // In principle outer pixels seem to have a different average light
226 // collection efficiency than outer ones. We set here the factor between
227 // the two.
228
229 fOuterPixelsLightCollection = outerlightcoll[90] / innerlightcoll[90];
230 // (at angle = 90 deg)
231
232
233 return kTRUE;
234}
235
236
237// --------------------------------------------------------------------------
238//
239// Obtain average ratio of photons in camera to image Size.
240//
241Int_t MMcCalibrationCalc::Process()
242{
243 //
244 // Exclude events with some high-gain saturated pixel
245 //
246 if (fPar->GetNumSatPixelsHG()>0)
247 return kTRUE;
248
249 const Float_t size = fHillas->GetSize();
250 const Float_t innersize = fNew->GetInnerSize();
251
252 // Size will at this point be in ADC counts (still uncalibrated)
253 //
254 // Exclude events with low Size (larger fluctuations)
255 // FIXME? The present cut (1000 "inner-pixel-counts") is somehow
256 // arbitrary. Might it be optimized?
257 //
258
259 if (size < 1000)
260 return kTRUE;
261
262 //
263 // PATCH: Convert number of photoelectrons in camera to the approximate number
264 // of photoelectrons that would have been registered if all pixels had the same
265 // light collection efficiency as inner ones (called here "corrected_photel").
266 //
267
268 const Float_t corrected_photel = (Float_t) fMcEvt->GetPhotElfromShower() /
269 (fOuterPixelsLightCollection + innersize / size * (1. - fOuterPixelsLightCollection));
270
271 fHistADC2PhotEl->Fill(TMath::Log10(corrected_photel/size));
272 fHistPhot2PhotEl->Fill( corrected_photel / (Float_t) fMcEvt->GetPassPhotCone() );
273
274 return kTRUE;
275}
276
277// --------------------------------------------------------------------------
278//
279// Fill the MCalibrationCam object
280//
281Int_t MMcCalibrationCalc::PostProcess()
282{
283 const Stat_t n = fHistADC2PhotEl->GetEntries();
284 if (n<1)
285 {
286 *fLog << err << "No events read... aborting." << endl;
287 return kFALSE;
288 }
289
290 fPhot2PhotEl = fHistPhot2PhotEl->GetMean(); // Average quantum efficiency
291
292 //
293 // Find peak of log10(photel/Size) histogram:
294 //
295 const Int_t reach = 2;
296 Stat_t summax = 0;
297 Int_t mode = 0;
298 for (Int_t ibin = 1+reach; ibin <= fHistADC2PhotEl->GetNbinsX()-reach; ibin++)
299 {
300 const Stat_t sum = fHistADC2PhotEl->Integral(ibin-reach, ibin+reach);
301
302 if (sum <= summax)
303 continue;
304
305 summax = sum;
306 mode = ibin;
307 }
308
309 fADC2PhotEl = TMath::Power(10, fHistADC2PhotEl->GetBinCenter(mode));
310
311 const Int_t num = fCalCam->GetSize();
312
313 for (int i=0; i<num; i++)
314 {
315 MCalibrationChargePix &calpix = (MCalibrationChargePix&)(*fCalCam)[i];
316 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam) [i];
317
318 Float_t qe = fPhot2PhotEl;
319 if (fGeom->GetPixRatio(i) < 1.)
320 qe *= fOuterPixelsLightCollection;
321 qepix.SetAverageQE(qe);
322
323 qepix.SetAvNormFFactor(1.);
324 // This factor should convert the default average QE for different colors to
325 // average QE for a spectrum like that of Cherenkov light (see the documentration
326 // of MCalibrationQEPix).
327 // Here we obtain average QE using already a Cherenkov spectrum so AvNormFFactor
328 // must be 1.
329
330
331 Float_t factor = fADC2PhotEl;
332
333 //
334 // We take into account the (possibly) different gain of outer pixels:
335 // FIXME: we are now assuming that all inner pixels have the same gain, and all
336 // outer pixels have the same gain (different from inner ones though). This can
337 // only be like this in camera 0.7, but may change in future versions of camera.
338 //
339
340 if (fGeom->GetPixRatio(i) < 1.)
341 factor *= fHeaderFadc->GetAmplitud()/fHeaderFadc->GetAmplitudOuter();
342
343 calpix.SetMeanConvFADC2Phe(factor);
344 calpix.SetMeanConvFADC2PheVar(0.);
345
346 calpix.SetMeanFFactorFADC2Phot(0.);
347 }
348
349 return kTRUE;
350}
Note: See TracBrowser for help on using the repository browser.