source: trunk/MagicSoft/Mars/macros/pedphotcalc.C@ 2860

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