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/Miscellaneous/2003_12_19/20031218_03522_P_Park_E.root",
|
---|
26 | TString calname="/mnt/Data/rootdata/Miscellaneous/2003_12_19/20031218_03527_C_Park_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(1);
|
---|
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 | // Create a empty Parameter List and an empty Task List
|
---|
86 | //
|
---|
87 | MParList plist2;
|
---|
88 | MTaskList tlist2;
|
---|
89 | plist2.AddToList(&tlist2);
|
---|
90 |
|
---|
91 | MExtractedSignalCam sigcam;
|
---|
92 | MCalibrationCam calcam;
|
---|
93 | //
|
---|
94 | // Get the previously created MPedestalCam into the new Parameter List
|
---|
95 | //
|
---|
96 | plist2.AddToList(&pedcam);
|
---|
97 | plist2.AddToList(&sigcam);
|
---|
98 | plist2.AddToList(&calcam);
|
---|
99 |
|
---|
100 | //
|
---|
101 | // Get the MAGIC geometry
|
---|
102 | //
|
---|
103 | tlist2.AddToList(&geomapl);
|
---|
104 | plist2.AddToList(&geomcam);
|
---|
105 |
|
---|
106 | //
|
---|
107 | // Now setup the new tasks and tasklist for the calibration
|
---|
108 | // ---------------------------------------------------
|
---|
109 | //
|
---|
110 |
|
---|
111 | MReadMarsFile read2("Events", calname);
|
---|
112 | read2.DisableAutoScheme();
|
---|
113 |
|
---|
114 | MExtractSignal sigsig;
|
---|
115 | MCalibrationCalc calcalc;
|
---|
116 |
|
---|
117 | //
|
---|
118 | // In case, we want to exclude a pre-defined list of bad pixels:
|
---|
119 | // (This is a preliminary feature)
|
---|
120 | //
|
---|
121 | // calcalc.ExcludePixelsFromAsciiFile("badpixels.dat");
|
---|
122 |
|
---|
123 | //
|
---|
124 | // As long, as we don't have digital modules,
|
---|
125 | // we have to set the color of the pulser LED by hand
|
---|
126 | //
|
---|
127 | calcalc.SetPulserColor(MCalibrationCalc::kECT1);
|
---|
128 |
|
---|
129 | //
|
---|
130 | // In case, we want to apply another fit function to the
|
---|
131 | // blind pixel
|
---|
132 | //
|
---|
133 | MCalibrationBlindPix *bp = calcam.GetBlindPixel();
|
---|
134 | bp->ChangeFitFunc(MHCalibrationBlindPixel::kEPoisson5);
|
---|
135 |
|
---|
136 | tlist2.AddToList(&read2);
|
---|
137 | tlist2.AddToList(&sigsig);
|
---|
138 | tlist2.AddToList(&calcalc);
|
---|
139 |
|
---|
140 | //
|
---|
141 | // Create and setup the eventloop
|
---|
142 | //
|
---|
143 | MEvtLoop evtloop2;
|
---|
144 | evtloop2.SetParList(&plist2);
|
---|
145 |
|
---|
146 | //
|
---|
147 | // Execute second analysis
|
---|
148 | //
|
---|
149 | if (!evtloop2.Eventloop())
|
---|
150 | return;
|
---|
151 |
|
---|
152 | tlist2.PrintStatistics();
|
---|
153 |
|
---|
154 | //
|
---|
155 | // just one example how to get the plots of individual pixels
|
---|
156 | //
|
---|
157 | MCalibrationCam *cam = plist2.FindObject("MCalibrationCam");
|
---|
158 | cam.Print();
|
---|
159 |
|
---|
160 | MHCamEvent camevt;
|
---|
161 | MHCamera disp1 (geomcam, "MCalibrationPix;Charge", "Fitted Mean Charges");
|
---|
162 | MHCamera disp3 (geomcam, "MCalibrationPix;SigmaCharge", "Sigma of Fitted Charges");
|
---|
163 | MHCamera disp5 (geomcam, "MCalibrationPix;ChargeProb", "Probability of Fit");
|
---|
164 | MHCamera disp6 (geomcam, "MCalibrationPix;Time", "Arrival Times");
|
---|
165 | MHCamera disp7 (geomcam, "MCalibrationPix;SigmaTime", "Sigma of Arrival Times");
|
---|
166 | MHCamera disp8 (geomcam, "MCalibrationPix;TimeChiSquare", "Chi Square of Time Fit");
|
---|
167 | MHCamera disp9 (geomcam, "MCalibrationPix;Ped", "Pedestals");
|
---|
168 | MHCamera disp10 (geomcam, "MCalibrationPix;PedRms", "Pedestal RMS");
|
---|
169 | MHCamera disp11 (geomcam, "MCalibrationPix;RSigma", "Reduced Sigmas");
|
---|
170 | MHCamera disp12 (geomcam, "MCalibrationPix;PheFFactorMethod", "Nr. of Phe's (F-Factor Method)");
|
---|
171 | MHCamera disp13 (geomcam, "MCalibrationPix;MeanConversionFFactorMethod",
|
---|
172 | "Conversion Factor (F-Factor Method)");
|
---|
173 | MHCamera disp14 (geomcam, "MCalibrationPix;MeanPhotInsidePlexiglass",
|
---|
174 | "Nr. of Photons (Blind Pixel Method)");
|
---|
175 | MHCamera disp15 (geomcam, "MCalibrationPix;MeanConversionBlindPixelMethod",
|
---|
176 | "Conversion Factor (Blind Pixel Method)");
|
---|
177 | MHCamera disp16 (geomcam, "MCalibrationPix;RSigma/Charge", "Reduced Sigma per Charge");
|
---|
178 |
|
---|
179 | disp1.SetCamContent(*cam, 0);
|
---|
180 | disp1.SetCamError(*cam,1);
|
---|
181 |
|
---|
182 | disp3.SetCamContent(*cam, 2);
|
---|
183 | disp3.SetCamError(*cam,3);
|
---|
184 |
|
---|
185 | disp5.SetCamContent(*cam, 4);
|
---|
186 |
|
---|
187 | disp6.SetCamContent(*cam, 5);
|
---|
188 | disp6.SetCamError(*cam, 6);
|
---|
189 | disp7.SetCamContent(*cam, 6);
|
---|
190 | disp8.SetCamContent(*cam, 7);
|
---|
191 |
|
---|
192 | disp9.SetCamContent(*cam, 8);
|
---|
193 | disp9.SetCamError(*cam, 9);
|
---|
194 |
|
---|
195 | disp10.SetCamContent(*cam, 9);
|
---|
196 | disp11.SetCamContent(*cam, 10);
|
---|
197 |
|
---|
198 | disp12.SetCamContent(*cam, 11);
|
---|
199 | disp12.SetCamError(*cam, 12);
|
---|
200 |
|
---|
201 | disp13.SetCamContent(*cam, 13);
|
---|
202 | disp13.SetCamError(*cam, 14);
|
---|
203 |
|
---|
204 | disp14.SetCamContent(*cam, 15);
|
---|
205 | disp15.SetCamContent(*cam, 16);
|
---|
206 | disp16.SetCamContent(*cam, 17);
|
---|
207 |
|
---|
208 |
|
---|
209 | disp1.SetYTitle("Charge [FADC counts]");
|
---|
210 | disp3.SetYTitle("\\sigma_{Charge} [FADC counts]");
|
---|
211 | disp5.SetYTitle("P_{Charge} [1]");
|
---|
212 | disp6.SetYTitle("Arr. Time [Time Slice Nr.]");
|
---|
213 | disp7.SetYTitle("\\sigma_{Time} [Time Slices]");
|
---|
214 | disp8.SetYTitle("\\chi^{2}_{Time} [1]");
|
---|
215 | disp9.SetYTitle("Ped [FADC Counts ]");
|
---|
216 | disp10.SetYTitle("RMS_{Ped} [FADC Counts ]");
|
---|
217 | disp11.SetYTitle("\\sqrt{\\sigma^{2}_{Charge} - RMS^{2}_{Ped}} [FADC Counts]");
|
---|
218 | disp12.SetYTitle("Nr. Photo-Electrons [1]");
|
---|
219 | disp13.SetYTitle("Conversion Factor [PhE/FADC Count]");
|
---|
220 | disp14.SetYTitle("Nr. Photons [1]");
|
---|
221 | disp15.SetYTitle("Conversion Factor [Phot/FADC Count]");
|
---|
222 | disp16.SetYTitle("Reduced Sigma / Charge [1]");
|
---|
223 |
|
---|
224 | MStatusDisplay *d3 = new MStatusDisplay;
|
---|
225 | d3->SetUpdateTime(3000);
|
---|
226 | d3->Resize(1180,900);
|
---|
227 |
|
---|
228 | gStyle->SetOptStat(1111);
|
---|
229 | gStyle->SetOptFit();
|
---|
230 |
|
---|
231 | // Charges
|
---|
232 | TCanvas &c1 = d3->AddTab("Fitted Charges");
|
---|
233 | c1.Divide(2,3);
|
---|
234 |
|
---|
235 | CamDraw(c1,disp1,cam,1,2,1);
|
---|
236 | CamDraw(c1,disp3,cam,2,2,2);
|
---|
237 |
|
---|
238 | // Fit Probability
|
---|
239 | TCanvas &c2 = d3->AddTab("Fit Prob.");
|
---|
240 | c2.Divide(1,3);
|
---|
241 |
|
---|
242 | CamDraw(c2,disp5,cam,1,1,3);
|
---|
243 |
|
---|
244 | // Times
|
---|
245 | TCanvas &c3 = d3->AddTab("Fitted Times");
|
---|
246 | c3.Divide(3,3);
|
---|
247 |
|
---|
248 | CamDraw(c3,disp6,cam,1,3,1);
|
---|
249 | CamDraw(c3,disp7,cam,2,3,0);
|
---|
250 | CamDraw(c3,disp8,cam,3,3,0);
|
---|
251 |
|
---|
252 | // Pedestals
|
---|
253 | TCanvas &c4 = d3->AddTab("Pedestals");
|
---|
254 | c4.Divide(2,3);
|
---|
255 |
|
---|
256 | CamDraw(c4,disp9,cam,1,2,0);
|
---|
257 | CamDraw(c4,disp10,cam,2,2,1);
|
---|
258 |
|
---|
259 | // Reduced Sigmas
|
---|
260 | TCanvas &c5 = d3->AddTab("Reduced Sigmas");
|
---|
261 | c5.Divide(2,3);
|
---|
262 |
|
---|
263 | CamDraw(c5,disp11,cam,1,2,2);
|
---|
264 | CamDraw(c5,disp16,cam,2,2,2);
|
---|
265 |
|
---|
266 | // F-Factor Method
|
---|
267 | TCanvas &c6 = d3->AddTab("F-Factor Method");
|
---|
268 | c6.Divide(2,3);
|
---|
269 |
|
---|
270 | CamDraw(c6,disp12,cam,1,2,1);
|
---|
271 | CamDraw(c6,disp13,cam,2,2,1);
|
---|
272 |
|
---|
273 | // Blind Pixel Method
|
---|
274 | TCanvas &c7 = d3->AddTab("Blind Pixel Method");
|
---|
275 | c7.Divide(2, 3);
|
---|
276 |
|
---|
277 | CamDraw(c7,disp14,cam,1,2,9);
|
---|
278 | CamDraw(c7,disp15,cam,2,2,1);
|
---|
279 |
|
---|
280 | }
|
---|
281 |
|
---|
282 | void CamDraw(TCanvas &c, MHCamera &cam, MCamEvent *evt, Int_t i, Int_t j, Int_t fit)
|
---|
283 | {
|
---|
284 |
|
---|
285 | c.cd(i);
|
---|
286 | gPad->SetBorderMode(0);
|
---|
287 | MHCamera *obj1=(MHCamera*)cam.DrawCopy("hist");
|
---|
288 | obj1->AddNotify(*evt);
|
---|
289 |
|
---|
290 | c.cd(i+j);
|
---|
291 | gPad->SetBorderMode(0);
|
---|
292 | obj1->Draw();
|
---|
293 | ((MHCamera*)obj1)->SetPrettyPalette();
|
---|
294 |
|
---|
295 | c.cd(i+2*j);
|
---|
296 | gPad->SetBorderMode(0);
|
---|
297 | Float_t he = gStyle->GetStatH();
|
---|
298 | Float_t wi = gStyle->GetStatH();
|
---|
299 | gStyle->SetStatH(0.4);
|
---|
300 | gStyle->SetStatW(0.25);
|
---|
301 |
|
---|
302 | TH1D *obj2 = (TH1D*)obj1->Projection();
|
---|
303 | obj2->Draw();
|
---|
304 | obj2->SetBit(kCanDelete);
|
---|
305 |
|
---|
306 | const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
|
---|
307 | const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
|
---|
308 | const Double_t integ = obj2->Integral("width")/2.5;
|
---|
309 | const Double_t mean = obj2->GetMean();
|
---|
310 | const Double_t rms = obj2->GetRMS();
|
---|
311 | const Double_t width = max-min;
|
---|
312 |
|
---|
313 | switch (fit)
|
---|
314 | {
|
---|
315 | case 0:
|
---|
316 | TF1 *sgaus = new TF1("sgaus","gaus(0)",min,max);
|
---|
317 | sgaus->SetParNames("Area","#mu","#sigma");
|
---|
318 | sgaus->SetParameters(integ/rms,mean,rms);
|
---|
319 | sgaus->SetParLimits(0,0.,integ);
|
---|
320 | sgaus->SetParLimits(1,min,max);
|
---|
321 | sgaus->SetParLimits(2,0,width/1.5);
|
---|
322 | obj2->Fit("sgaus","QLR");
|
---|
323 | obj2->GetFunction("sgaus")->SetLineColor(kYellow);
|
---|
324 | break;
|
---|
325 |
|
---|
326 | case 1:
|
---|
327 | TF1 *dgaus = new TF1("dgaus","gaus(0)+gaus(3)",min,max);
|
---|
328 | dgaus->SetParNames("A1","#mu1","#sigma1","A2","#mu2","#sigma2");
|
---|
329 | dgaus->SetParameters(integ/width,max-width/6.,width/4.,
|
---|
330 | integ/width,min+width/6.,width/4.);
|
---|
331 | dgaus->SetParLimits(0,0,integ);
|
---|
332 | dgaus->SetParLimits(1,min,max);
|
---|
333 | dgaus->SetParLimits(2,0,width/3.);
|
---|
334 | dgaus->SetParLimits(3,0,integ);
|
---|
335 | dgaus->SetParLimits(4,min,max);
|
---|
336 | dgaus->SetParLimits(5,0,width/3.);
|
---|
337 | obj2->Fit("dgaus","QLR");
|
---|
338 | obj2->GetFunction("dgaus")->SetLineColor(kYellow);
|
---|
339 | break;
|
---|
340 |
|
---|
341 | case 2:
|
---|
342 | TF1 *tgaus = new TF1("tgaus","gaus(0)+gaus(3)+gaus(6)",min,max);
|
---|
343 | tgaus->SetParNames("A1","#mu1","#sigma1","A2","#mu2","#sigma2","A3","#mu3","#sigma3");
|
---|
344 | tgaus->SetParameters(integ/width,max-width/6.,width/4.,
|
---|
345 | integ/width,min+width/6.,width/4.,
|
---|
346 | integ/width,min+width/6.,width/2.);
|
---|
347 | tgaus->SetParLimits(0,0,integ);
|
---|
348 | tgaus->SetParLimits(1,min,max);
|
---|
349 | tgaus->SetParLimits(2,0,width/4.);
|
---|
350 | tgaus->SetParLimits(3,0,integ);
|
---|
351 | tgaus->SetParLimits(4,min,max);
|
---|
352 | tgaus->SetParLimits(5,0,width/4.);
|
---|
353 | tgaus->SetParLimits(6,0,integ);
|
---|
354 | tgaus->SetParLimits(7,min,max);
|
---|
355 | tgaus->SetParLimits(8,0,width/2.);
|
---|
356 | obj2->Fit("tgaus","QLR");
|
---|
357 | obj2->GetFunction("tgaus")->SetLineColor(kYellow);
|
---|
358 | break;
|
---|
359 | case 3:
|
---|
360 | obj2->Fit("pol0","Q");
|
---|
361 | obj2->GetFunction("pol0")->SetLineColor(kYellow);
|
---|
362 | break;
|
---|
363 | case 9:
|
---|
364 | break;
|
---|
365 | default:
|
---|
366 | obj2->Fit("gaus","Q");
|
---|
367 | obj2->GetFunction("gaus")->SetLineColor(kYellow);
|
---|
368 | break;
|
---|
369 | }
|
---|
370 |
|
---|
371 | gStyle->SetStatH(he);
|
---|
372 | gStyle->SetStatW(wi);
|
---|
373 |
|
---|
374 | }
|
---|