source: trunk/MagicSoft/Mars/mcalib/MMcCalibrationCalc.cc@ 4592

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