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 |
|
---|
25 | void calibration(TString pedname="/mnt/Data/rootdata/2003_12_01/20031130_03340_P_CrabNebula_E.root",
|
---|
26 | TString calname="/mnt/Data/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 | tlist2.AddToList(&geomapl);
|
---|
97 |
|
---|
98 | //
|
---|
99 | // Now setup the new tasks and tasklist for the calibration
|
---|
100 | // ---------------------------------------------------
|
---|
101 | //
|
---|
102 |
|
---|
103 | MReadMarsFile read2("Events", calname);
|
---|
104 | read2.DisableAutoScheme();
|
---|
105 |
|
---|
106 | MExtractSignal sigsig;
|
---|
107 | MCalibrationCalc calcalc;
|
---|
108 | MExtractedSignalCam sigcam;
|
---|
109 |
|
---|
110 | plist2.AddToList(&geomcam);
|
---|
111 | plist2.AddToList(&sigcam);
|
---|
112 |
|
---|
113 | //
|
---|
114 | // As long, as we don't have digital modules,
|
---|
115 | // we have to set the color of the pulser LED by hand
|
---|
116 | //
|
---|
117 | calcalc.SetPulserColor(MCalibrationCalc::kECT1);
|
---|
118 |
|
---|
119 | tlist2.AddToList(&read2);
|
---|
120 | tlist2.AddToList(&sigsig);
|
---|
121 | tlist2.AddToList(&calcalc);
|
---|
122 |
|
---|
123 | MHCamEvent hist2;
|
---|
124 | hist2.SetType(8);
|
---|
125 | plist2.AddToList(&hist2);
|
---|
126 | MFillH fill2("MHCamEvent", "MCalibrationCam");
|
---|
127 | tlist2.AddToList(&fill2);
|
---|
128 |
|
---|
129 | MStatusDisplay *d2 = new MStatusDisplay;
|
---|
130 | d2->SetUpdateTime(3000);
|
---|
131 |
|
---|
132 | //
|
---|
133 | // Create and setup the eventloop
|
---|
134 | //
|
---|
135 | MEvtLoop evtloop2;
|
---|
136 | evtloop2.SetParList(&plist2);
|
---|
137 | evtloop2.SetDisplay(d2);
|
---|
138 |
|
---|
139 | //
|
---|
140 | // Execute second analysis
|
---|
141 | //
|
---|
142 | if (!evtloop2.Eventloop())
|
---|
143 | return;
|
---|
144 |
|
---|
145 | tlist2.PrintStatistics();
|
---|
146 |
|
---|
147 | //
|
---|
148 | // just one example how to get the plots of individual pixels
|
---|
149 | //
|
---|
150 | MCalibrationCam *cam = plist2.FindObject("MCalibrationCam");
|
---|
151 | cam.Print();
|
---|
152 |
|
---|
153 | //
|
---|
154 | // Here we are confronted to a serious bug in ROOT:
|
---|
155 | // If we do not apply the next command, gPad will get
|
---|
156 | // screwed up completely: (Thanks to tbretz for finding out
|
---|
157 | // the reason during several hours!!!)
|
---|
158 | //
|
---|
159 | gROOT->GetListOfCanvases()->Delete();
|
---|
160 |
|
---|
161 | MHCamEvent &h = *(MHCamEvent*)plist2->FindObject("MHCamEvent");
|
---|
162 | MHCamera &disp0 = *h.GetHistByName();
|
---|
163 | MHCamera disp1 (geomcam, "MCalibrationCam;q", "Fitted Mean Charges");
|
---|
164 | MHCamera disp3 (geomcam, "MCalibrationCam;sigmaq", "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 | MHCamera disp16 (geomcam, "MCalibrationCam;sigsqchargesq", "Sigma Square per Charges Square");
|
---|
177 |
|
---|
178 | disp1.SetCamContent(*cam, 0);
|
---|
179 | disp1.SetCamError(*cam,1);
|
---|
180 |
|
---|
181 | disp3.SetCamContent(*cam, 2);
|
---|
182 | disp3.SetCamError(*cam,3);
|
---|
183 |
|
---|
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 | disp12.SetCamError(*cam, 12);
|
---|
199 |
|
---|
200 | disp13.SetCamContent(*cam, 13);
|
---|
201 | disp13.SetCamError(*cam, 14);
|
---|
202 |
|
---|
203 | disp14.SetCamContent(*cam, 15);
|
---|
204 | disp15.SetCamContent(*cam, 16);
|
---|
205 | disp16.SetCamContent(*cam, 17);
|
---|
206 |
|
---|
207 |
|
---|
208 | disp1.SetYTitle("Q [FADC counts]");
|
---|
209 | // disp2.SetYTitle("\\Delta Q [FADC counts]");
|
---|
210 | disp3.SetYTitle("\\sigma_{Q} [FADC counts]");
|
---|
211 | // disp4.SetYTitle("\\Delta {\\sigma_{Q}} [FADC counts]");
|
---|
212 | disp5.SetYTitle("P_{Q} [au]");
|
---|
213 | disp6.SetYTitle("T [FADC slices]");
|
---|
214 | disp7.SetYTitle("\\sigma_{T} [FADC slices]");
|
---|
215 | disp8.SetYTitle("P_{T} [au]");
|
---|
216 | disp9.SetYTitle("P [Total FADC counts ]");
|
---|
217 | disp10.SetYTitle("RMS_{P} [Total FADC counts ]");
|
---|
218 | disp11.SetYTitle("\\sigma^{2}_{Q} - RMS^{2}_{P} [FADC counts^2]");
|
---|
219 | disp12.SetYTitle("Nr Phe's");
|
---|
220 | disp13.SetYTitle("Conversion Factor [Phe/FADC count]");
|
---|
221 | disp14.SetYTitle("Nr Photons");
|
---|
222 | disp15.SetYTitle("Conversion Factor [Ph/FADC count]");
|
---|
223 | disp16.SetYTitle("Sigma^2 per Charge^2 [1]");
|
---|
224 |
|
---|
225 | MStatusDisplay *d3 = new MStatusDisplay;
|
---|
226 | d3->SetUpdateTime(3000);
|
---|
227 | d3->Resize(1180,900);
|
---|
228 |
|
---|
229 | gStyle->SetOptStat(1111);
|
---|
230 | gStyle->SetOptFit();
|
---|
231 |
|
---|
232 | // Charges
|
---|
233 | TCanvas *c1 = &d3->AddTab("Fitted Charges");
|
---|
234 | c1->Divide(2,3);
|
---|
235 |
|
---|
236 | TObject *obj1;
|
---|
237 | TObject *obj2;
|
---|
238 |
|
---|
239 | c1->cd(1);
|
---|
240 | obj1=disp1.DrawCopy("hist");
|
---|
241 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
242 |
|
---|
243 | c1->cd(3);
|
---|
244 | gPad->SetBorderMode(0);
|
---|
245 | obj1->Draw();
|
---|
246 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
247 |
|
---|
248 | c1->cd(5);
|
---|
249 | gPad->SetBorderMode(0);
|
---|
250 | obj2 = obj1->DrawClone("proj");
|
---|
251 | gPad->Update();
|
---|
252 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
253 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
254 |
|
---|
255 | c1->cd(2);
|
---|
256 | obj1=disp3.DrawCopy("hist");
|
---|
257 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
258 |
|
---|
259 | c1->cd(4);
|
---|
260 | gPad->SetBorderMode(0);
|
---|
261 | obj1->Draw();
|
---|
262 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
263 |
|
---|
264 | c1->cd(6);
|
---|
265 | gPad->SetBorderMode(0);
|
---|
266 | obj2 = obj1->DrawClone("proj");
|
---|
267 | gPad->Update();
|
---|
268 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
269 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
270 |
|
---|
271 | // Fit Probability
|
---|
272 | TCanvas *c12 = &d3->AddTab("Fit Prob.");
|
---|
273 | c12->Divide(1, 3);
|
---|
274 |
|
---|
275 | c12->cd(1);
|
---|
276 | obj1=disp5.DrawCopy("hist");
|
---|
277 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
278 |
|
---|
279 | c12->cd(2);
|
---|
280 | gPad->SetBorderMode(0);
|
---|
281 | obj1->Draw();
|
---|
282 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
283 |
|
---|
284 | c12->cd(3);
|
---|
285 | gPad->SetBorderMode(0);
|
---|
286 | obj2 = obj1->DrawClone("proj");
|
---|
287 | gPad->Update();
|
---|
288 | ((MHCamera*)obj2)->GetYProj()->Fit("pol0","Q");
|
---|
289 | ((MHCamera*)obj2)->GetYProj()->GetFunction("pol0")->SetLineColor(kYellow);
|
---|
290 |
|
---|
291 | // Times
|
---|
292 | TCanvas *c2 = &d3->AddTab("Fitted Times");
|
---|
293 | c2->Divide(3, 3);
|
---|
294 |
|
---|
295 | c2->cd(1);
|
---|
296 | obj1=disp6.DrawCopy("hist");
|
---|
297 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
298 |
|
---|
299 | c2->cd(4);
|
---|
300 | obj1->Draw();
|
---|
301 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
302 |
|
---|
303 | c2->cd(7);
|
---|
304 | gPad->SetBorderMode(0);
|
---|
305 | obj2 = obj1->DrawClone("proj");
|
---|
306 | gPad->Update();
|
---|
307 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
308 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
309 |
|
---|
310 | c2->cd(2);
|
---|
311 | obj1=disp7.DrawCopy("hist");
|
---|
312 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
313 |
|
---|
314 | c2->cd(5);
|
---|
315 | obj1->Draw();
|
---|
316 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
317 |
|
---|
318 | c2->cd(8);
|
---|
319 | gPad->SetBorderMode(0);
|
---|
320 | obj2 = obj1->DrawClone("proj");
|
---|
321 | gPad->Update();
|
---|
322 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
323 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
324 |
|
---|
325 | c2->cd(3);
|
---|
326 | obj1=disp8.DrawCopy("hist");
|
---|
327 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
328 |
|
---|
329 | c2->cd(6);
|
---|
330 | obj1->Draw();
|
---|
331 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
332 |
|
---|
333 | c2->cd(9);
|
---|
334 | gPad->SetBorderMode(0);
|
---|
335 | obj2 = obj1->DrawClone("proj");
|
---|
336 | gPad->Update();
|
---|
337 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
338 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
339 |
|
---|
340 | // Pedestals
|
---|
341 | TCanvas *c3 = &d3->AddTab("Pedestals");
|
---|
342 | c3->Divide(2, 3);
|
---|
343 |
|
---|
344 | c3->cd(1);
|
---|
345 | obj1=disp9.DrawCopy("hist");
|
---|
346 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
347 |
|
---|
348 | c3->cd(3);
|
---|
349 | obj1->Draw();
|
---|
350 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
351 |
|
---|
352 | c3->cd(5);
|
---|
353 | gPad->SetBorderMode(0);
|
---|
354 | obj2 = obj1->DrawClone("proj");
|
---|
355 | gPad->Update();
|
---|
356 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
357 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
358 |
|
---|
359 | c3->cd(2);
|
---|
360 | obj1=disp10.DrawCopy("hist");
|
---|
361 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
362 |
|
---|
363 | c3->cd(4);
|
---|
364 | obj1->Draw();
|
---|
365 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
366 |
|
---|
367 | c3->cd(6);
|
---|
368 | gPad->SetBorderMode(0);
|
---|
369 | obj2 = obj1->DrawClone("proj");
|
---|
370 | gPad->Update();
|
---|
371 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
372 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
373 |
|
---|
374 | // Reduced Sigmas
|
---|
375 | TCanvas *c4 = &d3->AddTab("Reduced Sigmas");
|
---|
376 | c4->Divide(2,3);
|
---|
377 |
|
---|
378 | c4->cd(1);
|
---|
379 | obj1=disp11.DrawCopy("hist");
|
---|
380 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
381 |
|
---|
382 | c4->cd(3);
|
---|
383 | obj1->Draw();
|
---|
384 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
385 |
|
---|
386 | c4->cd(5);
|
---|
387 | gPad->SetBorderMode(0);
|
---|
388 | obj2 = obj1->DrawClone("proj");
|
---|
389 | gPad->Update();
|
---|
390 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
391 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
392 |
|
---|
393 | c4->cd(2);
|
---|
394 | obj1=disp16.DrawCopy("hist");
|
---|
395 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
396 |
|
---|
397 | c4->cd(4);
|
---|
398 | obj1->Draw();
|
---|
399 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
400 |
|
---|
401 | c4->cd(6);
|
---|
402 | gPad->SetBorderMode(0);
|
---|
403 | obj2 = obj1->DrawClone("proj");
|
---|
404 | gPad->Update();
|
---|
405 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
406 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
407 |
|
---|
408 | // F-Factor Method
|
---|
409 | TCanvas *c5 = &d3->AddTab("F-Factor Method");
|
---|
410 | c5->Divide(2, 3);
|
---|
411 |
|
---|
412 | c5->cd(1);
|
---|
413 | obj1=disp12.DrawCopy("hist");
|
---|
414 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
415 |
|
---|
416 | c5->cd(3);
|
---|
417 | obj1->Draw();
|
---|
418 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
419 |
|
---|
420 | c5->cd(5);
|
---|
421 | gPad->SetBorderMode(0);
|
---|
422 | obj2 = obj1->DrawClone("proj");
|
---|
423 | gPad->Update();
|
---|
424 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
425 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
426 |
|
---|
427 | c5->cd(2);
|
---|
428 | obj1=disp13.DrawCopy("hist");
|
---|
429 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
430 |
|
---|
431 | c5->cd(4);
|
---|
432 | obj1->Draw();
|
---|
433 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
434 |
|
---|
435 | c5->cd(6);
|
---|
436 | gPad->SetBorderMode(0);
|
---|
437 | obj2 = obj1->DrawClone("proj");
|
---|
438 | gPad->Update();
|
---|
439 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
440 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
441 |
|
---|
442 | // Blind Pixel Method
|
---|
443 | TCanvas *c6 = &d3->AddTab("Blind Pixel Method");
|
---|
444 | c6->Divide(2, 3);
|
---|
445 |
|
---|
446 | c6->cd(1);
|
---|
447 | obj1=disp14.DrawCopy("hist");
|
---|
448 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
449 |
|
---|
450 | c6->cd(3);
|
---|
451 | obj1->Draw();
|
---|
452 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
453 |
|
---|
454 |
|
---|
455 | c6->cd(2);
|
---|
456 | obj1=disp15.DrawCopy("hist");
|
---|
457 | ((MHCamera*)obj1)->AddNotify(*cam);
|
---|
458 |
|
---|
459 | c6->cd(4);
|
---|
460 | obj1->Draw();
|
---|
461 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
462 |
|
---|
463 | c6->cd(6);
|
---|
464 | gPad->SetBorderMode(0);
|
---|
465 | obj2 = obj1->DrawClone("proj");
|
---|
466 | gPad->Update();
|
---|
467 | ((MHCamera*)obj2)->GetYProj()->Fit("gaus","Q");
|
---|
468 | ((MHCamera*)obj2)->GetYProj()->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
469 | }
|
---|
470 |
|
---|