source: trunk/MagicSoft/Mars/macros/bootcampstandardanalysis.C@ 3784

Last change on this file since 3784 was 3754, checked in by gaug, 21 years ago
*** empty log message ***
File size: 11.5 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): Javier López, 12/2003 <mailto:jlopez@ifae.es>
19! Markus Gaug , 04/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25/////////////////////////////////////////////////////////////////////////////
26//
27// bootcampstandardanalysis.C
28//
29// Updated version of the macro designed at the Wuerzburg bootcamp and
30// compatible with the latest changes in Mars for general usage at the
31// Udine bootcamp.
32//
33// Needs as arguments the run number of a pedestal file ("*_P_*.root"),
34// one of a calibration file ("*_C_*.root") and one of a data file
35// ("*_D_*.root"). Performs the pedestal calculation, the calibration
36/// constants calculation and the calibration of the data.
37//
38// The TString inpath has to be set correctly.
39//
40// The macro searches for the pulser colour which corresponds to the calibration
41// run number. If the run number is smaller than 20000, pulser colour "CT1"
42// is assumed, otherwise, it searches for the strings "green", "blue", "uv" or
43// "ct1" in the filenames. If no colour or multiple colours are found, the
44// execution is aborted.
45//
46//////////////////////////////////////////////////////////////////////////////////
47const TString inpath = "/mnt/Data/rootdata/CrabNebula/2004_02_10/";
48const Int_t dpedrun = 14607;
49const Int_t dcalrun1 = 14608;
50const Int_t dcalrun2 = 0;
51const Int_t ddatrun1 = 14609;
52const Int_t ddatrun2 = 14614;
53const Bool_t usedisplay = kTRUE;
54
55void bootcampstandardanalysis(const Int_t prun=dpedrun, // pedestal file
56 const Int_t crun1=dcalrun1, const Int_t crun2=dcalrun2, // calibration file(s)
57 const Int_t drun1=ddatrun1, const Int_t drun2=ddatrun2 // data files(s)
58 )
59
60{
61
62 MRunIter pruns;
63 MRunIter cruns;
64 MRunIter druns;
65
66 pruns.AddRun(prun,inpath);
67
68 if (crun2==0)
69 cruns.AddRun(crun1,inpath);
70 else
71 cruns.AddRuns(crun1,crun2,inpath);
72
73 MCalibrationCam::PulserColor_t color;
74
75 if (crun1 < 20000)
76 color = MCalibrationCam::kCT1;
77 else
78 color = FindColor((MDirIter*)&cruns);
79
80 if (color == MCalibrationCam::kNONE)
81 {
82 TTimer timer("gSystem->ProcessEvents();", 50, kFALSE);
83
84 while (1)
85 {
86 timer.TurnOn();
87 TString input = Getline("Could not find the correct colour: Type 'q' to exit, "
88 "green, blue, uv or ct1 to go on: ");
89 timer.TurnOff();
90
91 if (input=="q\n")
92 return ;
93
94 if (input=="green")
95 color = MCalibrationCam::kGREEN;
96 if (input=="blue")
97 color = MCalibrationCam::kBLUE;
98 if (input=="uv")
99 color = MCalibrationCam::kUV;
100 if (input=="ct1")
101 color = MCalibrationCam::kCT1;
102 }
103 }
104
105 if (drun2==0)
106 druns.AddRun(drun1,inpath);
107 else
108 druns.AddRuns(drun1,drun2,inpath);
109
110 //
111 // Now setup the tasks and tasklist for the pedestals:
112 // ---------------------------------------------------
113 //
114 MBadPixelsCam badcam;
115 MGeomCamMagic geomcam;
116 MGeomApply geomapl;
117 //
118 // If you want to exclude pixels from the beginning, read
119 // an ascii-file with the corr. pixel numbers (see MBadPixelsCam)
120 //
121 // badcam.AsciiRead("badpixels.dat");
122
123 /************************************/
124 /* FIRST LOOP: PEDESTAL COMPUTATION */
125 /************************************/
126
127 MJPedestal pedloop;
128 pedloop.SetInput(&pruns);
129 if (usedisplay)
130 {
131 MStatusDisplay *display = new MStatusDisplay;
132 display->SetUpdateTime(3000);
133 display->Resize(850,700);
134 pedloop.SetDisplay(display);
135 }
136 pedloop.SetBadPixels(badcam);
137
138 if (!pedloop.Process())
139 return;
140
141 /****************************************/
142 /* SECOND LOOP: CALIBRATION COMPUTATION */
143 /****************************************/
144
145 //
146 // Now setup the new tasks for the calibration:
147 // ---------------------------------------------------
148 //
149 MCalibrationQECam qecam;
150 MJCalibration calloop;
151 calloop.SetColor(color);
152 calloop.SetInput(&cruns);
153 //
154 // Use as signal extractor MExtractSignal:
155 //
156 calloop.SetExtractorLevel(1);
157 //
158 // Set the corr. cams:
159 //
160 calloop.SetQECam(qecam);
161 calloop.SetBadPixels(pedloop.GetBadPixels());
162 //
163 // The next two commands are for the display:
164 //
165 if (usedisplay)
166 {
167 calloop.SetDisplay(display);
168 calloop.SetDataCheck();
169 }
170
171 //
172 // Apply rel. time calibration:
173 //
174 calloop.SetRelTimeCalibration();
175 //
176 // Use as arrival time extractor MArrivalTimeCalc2:
177 //
178 calloop.SetArrivalTimeLevel(2);
179
180 //
181 // Do the event-loop:
182 //
183 if (!calloop.Process(pedloop.GetPedestalCam()))
184 return;
185
186
187 /*************************************/
188 /* THIRD LOOP: PEDESTAL CALIBRATION */
189 /*************************************/
190
191 //
192 // Create a empty Parameter List and an empty Task List
193 //
194 MParList plist3;
195 MTaskList tlist3;
196 plist3.AddToList(&tlist3);
197
198 //
199 // Now setup the tasks and tasklist to calculate the pedestal rms in number of photons
200 // -----------------------------------------------------------------------------------
201 //
202
203 MCerPhotEvt nphot;
204 MPedPhotCam nphotrms;
205
206 plist3.AddToList(&geomcam);
207
208 //
209 // Retrieve the cameras from the previous runs:
210 //
211 plist3.AddToList(&pedloop.GetPedestalCam());
212 plist3.AddToList(&calloop.GetCalibrationCam());
213 plist3.AddToList(&calloop.GetQECam());
214 plist3.AddToList(&calloop.GetRelTimeCam());
215 plist3.AddToList(&calloop.GetBadPixels());
216 plist3.AddToList(&nphot);
217 plist3.AddToList(&nphotrms);
218
219 //tasks
220 MReadMarsFile read3("Events");
221 read3.DisableAutoScheme();
222 static_cast<MRead&>(read3).AddFiles(pruns);
223
224 MExtractSignal extsig;
225 MCalibrate photcalc;
226 photcalc.SetCalibrationMode(MCalibrate::kFfactor);
227 // MPedPhotCalc photrmscalc; //It doesn't exist yet
228
229 tlist3.AddToList(&read3);
230 tlist3.AddToList(&geomapl);
231 tlist3.AddToList(&extsig);
232 tlist3.AddToList(&photcalc);
233 // tlist3.AddToList(&photrmscalc);
234
235 //
236 // Create and setup the eventloop
237 //
238 MEvtLoop evtloop3;
239 evtloop3.SetParList(&plist3);
240
241 //
242 // Execute first analysis
243 //
244 if (!evtloop3.Eventloop())
245 return;
246
247 tlist3.PrintStatistics();
248
249 /*************************************/
250 /* FOURTH LOOP: DATA CALIBRATION */
251 /*************************************/
252
253 //
254 // Create a empty Parameter List and an empty Task List
255 //
256 MParList plist4;
257 MTaskList tlist4;
258 plist4.AddToList(&tlist4);
259
260 //
261 // Now setup the tasks and tasklist to analize the data
262 // -----------------------------------------------------
263 //
264
265 plist4.AddToList(&geomcam);
266 //
267 // Retrieve the cameras from the previous runs:
268 //
269 plist4.AddToList(&pedloop.GetPedestalCam());
270 plist4.AddToList(&calloop.GetCalibrationCam());
271 plist4.AddToList(&calloop.GetQECam());
272 plist4.AddToList(&calloop.GetRelTimeCam());
273 plist4.AddToList(&calloop.GetBadPixels());
274 plist4.AddToList(&nphot);
275 plist4.AddToList(&nphotrms);
276
277 MArrivalTime times;
278 plist4.AddToList(&times);
279
280 //tasks
281 MReadMarsFile read4("Events");
282 read4.DisableAutoScheme();
283 static_cast<MRead&>(read4).AddFiles(druns);
284
285 MArrivalTimeCalc2 timecalc;
286 MCalibrateRelTimes timecal;
287
288 tlist4.AddToList(&read4);
289 tlist4.AddToList(&geomapl);
290 tlist4.AddToList(&extsig);
291 tlist4.AddToList(&timecalc);
292 tlist4.AddToList(&photcalc);
293 tlist4.AddToList(&timecal);
294
295 //
296 // Create and setup the eventloop
297 //
298 MEvtLoop evtloop4;
299 evtloop4.SetParList(&plist4);
300
301 if (!evtloop4.PreProcess())
302
303 return;
304
305 TCanvas *c1 = new TCanvas;
306 MHCamera disp1(geomcam);
307 disp1.SetPrettyPalette();
308 //disp1.SetInvDeepBlueSeaPalette()
309 disp1.Draw();
310 gPad->SetLogy();
311 gPad->cd(1);
312
313 /*
314 TCanvas *c2 = new TCanvas;
315 MHCamera disp2(geomcam);
316 disp2.SetPrettyPalette();
317 //disp2.SetInvDeepBlueSeaPalette()
318 disp2.Draw();
319 gPad->SetLogy();
320 gPad->cd(1);
321 */
322 while (tlist4.Process())
323 {
324 disp1.SetCamContent(nphot);
325
326 gPad->Modified();
327 gPad->Update();
328
329 /*
330 disp2.SetCamContent(times);
331
332 gPad->Modified();
333 gPad->Update();
334 */
335
336 // Remove the comments if you want to go through the file
337 // event-by-event:
338 if (!HandleInput())
339 break;
340 }
341
342 evtloop4.PostProcess();
343
344}
345
346Bool_t HandleInput()
347{
348 TTimer timer("gSystem->ProcessEvents();", 50, kFALSE);
349
350 while (1)
351 {
352 //
353 // While reading the input process gui events asynchronously
354 //
355 timer.TurnOn();
356 TString input = Getline("Type 'q' to exit, <return> to go on: ");
357 timer.TurnOff();
358
359 if (input=="q\n")
360 return kFALSE;
361
362 if (input=="\n")
363 return kTRUE;
364 };
365
366 return kFALSE;
367}
368
369MCalibrationCam::PulserColor_t FindColor(MDirIter* run)
370{
371
372 MCalibrationCam::PulserColor_t col = MCalibrationCam::kNONE;
373
374 TString filenames;
375
376 while (!(filenames=run->Next()).IsNull())
377 {
378
379 filenames.ToLower();
380
381 if (filenames.Contains("green"))
382 if (col == MCalibrationCam::kNONE)
383 {
384 cout << "Found colour: Green in " << filenames << endl;
385 col = MCalibrationCam::kGREEN;
386 }
387 else if (col != MCalibrationCam::kGREEN)
388 {
389 cout << "Different colour found in " << filenames << "... abort" << endl;
390 return MCalibrationCam::kNONE;
391 }
392
393 if (filenames.Contains("blue"))
394 if (col == MCalibrationCam::kNONE)
395 {
396 cout << "Found colour: Blue in " << filenames << endl;
397 col = MCalibrationCam::kBLUE;
398 }
399 else if (col != MCalibrationCam::kBLUE)
400 {
401 cout << "Different colour found in " << filenames << "... abort" << endl;
402 return MCalibrationCam::kNONE;
403 }
404
405 if (filenames.Contains("uv"))
406 if (col == MCalibrationCam::kNONE)
407 {
408 cout << "Found colour: Uv in " << filenames << endl;
409 col = MCalibrationCam::kUV;
410 }
411 else if (col != MCalibrationCam::kUV)
412 {
413 cout << "Different colour found in " << filenames << "... abort" << endl;
414 return MCalibrationCam::kNONE;
415 }
416
417 if (filenames.Contains("ct1"))
418 if (col == MCalibrationCam::kNONE)
419 {
420 cout << "Found colour: Ct1 in " << filenames << endl;
421 col = MCalibrationCam::kCT1;
422 }
423 else if (col != MCalibrationCam::kCT1)
424 {
425 cout << "Different colour found in " << filenames << "... abort" << endl;
426 return MCalibrationCam::kNONE;
427 }
428
429 }
430
431
432
433 if (col == MCalibrationCam::kNONE)
434 cout << "No colour found in filenames of runs: " << ((MRunIter*)run)->GetRunsAsString()
435 << "... abort" << endl;
436
437 return col;
438}
Note: See TracBrowser for help on using the repository browser.