source: trunk/MagicSoft/Mars/macros/calibration.C@ 2629

Last change on this file since 2629 was 2629, checked in by gaug, 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!
18! Author(s): Markus Gaug, 11/2003 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25void calibration(TString pedname="/data/MAGIC/rootdata/2003_12_01/20031130_03340_P_CrabNebula_E.root",
26 TString calname="/data/MAGIC/rootdata/2003_12_01/20031130_03341_C_CrabNebula_E.root")
27{
28
29 //
30 // Create a empty Parameter List and an empty Task List
31 // The tasklist is identified in the eventloop by its name
32 //
33 MParList plist;
34
35 MTaskList tlist;
36 plist.AddToList(&tlist);
37
38 //
39 // Now setup the tasks and tasklist for the pedestals:
40 // ---------------------------------------------------
41 //
42
43 MReadMarsFile read("Events", pedname);
44 read.DisableAutoScheme();
45
46 MGeomApply geomapl;
47 MPedCalcPedRun pedcalc;
48 MGeomCamMagic geomcam;
49 MPedestalCam pedcam;
50
51 tlist.AddToList(&read);
52 tlist.AddToList(&geomapl);
53 tlist.AddToList(&pedcalc);
54
55 plist.AddToList(&pedcam);
56
57 MHCamEvent hist("Pedestal");
58 hist.SetType(0);
59 plist.AddToList(&hist);
60 MFillH fill(&hist, "MPedestalCam");
61
62 tlist.AddToList(&fill);
63
64 MStatusDisplay *d1 = new MStatusDisplay;
65
66 // Set update time to 3s
67 d1->SetUpdateTime(3000);
68
69 //
70 // Create and setup the eventloop
71 //
72 MEvtLoop evtloop;
73 evtloop.SetParList(&plist);
74 evtloop.SetDisplay(d1);
75
76 //
77 // Execute first analysis
78 //
79 if (!evtloop.Eventloop())
80 return;
81
82 tlist.PrintStatistics();
83
84
85 //
86 // Create a empty Parameter List and an empty Task List
87 //
88 MParList plist2;
89
90 MTaskList tlist2;
91 plist2.AddToList(&tlist2);
92
93
94 plist2.AddToList((MPedestalCam*)plist.FindObject("MPedestalCam"));
95
96// MGeomApply geomapl2;
97 tlist2.AddToList(&geomapl);
98
99 //
100 // Now setup the new tasks and tasklist for the calibration
101 // ---------------------------------------------------
102 //
103
104 MReadMarsFile read2("Events", calname);
105 read2.DisableAutoScheme();
106
107 MExtractSignal sigsig;
108 MCalibrationCalc calcalc;
109// calcalc.SetSkipTFits();
110 MExtractedSignalCam sigcam;
111
112 plist2.AddToList(&geomcam);
113 plist2.AddToList(&sigcam);
114
115 //
116 // As long, as we don't have digital modules,
117 // we have to set the color of the pulser LED by hand
118 //
119 calcalc.SetPulserColor(MCalibrationCalc::kECT1);
120
121 tlist2.AddToList(&read2);
122 tlist2.AddToList(&calcalc);
123 tlist2.AddToList(&sigsig);
124
125 MHCamEvent hist2;
126 hist2.SetType(8);
127 plist2.AddToList(&hist2);
128 MFillH fill2("MHCamEvent", "MCalibrationCam");
129 tlist2.AddToList(&fill2);
130
131 //
132 // Create and setup the eventloop
133 //
134 MEvtLoop evtloop2;
135 evtloop2.SetParList(&plist2);
136
137 //
138 // Execute second analysis
139 //
140 if (!evtloop2.Eventloop())
141 return;
142
143 tlist2.PrintStatistics();
144
145 //
146 // just one example how to get the plots of individual pixels
147 //
148 MCalibrationCam *cam = plist2.FindObject("MCalibrationCam");
149 cam.Print();
150
151 //
152 // Here we are confronted to a serious bug in ROOT:
153 // If we do not apply the next command, gPad will get
154 // screwed up completely: (Thanks to tbretz for finding out
155 // the reason during several hours!!!)
156 //
157 gROOT->GetListOfCanvases()->Delete();
158
159 MHCamEvent &h = *(MHCamEvent*)plist2->FindObject("MHCamEvent");
160 MHCamera &disp0 = *h.GetHistByName();
161 MHCamera disp1 (geomcam, "MCalibrationCam;q", "Fitted Mean Charges");
162 // MHCamera disp2 (geomcam, "MCalibrationCam;errq", "Error of Fitted Mean Charges");
163 MHCamera disp3 (geomcam, "MCalibrationCam;sigmaq", "Sigma of Fitted Mean Charges");
164 // MHCamera disp4 (geomcam, "MCalibrationCam;errsigmaq", "Error of Sigma of Fitted Mean Charges");
165 MHCamera disp5 (geomcam, "MCalibrationCam;probq", "Probability of Fit");
166 MHCamera disp6 (geomcam, "MCalibrationCam;t", "Arrival Times");
167 MHCamera disp7 (geomcam, "MCalibrationCam;sigmat", "Sigma of Arrival Times");
168 MHCamera disp8 (geomcam, "MCalibrationCam;probt", "Probability of Time Fit");
169 MHCamera disp9 (geomcam, "MCalibrationCam;ped", "Pedestals");
170 MHCamera disp10 (geomcam, "MCalibrationCam;pedrms", "Pedestal RMS");
171 MHCamera disp11 (geomcam, "MCalibrationCam;rsigma", "Reduced Sigmas");
172 MHCamera disp12 (geomcam, "MCalibrationCam;phe", "Nr. of Phe's (F-Factor Method)");
173 MHCamera disp13 (geomcam, "MCalibrationCam;convphe", "Conversion Factor (F-Factor Method)");
174 MHCamera disp14 (geomcam, "MCalibrationCam;photons", "Nr. of Photons (Blind Pixel Method)");
175 MHCamera disp15 (geomcam, "MCalibrationCam;convphot", "Conversion Factor (Blind Pixel Method)");
176
177 disp1.SetCamContent(*cam, 0);
178 disp1.SetCamError(*cam,1);
179 // disp2.SetCamContent(*cam, 1);
180
181 disp3.SetCamContent(*cam, 2);
182 disp3.SetCamError(*cam,3);
183 // disp4.SetCamContent(*cam, 3);
184 disp5.SetCamContent(*cam, 4);
185
186 disp6.SetCamContent(*cam, 5);
187 disp6.SetCamError(*cam,6);
188 disp7.SetCamContent(*cam, 6);
189 disp8.SetCamContent(*cam, 7);
190
191 disp9.SetCamContent(*cam, 8);
192 disp9.SetCamError(*cam,9);
193 disp10.SetCamContent(*cam, 9);
194
195 disp11.SetCamContent(*cam, 10);
196
197 disp12.SetCamContent(*cam, 11);
198 disp13.SetCamContent(*cam, 12);
199 disp14.SetCamContent(*cam, 13);
200 disp15.SetCamContent(*cam, 14);
201 // disp16.SetCamError(*cam, 16);
202
203
204 disp1.SetYTitle("Q [FADC counts]");
205 // disp2.SetYTitle("\\Delta Q [FADC counts]");
206 disp3.SetYTitle("\\sigma_{Q} [FADC counts]");
207 // disp4.SetYTitle("\\Delta {\\sigma_{Q}} [FADC counts]");
208 disp5.SetYTitle("P [au]");
209 disp6.SetYTitle("T [FADC slices]");
210 disp7.SetYTitle("\\Delta T [FADC slices]");
211 disp8.SetYTitle("P [au]");
212 disp9.SetYTitle("P [FADC counts/ slice ]");
213 disp10.SetYTitle("RMS_{P} [FADC counts / slice ]");
214 disp11.SetYTitle("\\sigma^2_{Q} - RMS^2_{P} [FADC counts^2]");
215 disp12.SetYTitle("Nr Phe's");
216 disp13.SetYTitle("Conversion Factor [Phe/FADC count]");
217 disp14.SetYTitle("Nr Photons");
218 disp15.SetYTitle("Conversion Factor [Ph/FADC count]");
219
220 disp1.SetPrettyPalette();
221 disp3.SetPrettyPalette();
222 disp5.SetPrettyPalette();
223 disp6.SetPrettyPalette();
224 disp7.SetPrettyPalette();
225 disp8.SetPrettyPalette();
226 disp9.SetPrettyPalette();
227 disp10.SetPrettyPalette();
228 disp11.SetPrettyPalette();
229 disp12.SetPrettyPalette();
230 disp13.SetPrettyPalette();
231 disp14.SetPrettyPalette();
232 disp15.SetPrettyPalette();
233
234
235
236 MStatusDisplay *d2 = new MStatusDisplay;
237
238 // Set update time to 1s
239 d2->SetUpdateTime(1000);
240
241 TCanvas *c1 = &d2->AddTab("Fitted Charges");
242 c1->Divide(2, 2);
243
244 TObject *obj;
245
246 c1->cd(1);
247 gStyle->SetOptStat(1111);
248 obj=disp1.DrawCopy("hist");
249 ((MHCamera*)obj)->AddNotify(*cam);
250
251 c1->cd(3);
252 gPad->SetBorderMode(0);
253 obj->Draw();
254
255 c1->cd(2);
256 gStyle->SetOptStat(1101);
257 obj=disp3.DrawCopy("hist");
258 ((MHCamera*)obj)->AddNotify(*cam);
259
260 c1->cd(4);
261 gPad->SetBorderMode(0);
262 obj->Draw();
263
264
265 TCanvas *c12 = &d2->AddTab("Fit Prob.");
266 c12->Divide(1, 2);
267
268 c12->cd(1);
269 gStyle->SetOptStat(1101);
270 obj=disp5.DrawCopy("hist");
271 ((MHCamera*)obj)->AddNotify(*cam);
272
273 c12->cd(2);
274 gPad->SetBorderMode(0);
275 obj->Draw();
276
277 TCanvas *c2 = &d2->AddTab("Fitted Times");
278 c2->Divide(3, 2);
279
280 c2->cd(1);
281 gStyle->SetOptStat(1111);
282 obj=disp6.DrawCopy("hist");
283 ((MHCamera*)obj)->AddNotify(*cam);
284
285 c2->cd(4);
286 obj->Draw();
287
288 c2->cd(2);
289 gStyle->SetOptStat(1101);
290 obj=disp7.DrawCopy("hist");
291 ((MHCamera*)obj)->AddNotify(*cam);
292
293 c2->cd(5);
294 obj->Draw();
295
296 c2->cd(3);
297 gStyle->SetOptStat(1101);
298 obj=disp8.DrawCopy("hist");
299 ((MHCamera*)obj)->AddNotify(*cam);
300
301 c2->cd(6);
302 obj->Draw();
303
304 TCanvas *c3 = &d2->AddTab("Pedestals");
305 c3->Divide(2, 2);
306
307 c3->cd(1);
308 gStyle->SetOptStat(1111);
309 obj=disp9.DrawCopy("hist");
310 ((MHCamera*)obj)->AddNotify(*cam);
311
312 c3->cd(3);
313 obj->Draw();
314
315 c3->cd(2);
316 gStyle->SetOptStat(1111);
317 obj=disp10.DrawCopy("hist");
318 ((MHCamera*)obj)->AddNotify(*cam);
319
320 c3->cd(4);
321 obj->Draw();
322
323 TCanvas *c4 = &d2->AddTab("Reduced Charges");
324 c4->Divide(2,1);
325
326 c4->cd(1);
327 gStyle->SetOptStat(1111);
328 obj=disp11.DrawCopy("hist");
329 ((MHCamera*)obj)->AddNotify(*cam);
330
331 c4->cd(3);
332 obj->Draw();
333
334 TCanvas *c5 = &d2->AddTab("F-Factor Method");
335 c5->Divide(2, 2);
336
337 c5->cd(1);
338 gStyle->SetOptStat(1111);
339 obj=disp12.DrawCopy("hist");
340 ((MHCamera*)obj)->AddNotify(*cam);
341
342 c5->cd(3);
343 obj->Draw();
344
345 c5->cd(2);
346 gStyle->SetOptStat(1101);
347 obj=disp13.DrawCopy("hist");
348 ((MHCamera*)obj)->AddNotify(*cam);
349
350 c5->cd(4);
351 obj->Draw();
352
353 TCanvas *c6 = &d2->AddTab("Blind Pixel Method");
354 c6->Divide(2, 2);
355
356 c6->cd(1);
357 gStyle->SetOptStat(1111);
358 obj=disp14.DrawCopy("hist");
359 ((MHCamera*)obj)->AddNotify(*cam);
360
361 c6->cd(3);
362 obj->Draw();
363
364 c6->cd(2);
365 gStyle->SetOptStat(1101);
366 obj=disp15.DrawCopy("hist");
367 ((MHCamera*)obj)->AddNotify(*cam);
368
369 c6->cd(4);
370 obj->Draw();
371
372}
373
Note: See TracBrowser for help on using the repository browser.