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

Last change on this file since 4166 was 3852, checked in by moralejo, 21 years ago
*** empty log message ***
File size: 8.4 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// MRawRunHeader
31// MMcFadcHeader
32// MHillas
33// MNewImagePar
34// MMcEvt
35//
36// Output Containers:
37// (containers mus exist already, they are filled with new values).
38// MCalibrationChargeCam
39// MCalibrationQECam
40//
41/////////////////////////////////////////////////////////////////////////////
42#include "MMcCalibrationCalc.h"
43
44#include <TH1.h>
45
46#include "MLog.h"
47#include "MLogManip.h"
48
49#include "MParList.h"
50
51#include "MCalibrationChargePix.h"
52#include "MCalibrationChargeCam.h"
53
54#include "MCalibrationQEPix.h"
55#include "MCalibrationQECam.h"
56
57#include "MGeomCam.h"
58#include "MRawRunHeader.h"
59
60#include "MHillas.h"
61#include "MNewImagePar.h"
62
63#include "MMcEvt.hxx"
64#include "MMcFadcHeader.hxx"
65
66ClassImp(MMcCalibrationCalc);
67
68using namespace std;
69
70MMcCalibrationCalc::MMcCalibrationCalc(const char *name, const char *title)
71{
72 fName = name ? name : "MMcCalibrationCalc";
73 fTitle = title ? title : "Calculate and write conversion factors into MCalibrationChargeCam and MCalibrationQECam containers";
74
75 fHistADC2PhotEl = new TH1F(AddSerialNumber("ADC2PhotEl"), "log10(fPhotElfromShower/fSize)", 1500, -3., 3.);
76 fHistADC2PhotEl->SetXTitle("log_{10}(fPhotElfromShower / fSize) [photel/ADC count]");
77
78
79 fHistPhot2PhotEl = new TH1F(AddSerialNumber("Phot2PhotEl"), "Photon conversion efficiency", 1000, 0., 1.);
80 fHistPhot2PhotEl->SetXTitle("Overall photon conversion efficiency [photoelectron/photon]");
81
82}
83
84// --------------------------------------------------------------------------
85//
86// Check for the run type. Return kTRUE if it is a MC run or if there
87// is no MC run header (old camera files) kFALSE in case of a different
88// run type
89//
90Bool_t MMcCalibrationCalc::CheckRunType(MParList *pList) const
91{
92 const MRawRunHeader *run = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
93 if (!run)
94 {
95 *fLog << warn << "Warning - cannot check file type, MRawRunHeader not found." << endl;
96 return kTRUE;
97 }
98
99 return run->IsMonteCarloRun();
100}
101
102// --------------------------------------------------------------------------
103//
104// Make sure, that there is an MCalibrationCam Object in the Parameter List.
105//
106Int_t MMcCalibrationCalc::PreProcess(MParList *pList)
107{
108 fHistADC2PhotEl->Reset();
109 fHistPhot2PhotEl->Reset();
110
111 fADC2PhotEl = 0;
112 fPhot2PhotEl = 0;
113
114 fCalCam = (MCalibrationChargeCam*) pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
115 if (!fCalCam)
116 {
117 *fLog << err << AddSerialNumber("MCalibrationChargeCam") << "not found... aborting." << endl;
118 return kFALSE;
119 }
120
121 fQECam = (MCalibrationQECam*) pList->FindObject(AddSerialNumber("MCalibrationQECam"));
122 if (!fQECam)
123 {
124 *fLog << err << AddSerialNumber("MCalibrationQECam") << "not found... aborting." << endl;
125 return kFALSE;
126 }
127
128 fHillas = (MHillas*) pList->FindObject(AddSerialNumber("MHillas"));
129 if ( !fHillas)
130 {
131 *fLog << err << AddSerialNumber("MHillas") << "not found... aborting." << endl;
132 return kFALSE;
133 }
134
135 fNew = (MNewImagePar*)pList->FindObject(AddSerialNumber("MNewImagePar"));
136 if (!fNew)
137 {
138 *fLog << err << AddSerialNumber("MNewImagePar") << "not found... aborting." << endl;
139 return kFALSE;
140 }
141
142 fMcEvt = (MMcEvt*) pList->FindObject(AddSerialNumber("MMcEvt"));
143 if (!fMcEvt)
144 {
145 *fLog << err << AddSerialNumber("MMcEvt") << "not found... aborting." << endl;
146 return kFALSE;
147 }
148
149 return kTRUE;
150}
151
152// --------------------------------------------------------------------------
153//
154// Check for the runtype.
155// Search for MGeomCam and MMcFadcHeader.
156//
157Bool_t MMcCalibrationCalc::ReInit(MParList *pList)
158{
159 //
160 // If it is no MC file display error and exit
161 //
162 if (!CheckRunType(pList))
163 {
164 *fLog << err << "MMcCalibrationCalc can only be used with MC files... aborting." << endl;
165 return kFALSE;
166 }
167
168 //
169 // Now check the existence of all necessary containers.
170 //
171 fGeom = (MGeomCam*) pList->FindObject(AddSerialNumber("MGeomCam"));
172 if (!fGeom)
173 {
174 *fLog << err << AddSerialNumber("MGeomCam") << " not found... aborting." << endl;
175 return kFALSE;
176 }
177
178 fHeaderFadc = (MMcFadcHeader*)pList->FindObject(AddSerialNumber("MMcFadcHeader"));
179 if (!fHeaderFadc)
180 {
181 *fLog << err << AddSerialNumber("MMcFadcHeader") << " not found... aborting." << endl;
182 return kFALSE;
183 }
184
185 for (UInt_t ipix = 0; ipix < fGeom->GetNumPixels(); ipix++)
186 {
187 if (fHeaderFadc->GetPedestalRmsHigh(ipix) > 0 ||
188 fHeaderFadc->GetPedestalRmsLow(ipix) > 0 )
189 {
190 *fLog << err << "Trying to calibrate the data using a Camera file produced with added noise." << endl;
191 *fLog << "Please use a noiseless file for calibration... aborting." << endl << endl;
192 return kFALSE;
193 }
194 }
195
196 return kTRUE;
197}
198
199
200// --------------------------------------------------------------------------
201//
202// Obtain average ratio of photons in camera to image Size.
203//
204Int_t MMcCalibrationCalc::Process()
205{
206 //
207 // Exclude events with some saturated pixel
208 //
209 if (fNew->GetNumSaturatedPixels()>0)
210 return kTRUE;
211
212 //
213 // Exclude events with low Size (larger fluctuations)
214 // FIXME? The present cut (1000 "inner-pixel-counts") is somehow
215 // arbitrary. Might it be optimized?
216 //
217
218 const Float_t size = fHillas->GetSize();
219 // Size will at this point be in ADC counts (still uncalibrated)
220
221 if (size < 1000)
222 return kTRUE;
223
224 fHistADC2PhotEl->Fill(TMath::Log10(fMcEvt->GetPhotElfromShower()/size));
225 fHistPhot2PhotEl->Fill( (Float_t) fMcEvt->GetPhotElfromShower() /
226 (Float_t) fMcEvt->GetPassPhotCone() );
227
228 return kTRUE;
229}
230
231// --------------------------------------------------------------------------
232//
233// Fill the MCalibrationCam object
234//
235Int_t MMcCalibrationCalc::PostProcess()
236{
237 const Stat_t n = fHistADC2PhotEl->GetEntries();
238 if (n<1)
239 {
240 *fLog << err << "No events read... aborting." << endl;
241 return kFALSE;
242 }
243
244 fPhot2PhotEl = fHistPhot2PhotEl->GetMean(); // Average quantum efficiency
245
246 //
247 // Find peak of log10(photel/Size) histogram:
248 //
249 const Int_t reach = 2;
250 Stat_t summax = 0;
251 Int_t mode = 0;
252 for (Int_t ibin = 1+reach; ibin <= fHistADC2PhotEl->GetNbinsX()-reach; ibin++)
253 {
254 const Stat_t sum = fHistADC2PhotEl->Integral(ibin-reach, ibin+reach);
255
256 if (sum <= summax)
257 continue;
258
259 summax = sum;
260 mode = ibin;
261 }
262
263 fADC2PhotEl = TMath::Power(10, fHistADC2PhotEl->GetBinCenter(mode));
264
265 const Int_t num = fCalCam->GetSize();
266
267 for (int i=0; i<num; i++)
268 {
269 MCalibrationChargePix &calpix = (MCalibrationChargePix&)(*fCalCam)[i];
270 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam) [i];
271
272 qepix.SetAverageQE(fPhot2PhotEl);
273
274 qepix.SetAvNormFFactor(1.);
275 // This factor should convert the default average QE for different colors to
276 // average QE for a spectrum like that of Cherenkov light (see the documentration
277 // of MCalibrationQEPix).
278 // Here we obtain average QE using already a Cherenkov spectrum so AvNormFFactor
279 // must be 1.
280
281
282 Float_t factor = fADC2PhotEl;
283
284 //
285 // We take into account the (possibly) different gain of outer pixels:
286 //
287 if (fGeom->GetPixRatio(i) < 1.)
288 factor *= fHeaderFadc->GetAmplitud()/fHeaderFadc->GetAmplitudOuter();
289
290 calpix.SetMeanConvFADC2Phe(factor);
291 calpix.SetMeanConvFADC2PheVar(0.);
292
293 calpix.SetMeanFFactorFADC2Phot(0.);
294 }
295
296 return kTRUE;
297}
Note: See TracBrowser for help on using the repository browser.