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 | ! Author(s): Josep Flix, 01/2004 <mailto:jflix@ifae.es>
|
---|
18 | ! Javier Rico, 01/2004 <mailto:jrico@ifae.es>
|
---|
19 | ! Markus Gaug, 03/2004 <mailto:markus@ifae.es>
|
---|
20 | !
|
---|
21 | ! (based on bootcampstandardanalysis.C by Javier López)
|
---|
22 | !
|
---|
23 | !
|
---|
24 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
25 | !
|
---|
26 | !
|
---|
27 | \* ======================================================================== */
|
---|
28 |
|
---|
29 | #include "MParList.h"
|
---|
30 | #include "MTaskList.h"
|
---|
31 | #include "MPedestalCam.h"
|
---|
32 | #include "MReadMarsFile.h"
|
---|
33 | #include "MGeomApply.h"
|
---|
34 | #include "MPedCalcPedRun.h"
|
---|
35 | #include "MGeomCamMagic.h"
|
---|
36 | #include "MEvtLoop.h"
|
---|
37 | #include "MCalibrationChargeCam.h"
|
---|
38 | #include "MHCalibrationChargeCam.h"
|
---|
39 | #include "MCalibrationChargePINDiode.h"
|
---|
40 | #include "MHCalibrationChargePINDiode.h"
|
---|
41 | #include "MCalibrationChargeBlindPix.h"
|
---|
42 | #include "MHCalibrationChargeBlindPix.h"
|
---|
43 | #include "MCalibrationChargeCalc.h"
|
---|
44 | #include "MExtractedSignalCam.h"
|
---|
45 | #include "MExtractSignal.h"
|
---|
46 | #include "MExtractedSignalPINDiode.h"
|
---|
47 | #include "MExtractPINDiode.h"
|
---|
48 | #include "MExtractedSignalBlindPixel.h"
|
---|
49 | #include "MExtractBlindPixel.h"
|
---|
50 | #include "MCerPhotEvt.h"
|
---|
51 | #include "MCalibrate.h"
|
---|
52 | #include "MPedPhotCalc.h"
|
---|
53 | #include "MPedPhotCam.h"
|
---|
54 | #include "MBadPixelsCam.h"
|
---|
55 | #include "MHCamera.h"
|
---|
56 |
|
---|
57 | #include "TTimer.h"
|
---|
58 | #include "TString.h"
|
---|
59 | #include "TCanvas.h"
|
---|
60 |
|
---|
61 | #include <iostream>
|
---|
62 |
|
---|
63 |
|
---|
64 | const TString pedfile="/disc02/Data/rootdata/Miscellaneous/2003_12_19/20031218_03522_P_Park_E.root";
|
---|
65 | const TString calfile="/disc02/Data/rootdata/Miscellaneous/2003_12_19/20031218_03527_C_Park_E.root";
|
---|
66 |
|
---|
67 | void pedphotcalc(TString pedname=pedfile,TString calname=calfile)
|
---|
68 | {
|
---|
69 |
|
---|
70 | /*
|
---|
71 | This macro is an example of the use of MPedPhotCalc. It computes
|
---|
72 | the pedestal mean and rms from pedestal files undergoing the
|
---|
73 | signal extraction + calibration chain, in units of photons. It
|
---|
74 | produces plots of the relevant computed quantities.
|
---|
75 | */
|
---|
76 |
|
---|
77 | /************************************/
|
---|
78 | /* FIRST LOOP: PEDESTAL COMPUTATION */
|
---|
79 | /************************************/
|
---|
80 |
|
---|
81 | MParList plist;
|
---|
82 | MTaskList tlist;
|
---|
83 | plist.AddToList(&tlist);
|
---|
84 |
|
---|
85 | // containers
|
---|
86 | MPedestalCam pedcam;
|
---|
87 | plist.AddToList(&pedcam);
|
---|
88 |
|
---|
89 | //tasks
|
---|
90 | MReadMarsFile read("Events", pedname);
|
---|
91 | MGeomApply geomapl;
|
---|
92 | MPedCalcPedRun pedcalc;
|
---|
93 | MGeomCamMagic geomcam;
|
---|
94 |
|
---|
95 | read.DisableAutoScheme();
|
---|
96 |
|
---|
97 | tlist.AddToList(&read);
|
---|
98 | tlist.AddToList(&geomapl);
|
---|
99 | tlist.AddToList(&pedcalc);
|
---|
100 |
|
---|
101 | // Create and execute the event looper
|
---|
102 | MEvtLoop evtloop;
|
---|
103 | evtloop.SetParList(&plist);
|
---|
104 |
|
---|
105 | cout << "*************************" << endl;
|
---|
106 | cout << "** COMPUTING PEDESTALS **" << endl;
|
---|
107 | cout << "*************************" << endl;
|
---|
108 | if (!evtloop.Eventloop())
|
---|
109 | return;
|
---|
110 |
|
---|
111 | tlist.PrintStatistics();
|
---|
112 |
|
---|
113 | /****************************/
|
---|
114 | /* SECOND LOOP: CALIBRATION */
|
---|
115 | /****************************/
|
---|
116 |
|
---|
117 | MParList plist2;
|
---|
118 | MTaskList tlist2;
|
---|
119 | plist2.AddToList(&tlist2);
|
---|
120 |
|
---|
121 | // containers
|
---|
122 | MCalibrationChargeCam calcam;
|
---|
123 | MExtractedSignalCam sigcam;
|
---|
124 | MBadPixelsCam badcam;
|
---|
125 |
|
---|
126 | MCalibrationChargeCam chargecam;
|
---|
127 | MCalibrationChargePINDiode chargepin;
|
---|
128 | MCalibrationChargeBlindPixel chargeblind;
|
---|
129 | chargeblind.SetColor(MCalibrationChargeBlindPix::kECT1);
|
---|
130 | chargepin.SetColor( MCalibrationChargePINDiode::kECT1);
|
---|
131 | //
|
---|
132 | plist2.AddToList(&geomcam);
|
---|
133 | plist2.AddToList(&pedcam);
|
---|
134 | plist2.AddToList(&calcam);
|
---|
135 | plist2.AddToList(&sigcam);
|
---|
136 | plist2.AddToList(&badcam);
|
---|
137 | plist2.AddToList(&chargecam);
|
---|
138 | plist2.AddToList(&chargepin);
|
---|
139 | plist2.AddToList(&chargeblind);
|
---|
140 |
|
---|
141 | //tasks
|
---|
142 | MReadMarsFile read2("Events", calname);
|
---|
143 | read2.DisableAutoScheme();
|
---|
144 |
|
---|
145 | MExtractPINDiode pincalc;
|
---|
146 | MExtractBlindPixel blindcalc;
|
---|
147 | MExtractSignal sigcalc;
|
---|
148 | MCalibrationChargeCalc calcalc;
|
---|
149 |
|
---|
150 | MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode");
|
---|
151 | MFillH fillblind("MHCalibrationChargeBlindPixel", "MExtractedSignalBlindPixel");
|
---|
152 | MFillH fillcam("MHCalibrationChargeCam" , "MExtractedSignalCam");
|
---|
153 | //
|
---|
154 | // Apply a filter against cosmics
|
---|
155 | // (was directly in MCalibrationChargeCalc in earlier versions)
|
---|
156 | //
|
---|
157 | MFCosmics cosmics;
|
---|
158 | MContinue cont(&cosmics);
|
---|
159 |
|
---|
160 | tlist2.AddToList(&read2);
|
---|
161 | tlist2.AddToList(&geomapl);
|
---|
162 | tlist2.AddToList(&blindcalc);
|
---|
163 | tlist2.AddToList(&pincalc);
|
---|
164 | tlist2.AddToList(&sigcalc);
|
---|
165 | //
|
---|
166 | // In case, you want to skip the cosmics rejection,
|
---|
167 | // uncomment the next line
|
---|
168 | //
|
---|
169 | tlist2.AddToList(&cont);
|
---|
170 | //
|
---|
171 | tlist2.AddToList(&fillpin);
|
---|
172 | tlist2.AddToList(&fillblind);
|
---|
173 | tlist2.AddToList(&fillcam);
|
---|
174 | tlist2.AddToList(&calcalc);
|
---|
175 |
|
---|
176 | // Execute second loop (calibration)
|
---|
177 | MEvtLoop evtloop2;
|
---|
178 | evtloop2.SetParList(&plist2);
|
---|
179 |
|
---|
180 | cout << "***********************************" << endl;
|
---|
181 | cout << "** COMPUTING CALIBRATION FACTORS **" << endl;
|
---|
182 | cout << "***********************************" << endl;
|
---|
183 | if (!evtloop2.Eventloop())
|
---|
184 | return;
|
---|
185 |
|
---|
186 | tlist2.PrintStatistics();
|
---|
187 |
|
---|
188 |
|
---|
189 | /************************************************************************/
|
---|
190 | /* THIRD LOOP: PEDESTAL COMPUTATION USING EXTRACTED SIGNAL (IN PHOTONS) */
|
---|
191 | /************************************************************************/
|
---|
192 |
|
---|
193 | // Create an empty Parameter List and an empty Task List
|
---|
194 | MParList plist3;
|
---|
195 | MTaskList tlist3;
|
---|
196 | plist3.AddToList(&tlist3);
|
---|
197 |
|
---|
198 | // containers
|
---|
199 | MCerPhotEvt photcam;
|
---|
200 | MPedPhotCam pedphotcam;
|
---|
201 | MExtractedSignalCam extedsig;
|
---|
202 |
|
---|
203 | plist3.AddToList(&geomcam);
|
---|
204 | plist3.AddToList(&pedcam);
|
---|
205 | plist3.AddToList(&calcam);
|
---|
206 | plist3.AddToList(&photcam);
|
---|
207 | plist3.AddToList(&extedsig);
|
---|
208 | plist3.AddToList(&pedphotcam);
|
---|
209 |
|
---|
210 | //tasks
|
---|
211 | MReadMarsFile read3("Events", pedname);
|
---|
212 | MCalibrate photcalc;
|
---|
213 | MPedPhotCalc pedphotcalc;
|
---|
214 |
|
---|
215 | read3.DisableAutoScheme();
|
---|
216 |
|
---|
217 | tlist3.AddToList(&read3);
|
---|
218 | tlist3.AddToList(&geomapl);
|
---|
219 | tlist3.AddToList(&sigcalc);
|
---|
220 | tlist3.AddToList(&photcalc);
|
---|
221 | tlist3.AddToList(&pedphotcalc);
|
---|
222 |
|
---|
223 | // Create and execute eventloop
|
---|
224 | MEvtLoop evtloop3;
|
---|
225 | evtloop3.SetParList(&plist3);
|
---|
226 |
|
---|
227 | cout << "*************************************************************" << endl;
|
---|
228 | cout << "** COMPUTING PEDESTALS USING EXTRACTED SIGNAL (IN PHOTONS) **" << endl;
|
---|
229 | cout << "*************************************************************" << endl;
|
---|
230 |
|
---|
231 | if (!evtloop3.Eventloop())
|
---|
232 | return;
|
---|
233 | tlist3.PrintStatistics();
|
---|
234 |
|
---|
235 |
|
---|
236 | /**********************/
|
---|
237 | /* PRODUCE NICE PLOTS */
|
---|
238 | /**********************/
|
---|
239 |
|
---|
240 | const UShort_t npix = 577;
|
---|
241 |
|
---|
242 | // declare histograms
|
---|
243 | // pedestals
|
---|
244 | TH1F* pedhist = new TH1F("pedhist","Pedestal",npix,0,npix);
|
---|
245 | TH1F* pedrmshist = new TH1F("pedrmshist","Pedestal RMS",npix,0,npix);
|
---|
246 | TH1F* peddist = new TH1F("peddist","Pedestal",100,0,20);
|
---|
247 | TH1F* pedrmsdist = new TH1F("pedrmsdist","Pedestal RMS",100,0,15);
|
---|
248 |
|
---|
249 | // calibration factors
|
---|
250 | TH1F* calhist = new TH1F("calhist","Calibration factors",npix,0,npix);
|
---|
251 | TH1F* caldist = new TH1F("caldist","Calibration factors",100,0,1);
|
---|
252 |
|
---|
253 | // pedestal signals
|
---|
254 | TH1F* pedphothist = new TH1F("pedphothist","Pedestal Signal",npix,0,npix);
|
---|
255 | TH1F* pedphotrmshist = new TH1F("pedphotrmshist","Pedestal Signal RMS",npix,0,npix);
|
---|
256 | TH1F* pedphotdist = new TH1F("pedphotdist","Pedestal Signal",100,-0.4,0.4);
|
---|
257 | TH1F* pedphotrmsdist = new TH1F("pedphotrmsdist","Pedestal Signal RMS",100,0,25);
|
---|
258 |
|
---|
259 | // fill histograms
|
---|
260 | for(int i=0;i<npix;i++)
|
---|
261 | {
|
---|
262 | const Float_t ped = pedcam[i].GetPedestal();
|
---|
263 | const Float_t pedrms = pedcam[i].GetPedestalRms();
|
---|
264 | const Float_t cal = calcam[i].GetMeanConversionBlindPixelMethod();
|
---|
265 | const Float_t calerr = calcam[i].GetErrorConversionBlindPixelMethod();
|
---|
266 | const Float_t pedphot = pedphotcam[i].GetMean();
|
---|
267 | const Float_t pedphotrms = pedphotcam[i].GetRms();
|
---|
268 |
|
---|
269 | pedhist->Fill(i,ped);
|
---|
270 | peddist->Fill(ped);
|
---|
271 | pedrmshist->Fill(i,pedrms);
|
---|
272 | pedrmsdist->Fill(pedrms);
|
---|
273 |
|
---|
274 | calhist->Fill(i,cal);
|
---|
275 | caldist->Fill(cal);
|
---|
276 |
|
---|
277 | pedphothist->Fill(i,pedphot);
|
---|
278 | pedphotdist->Fill(pedphot);
|
---|
279 | pedphotrmshist->Fill(i,pedphotrms);
|
---|
280 | pedphotrmsdist->Fill(pedphotrms);
|
---|
281 | }
|
---|
282 |
|
---|
283 | // Draw
|
---|
284 | gROOT->Reset();
|
---|
285 | gStyle->SetCanvasColor(0);
|
---|
286 | TCanvas* canvas = new TCanvas("canvas","pedphotcalc.C", 0, 100, 650, 800);
|
---|
287 | canvas->SetBorderMode(0); // Delete the Canvas' border line
|
---|
288 | canvas->cd();
|
---|
289 | canvas->SetBorderMode(0);
|
---|
290 |
|
---|
291 | canvas->Divide(2,4);
|
---|
292 |
|
---|
293 | // draw pedestal histo
|
---|
294 | canvas->cd(1);
|
---|
295 | gPad->cd();
|
---|
296 | gPad->SetBorderMode(0);
|
---|
297 |
|
---|
298 | pedhist->SetStats(kFALSE);
|
---|
299 | pedhist->GetXaxis()->SetTitle("Pixel SW Id");
|
---|
300 | pedhist->GetYaxis()->SetTitle("Pedestal (ADC counts)");
|
---|
301 | pedrmshist->SetStats(kFALSE);
|
---|
302 | pedrmshist->SetLineColor(2);
|
---|
303 | pedhist->Draw();
|
---|
304 | pedrmshist->Draw("same");
|
---|
305 |
|
---|
306 | TLegend* leg1 = new TLegend(.14,.68,.39,.88);
|
---|
307 | leg1->SetHeader("");
|
---|
308 | leg1->AddEntry(pedhist,"Pedestal","L");
|
---|
309 | leg1->AddEntry(pedrmshist,"Pedestal RMS","L");
|
---|
310 | leg1->SetFillColor(0);
|
---|
311 | leg1->SetLineColor(0);
|
---|
312 | leg1->SetBorderSize(0);
|
---|
313 | leg1->Draw();
|
---|
314 |
|
---|
315 | // draw pedestal distribution
|
---|
316 | canvas->cd(2);
|
---|
317 | gPad->cd();
|
---|
318 | gPad->SetBorderMode(0);
|
---|
319 | peddist->GetXaxis()->SetTitle("Pedestal (ADC counts)");
|
---|
320 | pedrmsdist->SetLineColor(2);
|
---|
321 | peddist->Draw();
|
---|
322 | pedrmsdist->Draw("same");
|
---|
323 |
|
---|
324 |
|
---|
325 | TLegend* leg2 = new TLegend(.14,.68,.39,.88);
|
---|
326 | leg2->SetHeader("");
|
---|
327 | leg2->AddEntry(pedhist,"Pedestal","L");
|
---|
328 | leg2->AddEntry(pedrmshist,"Pedestal RMS","L");
|
---|
329 | leg2->SetFillColor(0);
|
---|
330 | leg2->SetLineColor(0);
|
---|
331 | leg2->SetBorderSize(0);
|
---|
332 | leg2->Draw();
|
---|
333 |
|
---|
334 | // draw calibration histo
|
---|
335 | canvas->cd(3);
|
---|
336 | gPad->cd();
|
---|
337 | gPad->SetBorderMode(0);
|
---|
338 | calhist->GetXaxis()->SetTitle("Pixel SW Id");
|
---|
339 | calhist->SetMaximum(1);
|
---|
340 | calhist->SetMinimum(0);
|
---|
341 | calhist->GetYaxis()->SetTitle("Calibration factor (photon/ADC count)");
|
---|
342 | calhist->SetStats(kFALSE);
|
---|
343 | calhist->Draw();
|
---|
344 |
|
---|
345 | // draw calibration distribution
|
---|
346 | canvas->cd(4);
|
---|
347 | gPad->cd();
|
---|
348 | gPad->SetBorderMode(0);
|
---|
349 | caldist->GetXaxis()->SetTitle("Calibration factor (photon/ADC count)");
|
---|
350 | caldist->Draw();
|
---|
351 |
|
---|
352 | // draw pedestal signal histo
|
---|
353 | canvas->cd(5);
|
---|
354 | gPad->cd();
|
---|
355 | gPad->SetBorderMode(0);
|
---|
356 | pedphothist->GetXaxis()->SetTitle("Pixel SW Id");
|
---|
357 | pedphothist->GetYaxis()->SetTitle("Pedestal signal (photons)");
|
---|
358 | pedphothist->SetStats(kFALSE);
|
---|
359 | pedphothist->Draw();
|
---|
360 |
|
---|
361 | // draw pedestal signal distribution
|
---|
362 | canvas->cd(6);
|
---|
363 | gPad->cd();
|
---|
364 | gPad->SetBorderMode(0);
|
---|
365 | pedphotdist->GetXaxis()->SetTitle("Pedestal signal (photons)");
|
---|
366 | pedphotdist->Draw();
|
---|
367 |
|
---|
368 | // draw pedestal signal rms histo
|
---|
369 | canvas->cd(7);
|
---|
370 | gPad->cd();
|
---|
371 | gPad->SetBorderMode(0);
|
---|
372 | pedphotrmshist->GetXaxis()->SetTitle("Pixel SW Id");
|
---|
373 | pedphotrmshist->GetYaxis()->SetTitle("Pedestal signal rms (photons)");
|
---|
374 | pedphotrmshist->SetStats(kFALSE);
|
---|
375 | pedphotrmshist->Draw();
|
---|
376 |
|
---|
377 | // draw pedestal signal rms distribution
|
---|
378 | canvas->cd(8);
|
---|
379 | gPad->cd();
|
---|
380 | gPad->SetBorderMode(0);
|
---|
381 | pedphotrmsdist->GetXaxis()->SetTitle("Pedestal signal rms (photons)");
|
---|
382 | pedphotrmsdist->Draw();
|
---|
383 |
|
---|
384 | return;
|
---|
385 | }
|
---|
386 |
|
---|