source: trunk/MagicSoft/Mars/mjobs/MJCalibration.cc@ 3916

Last change on this file since 3916 was 3913, checked in by gaug, 21 years ago
*** empty log message ***
File size: 32.6 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): Thomas Bretz, 1/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Markus Gaug, 02/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25/////////////////////////////////////////////////////////////////////////////
26//
27// MJCalibration
28//
29// Do one calibration loop over serious of runs with the same pulser
30// colour and the same intensity. The following containers (rectangular) and
31// tasks (ellipses) are called to produce an MCalibrationChargeCam and to
32// update the MCalibrationQECam: (MCalibrate is not called from this class)
33//
34//Begin_Html
35/*
36<img src="images/CalibClasses.gif">
37*/
38//End_Html
39//
40// Different signal extractors can be set with the command SetExtractor()
41// Only extractors deriving from MExtractor can be set, default is MExtractSlidingWindow
42//
43// Different arrival time extractors can be set with the command SetTimeExtractor()
44// Only extractors deriving from MTimeExtractor can be set, default is MExtractTimeSpline
45//
46// At the end of the eventloop, plots and results are displayed, depending on
47// the flags set (see DisplayResult())
48//
49// If the flag SetFullDisplay() is set, all MHCameras will be displayed.
50// if the flag SetDataCheckDisplay() is set, only the most important ones are displayed
51// Otherwise, (default: SetNormalDisplay()), a good selection of plots is given
52//
53// If the flag SetDataCheck() is set, the calibration is used as in the data check at
54// La Palma, which mean especially running it on raw data files.
55//
56// See also: MHCalibrationChargePix, MHCalibrationChargeCam, MHGausEvents
57// MHCalibrationChargeBlindPix, MHCalibrationChargePINDiode
58// MCalibrationChargePix, MCalibrationChargeCam, MCalibrationChargeCalc
59// MCalibrationChargeBlindPix, MCalibrationChargePINDiode,
60// MCalibrationQECam, MBadPixelsPix, MBadPixelsCam
61//
62// If the flag RelTimeCalibration() is set, a calibration of the relative arrival
63// times is also performed. The following containers (rectangular) and
64// tasks (ellipses) are called to produce an MCalibrationRelTimeCam used by
65// MCalibrateTime to correct timing offset between pixels: (MCalibrateTime is not
66// called from this class)
67//
68//Begin_Html
69/*
70<img src="images/RelTimeClasses.gif">
71*/
72//End_Html
73//
74// Different arrival time extractors can be chosen via the command SetArrivalTimeLevel(UInt_t i)
75// Up to now, the following extractors are available:
76// i=1: Use MArrivalTimeCalc (using MCubicSpline, arrival time == position at half maximum)
77// i=2: Use MArrivalTimeCalc2 (mean time of fWindowSize time slices with the highest integral content: default)
78//
79// See also: MHCalibrationRelTimePix, MHCalibrationRelTimeCam, MHGausEvents
80// MCalibrationRelTimePix, MCalibrationRelTimeCam
81// MBadPixelsPix, MBadPixelsCam
82//
83/////////////////////////////////////////////////////////////////////////////
84#include "MJCalibration.h"
85
86#include <TFile.h>
87#include <TStyle.h>
88#include <TCanvas.h>
89#include <TSystem.h>
90
91#include "MLog.h"
92#include "MLogManip.h"
93
94#include "MRunIter.h"
95#include "MParList.h"
96#include "MTaskList.h"
97#include "MEvtLoop.h"
98
99#include "MHCamera.h"
100#include "MGeomCam.h"
101
102#include "MPedestalCam.h"
103#include "MCalibrationCam.h"
104#include "MCalibrationQECam.h"
105#include "MCalibrationChargeCam.h"
106#include "MCalibrationChargePINDiode.h"
107#include "MCalibrationChargeBlindPix.h"
108#include "MCalibrationChargeCalc.h"
109
110#include "MHGausEvents.h"
111#include "MHCalibrationCam.h"
112#include "MHCalibrationChargeCam.h"
113#include "MHCalibrationRelTimeCam.h"
114#include "MCalibrationRelTimeCam.h"
115
116#include "MReadMarsFile.h"
117#include "MRawFileRead.h"
118#include "MGeomApply.h"
119#include "MBadPixelsMerge.h"
120#include "MBadPixelsCam.h"
121#include "MTimeExtractor.h"
122#include "MExtractor.h"
123#include "MExtractPINDiode.h"
124#include "MExtractBlindPixel.h"
125#include "MExtractSlidingWindow.h"
126#include "MExtractTimeSpline.h"
127#include "MFCosmics.h"
128#include "MContinue.h"
129#include "MFillH.h"
130
131#include "MArrivalTimeCam.h"
132#include "MArrivalTimeCalc.h"
133#include "MArrivalTimeCalc2.h"
134
135#include "MStatusDisplay.h"
136
137ClassImp(MJCalibration);
138
139using namespace std;
140
141const Int_t MJCalibration::gkIFAEBoxInaugurationRun = 20113;
142// --------------------------------------------------------------------------
143//
144// Default constructor.
145//
146// Sets fRuns to 0, fExtractor to NULL, fTimeExtractor to NULL, fColor to kNONE,
147// fDisplay to kNormalDisplay, fRelTime to kFALSE, fDataCheck to kFALSE
148//
149MJCalibration::MJCalibration(const char *name, const char *title)
150 : fRuns(0), fExtractor(NULL), fTimeExtractor(NULL),
151 fColor(MCalibrationCam::kNONE), fDisplayType(kNormalDisplay),
152 fRelTimes(kFALSE), fDataCheck(kFALSE)
153{
154 fName = name ? name : "MJCalibration";
155 fTitle = title ? title : "Tool to create the calibration constants for one calibration run";
156
157}
158
159
160// --------------------------------------------------------------------------
161//
162// Display the results in MStatusDisplay:
163//
164// - Add "Calibration" to the MStatusDisplay title
165// - Retrieve the MGeomCam from MParList
166// - Initialize the following MHCamera's:
167// 1) MCalibrationPix::GetMean()
168// 2) MCalibrationPix::Sigma()
169// 3) MCalibrationChargePix::GetRSigma()
170// 4) MCalibrationChargePix::GetRSigmaPerCharge()
171// 5) MCalibrationChargePix::GetPheFFactorMethod()
172// 6) MCalibrationChargePix::GetMeanConvFADC2Phe()
173// 7) MCalibrationChargePix::GetMeanFFactorFADC2Phot()
174// 8) MCalibrationQEPix::GetQECascadesFFactor()
175// 9) MCalibrationQEPix::GetQECascadesBlindPixel()
176// 10) MCalibrationQEPix::GetQECascadesPINDiode()
177// 11) MCalibrationQEPix::GetQECascadesCombined()
178// 12) MCalibrationQEPix::IsAverageQEFFactorAvailable()
179// 13) MCalibrationQEPix::IsAverageQEBlindPixelAvailable()
180// 14) MCalibrationQEPix::IsAverageQEPINDiodeAvailable()
181// 15) MCalibrationQEPix::IsAverageQECombinedAvailable()
182// 16) MCalibrationChargePix::IsHiGainSaturation()
183// 17) MCalibrationPix::GetHiLoMeansDivided()
184// 18) MCalibrationPix::GetHiLoSigmasDivided()
185// 19) MCalibrationChargePix::GetHiGainPickup()
186// 20) MCalibrationChargePix::GetLoGainPickup()
187// 21) MCalibrationChargePix::GetHiGainBlackout()
188// 22) MCalibrationChargePix::GetLoGainBlackout()
189// 23) MCalibrationPix::IsExcluded()
190// 24) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnsuitableRun)
191// 25) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnreliableRun)
192// 26) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kHiGainOscillating)
193// 27) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kLoGainOscillating)
194// 28) MCalibrationChargePix::GetAbsTimeMean()
195// 29) MCalibrationChargePix::GetAbsTimeRms()
196//
197// If the flag SetFullDisplay() is set, all MHCameras will be displayed.
198// if the flag SetDataCheckDisplay() is set, only the most important ones are displayed
199// and otherwise, (default: SetNormalDisplay()), a good selection of plots is given
200//
201void MJCalibration::DisplayResult(MParList &plist)
202{
203 if (!fDisplay)
204 return;
205
206 //
207 // Update display
208 //
209 TString title = fDisplay->GetTitle();
210 title += "-- Calibration ";
211 title += fRuns->GetRunsAsString();
212 title += " --";
213 fDisplay->SetTitle(title);
214
215 //
216 // Get container from list
217 //
218 MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
219
220 // Create histograms to display
221 MHCamera disp1 (geomcam, "Charge", "Fitted Mean Charges");
222 MHCamera disp2 (geomcam, "SigmaCharge", "Sigma of Fitted Charges");
223 MHCamera disp3 (geomcam, "RSigma", "Reduced Sigmas");
224 MHCamera disp4 (geomcam, "RSigmaPerCharge", "Reduced Sigma per Charge");
225 MHCamera disp5 (geomcam, "NumPhes", "Nr. of Phe's (F-Factor Method)");
226 MHCamera disp6 (geomcam, "ConvFADC2Phes", "Conversion Factor (F-Factor Method)");
227 MHCamera disp7 (geomcam, "TotalFFactor", "Total F-Factor (F-Factor Method)");
228 MHCamera disp8 (geomcam, "CascadesQEFFactor", "Cascades QE (F-Factor Method)");
229 MHCamera disp9 (geomcam, "CascadesQEBlindPix","Cascades QE (Blind Pixel Method)");
230 MHCamera disp10(geomcam, "CascadesQEPINDiode","Cascades QE (PIN Diode Method)");
231 MHCamera disp11(geomcam, "CascadesQECombined","Cascades QE (Combined Method)");
232 MHCamera disp12(geomcam, "FFactorValid", "Pixels with valid F-Factor calibration");
233 MHCamera disp13(geomcam, "BlindPixelValid", "Pixels with valid BlindPixel calibration");
234 MHCamera disp14(geomcam, "PINdiodeValid", "Pixels with valid PINDiode calibration");
235 MHCamera disp15(geomcam, "CombinedValid", "Pixels with valid Combined calibration");
236 MHCamera disp16(geomcam, "Saturation", "Pixels with saturated Hi Gain");
237 MHCamera disp17(geomcam, "ConversionMeans", "Conversion HiGain.vs.LoGain Means");
238 MHCamera disp18(geomcam, "ConversionSigmas", "Conversion HiGain.vs.LoGain Sigmas");
239 MHCamera disp19(geomcam, "HiGainPickup", "Number Pickup events Hi Gain");
240 MHCamera disp20(geomcam, "LoGainPickup", "Number Pickup events Lo Gain");
241 MHCamera disp21(geomcam, "HiGainBlackout", "Number Blackout events Hi Gain");
242 MHCamera disp22(geomcam, "LoGainBlackout", "Number Blackout events Lo Gain");
243 MHCamera disp23(geomcam, "Excluded", "Pixels previously excluded");
244 MHCamera disp24(geomcam, "Bad;UnSuitable", "Pixels not suited for further analysis");
245 MHCamera disp25(geomcam, "Bad;UnReliable", "Pixels not reliable for further analysis");
246 MHCamera disp26(geomcam, "Bad;HiGainOscillating", "Oscillating Pixels High Gain");
247 MHCamera disp27(geomcam, "Bad;LoGainOscillating", "Oscillating Pixels Low Gain");
248 MHCamera disp28(geomcam, "Cal;AbsTimeMean", "Abs. Arrival Times");
249 MHCamera disp29(geomcam, "Cal;AbsTimeRms", "RMS of Arrival Times");
250 MHCamera disp30(geomcam, "time;MeanTime", "Mean Rel. Arrival Times");
251 MHCamera disp31(geomcam, "time;SigmaTime", "Sigma Rel. Arrival Times");
252 MHCamera disp32(geomcam, "time;TimeProb", "Probability of Time Fit");
253 MHCamera disp33(geomcam, "time;NotFitValid", "Pixels with not valid fit results");
254 MHCamera disp34(geomcam, "time;Oscillating", "Oscillating Pixels");
255
256
257
258 // Fitted charge means and sigmas
259 disp1.SetCamContent(fCalibrationCam, 0);
260 disp1.SetCamError( fCalibrationCam, 1);
261 disp2.SetCamContent(fCalibrationCam, 2);
262 disp2.SetCamError( fCalibrationCam, 3);
263
264 // Reduced Sigmas and reduced sigmas per charge
265 disp3.SetCamContent(fCalibrationCam, 5);
266 disp3.SetCamError( fCalibrationCam, 6);
267 disp4.SetCamContent(fCalibrationCam, 7);
268 disp4.SetCamError( fCalibrationCam, 8);
269
270 // F-Factor Method
271 disp5.SetCamContent(fCalibrationCam, 9);
272 disp5.SetCamError( fCalibrationCam, 10);
273 disp6.SetCamContent(fCalibrationCam, 11);
274 disp6.SetCamError( fCalibrationCam, 12);
275 disp7.SetCamContent(fCalibrationCam, 13);
276 disp7.SetCamError( fCalibrationCam, 14);
277
278 // Quantum Efficiencies
279 disp8.SetCamContent (fQECam, 0 );
280 disp8.SetCamError (fQECam, 1 );
281 disp9.SetCamContent (fQECam, 2 );
282 disp9.SetCamError (fQECam, 3 );
283 disp10.SetCamContent(fQECam, 4 );
284 disp10.SetCamError (fQECam, 5 );
285 disp11.SetCamContent(fQECam, 6 );
286 disp11.SetCamError (fQECam, 7 );
287
288 // Valid flags
289 disp12.SetCamContent(fQECam, 8 );
290 disp13.SetCamContent(fQECam, 9 );
291 disp14.SetCamContent(fQECam, 10);
292 disp15.SetCamContent(fQECam, 11);
293
294 // Conversion Hi-Lo
295 disp16.SetCamContent(fCalibrationCam, 25);
296 disp17.SetCamContent(fCalibrationCam, 16);
297 disp17.SetCamError (fCalibrationCam, 17);
298 disp18.SetCamContent(fCalibrationCam, 18);
299 disp18.SetCamError (fCalibrationCam, 19);
300
301 // Pickup and Blackout
302 disp19.SetCamContent(fCalibrationCam, 21);
303 disp20.SetCamContent(fCalibrationCam, 22);
304 disp21.SetCamContent(fCalibrationCam, 23);
305 disp22.SetCamContent(fCalibrationCam, 24);
306
307 // Pixels with defects
308 disp23.SetCamContent(fCalibrationCam, 20);
309 disp24.SetCamContent(fBadPixels, 1);
310 disp25.SetCamContent(fBadPixels, 3);
311
312 // Oscillations
313 disp26.SetCamContent(fBadPixels, 10);
314 disp27.SetCamContent(fBadPixels, 11);
315
316 // Arrival Times
317 disp28.SetCamContent(fCalibrationCam, 26);
318 disp28.SetCamError( fCalibrationCam, 27);
319 disp29.SetCamContent(fCalibrationCam, 27);
320
321 disp1.SetYTitle("Q [FADC counts]");
322 disp2.SetYTitle("\\sigma_{Q} [FADC counts]");
323
324 disp3.SetYTitle("\\sqrt{\\sigma^{2}_{Q} - RMS^{2}_{Ped}} [FADC Counts]");
325 disp4.SetYTitle("Red.Sigma/<Q> [1]");
326
327 disp5.SetYTitle("Nr. Phe's [1]");
328 disp6.SetYTitle("Conv.Factor [PhE/FADC counts]");
329 disp7.SetYTitle("Total F-Factor [1]");
330
331 disp8.SetYTitle("QE [1]");
332 disp9.SetYTitle("QE [1]");
333 disp10.SetYTitle("QE [1]");
334 disp11.SetYTitle("QE [1]");
335
336 disp12.SetYTitle("[1]");
337 disp13.SetYTitle("[1]");
338 disp14.SetYTitle("[1]");
339 disp15.SetYTitle("[1]");
340 disp16.SetYTitle("[1]");
341
342 disp17.SetYTitle("<Q>(High)/<Q>(Low) [1]");
343 disp18.SetYTitle("\\sigma_{Q}(High)/\\sigma_{Q}(Low) [1]");
344
345 disp19.SetYTitle("[1]");
346 disp20.SetYTitle("[1]");
347 disp21.SetYTitle("[1]");
348 disp22.SetYTitle("[1]");
349 disp23.SetYTitle("[1]");
350 disp24.SetYTitle("[1]");
351 disp25.SetYTitle("[1]");
352 disp26.SetYTitle("[1]");
353 disp27.SetYTitle("[1]");
354
355 disp28.SetYTitle("Mean Abs. Time [FADC slice]");
356 disp29.SetYTitle("RMS Abs. Time [FADC slices]");
357
358 if (fRelTimes)
359 {
360
361 disp30.SetCamContent(fRelTimeCam,0);
362 disp30.SetCamError( fRelTimeCam,1);
363 disp31.SetCamContent(fRelTimeCam,2);
364 disp31.SetCamError( fRelTimeCam,3);
365 disp32.SetCamContent(fRelTimeCam,4);
366 disp33.SetCamContent(fBadPixels,20);
367 disp34.SetCamContent(fBadPixels,21);
368
369 disp30.SetYTitle("Time Offset [FADC units]");
370 disp31.SetYTitle("Timing resolution [FADC units]");
371 disp32.SetYTitle("P_{Time} [1]");
372 disp33.SetYTitle("[1]");
373 disp34.SetYTitle("[1]");
374 }
375
376 if (fDisplayType == kDataCheckDisplay)
377 {
378 TCanvas &c1 = fDisplay->AddTab("Fit.Charge");
379 c1.Divide(3, 3);
380
381 CamDraw(c1, 1, 3, disp1, 2);
382 CamDraw(c1, 2, 3, disp4, 2);
383 CamDraw(c1, 3, 3, disp28, 2);
384
385 // F-Factor
386 TCanvas &c2 = fDisplay->AddTab("Phe's");
387 c2.Divide(3,4);
388
389 CamDraw(c2, 1, 3, disp6, 2, 1);
390 CamDraw(c2, 2, 3, disp7, 2, 1);
391 CamDraw(c2, 3, 3, disp8, 2, 1);
392
393 // QE's
394 TCanvas &c3 = fDisplay->AddTab("QE's");
395 c3.Divide(3,4);
396
397 CamDraw(c3, 1, 3, disp8, 2, 1);
398 CamDraw(c3, 2, 3, disp9, 2, 1);
399 CamDraw(c3, 3, 3, disp10, 2, 1);
400
401 // Defects
402 TCanvas &c4 = fDisplay->AddTab("Defect");
403 // c4.Divide(3,2);
404 c4.Divide(2,2);
405
406 // CamDraw(c4, 1, 3, disp23, 0);
407 // CamDraw(c4, 2, 3, disp24, 0);
408 // CamDraw(c4, 3, 3, disp25, 0);
409 CamDraw(c4, 1, 2, disp24, 0);
410 CamDraw(c4, 2, 2, disp25, 0);
411
412 if (fRelTimes)
413 {
414 // Rel. Times
415 TCanvas &c5 = fDisplay->AddTab("Rel. Times");
416 c5.Divide(2,4);
417
418 CamDraw(c5, 1, 2, disp30, 2);
419 CamDraw(c5, 2, 2, disp31, 2);
420 }
421
422
423 return;
424 }
425
426 if (fDisplayType == kNormalDisplay)
427 {
428
429 // Charges
430 TCanvas &c11 = fDisplay->AddTab("Fit.Charge");
431 c11.Divide(2, 4);
432
433 CamDraw(c11, 1, 2, disp1, 5, 1);
434 CamDraw(c11, 2, 2, disp2, 5, 1);
435
436 // Reduced Sigmas
437 TCanvas &c12 = fDisplay->AddTab("Red.Sigma");
438 c12.Divide(2,4);
439
440 CamDraw(c12, 1, 2, disp3, 5, 1);
441 CamDraw(c12, 2, 2, disp4, 5, 1);
442
443 // F-Factor
444 TCanvas &c13 = fDisplay->AddTab("Phe's");
445 c13.Divide(3,4);
446
447 CamDraw(c13, 1, 3, disp5, 5, 1);
448 CamDraw(c13, 2, 3, disp6, 5, 1);
449 CamDraw(c13, 3, 3, disp7, 5, 1);
450
451 // QE's
452 TCanvas &c14 = fDisplay->AddTab("QE's");
453 c14.Divide(4,4);
454
455 CamDraw(c14, 1, 4, disp8, 5, 1);
456 CamDraw(c14, 2, 4, disp9, 5, 1);
457 CamDraw(c14, 3, 4, disp10, 5, 1);
458 CamDraw(c14, 4, 4, disp11, 5, 1);
459
460 // Defects
461 TCanvas &c15 = fDisplay->AddTab("Defect");
462 // c15.Divide(5,2);
463 c15.Divide(4,2);
464
465 /*
466 CamDraw(c15, 1, 5, disp23, 0);
467 CamDraw(c15, 2, 5, disp24, 0);
468 CamDraw(c15, 3, 5, disp25, 0);
469 CamDraw(c15, 4, 5, disp26, 0);
470 CamDraw(c15, 5, 5, disp27, 0);
471 */
472 CamDraw(c15, 1, 4, disp24, 0);
473 CamDraw(c15, 2, 4, disp25, 0);
474 CamDraw(c15, 3, 4, disp26, 0);
475 CamDraw(c15, 4, 4, disp27, 0);
476
477 // Abs. Times
478 TCanvas &c16 = fDisplay->AddTab("Abs. Times");
479 c16.Divide(2,3);
480
481 CamDraw(c16, 1, 2, disp28, 5);
482 CamDraw(c16, 2, 2, disp29, 5);
483
484 if (fRelTimes)
485 {
486 // Rel. Times
487 TCanvas &c17 = fDisplay->AddTab("Rel. Times");
488 c17.Divide(2,4);
489
490 CamDraw(c17, 1, 2, disp30, 5, 1);
491 CamDraw(c17, 2, 2, disp31, 5, 1);
492 }
493
494 return;
495 }
496
497 if (fDisplayType == kFullDisplay)
498 {
499
500 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationChargeCam");
501
502 for (Int_t aidx=0;aidx<cam->GetAverageAreas();aidx++)
503 {
504 cam->GetAverageHiGainArea(aidx).DrawClone("all");
505 cam->GetAverageLoGainArea(aidx).DrawClone("all");
506 }
507
508 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
509 {
510 cam->GetAverageHiGainSector(sector).DrawClone("all");
511 cam->GetAverageLoGainSector(sector).DrawClone("all");
512 }
513
514 // Charges
515 TCanvas &c21 = fDisplay->AddTab("Fit.Charge");
516 c21.Divide(2, 4);
517
518 CamDraw(c21, 1, 2, disp1, 2, 1);
519 CamDraw(c21, 2, 2, disp2, 2, 1);
520
521 // Reduced Sigmas
522 TCanvas &c23 = fDisplay->AddTab("Red.Sigma");
523 c23.Divide(2,4);
524
525 CamDraw(c23, 1, 2, disp3, 2, 1);
526 CamDraw(c23, 2, 2, disp4, 2, 1);
527
528 // F-Factor
529 TCanvas &c24 = fDisplay->AddTab("Phe's");
530 c24.Divide(3,4);
531
532 CamDraw(c24, 1, 3, disp5, 2, 1);
533 CamDraw(c24, 2, 3, disp6, 2, 1);
534 CamDraw(c24, 3, 3, disp7, 2, 1);
535
536 // QE's
537 TCanvas &c25 = fDisplay->AddTab("QE's");
538 c25.Divide(4,4);
539
540 CamDraw(c25, 1, 4, disp8, 2, 1);
541 CamDraw(c25, 2, 4, disp9, 2, 1);
542 CamDraw(c25, 3, 4, disp10, 2, 1);
543 CamDraw(c25, 4, 4, disp11, 2, 1);
544
545 // Validity
546 TCanvas &c26 = fDisplay->AddTab("Valid");
547 c26.Divide(4,2);
548
549 CamDraw(c26, 1, 4, disp12, 0);
550 CamDraw(c26, 2, 4, disp13, 0);
551 CamDraw(c26, 3, 4, disp14, 0);
552 CamDraw(c26, 4, 4, disp15, 0);
553
554 // Other info
555 TCanvas &c27 = fDisplay->AddTab("HiLoGain");
556 c27.Divide(3,3);
557
558 CamDraw(c27, 1, 3, disp16, 0);
559 CamDraw(c27, 2, 3, disp17, 1);
560 CamDraw(c27, 3, 3, disp18, 1);
561
562 // Pickup
563 TCanvas &c28 = fDisplay->AddTab("Pickup");
564 c28.Divide(4,2);
565
566 CamDraw(c28, 1, 4, disp19, 0);
567 CamDraw(c28, 2, 4, disp20, 0);
568 CamDraw(c28, 3, 4, disp21, 0);
569 CamDraw(c28, 4, 4, disp22, 0);
570
571 // Defects
572 TCanvas &c29 = fDisplay->AddTab("Defect");
573 // c29.Divide(5,2);
574 c29.Divide(4,2);
575
576 CamDraw(c29, 1, 4, disp24, 0);
577 CamDraw(c29, 2, 4, disp25, 0);
578 CamDraw(c29, 3, 4, disp26, 0);
579 CamDraw(c29, 4, 4, disp27, 0);
580
581 // Abs. Times
582 TCanvas &c30 = fDisplay->AddTab("Abs. Times");
583 c30.Divide(2,3);
584
585 CamDraw(c30, 1, 2, disp28, 2);
586 CamDraw(c30, 2, 2, disp29, 1);
587
588 if (fRelTimes)
589 {
590 // Rel. Times
591 TCanvas &c31 = fDisplay->AddTab("Rel. Times");
592 c31.Divide(3,4);
593
594 CamDraw(c31, 1, 3, disp30, 2, 1);
595 CamDraw(c31, 2, 3, disp31, 2, 1);
596 CamDraw(c31, 3, 3, disp32, 4, 1);
597
598 // Time Defects
599 TCanvas &c32 = fDisplay->AddTab("Time Def.");
600 c32.Divide(2,2);
601
602 CamDraw(c32, 1, 2, disp33,0);
603 CamDraw(c32, 2, 2, disp34,0);
604
605 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationRelTimeCam");
606
607 for (Int_t aidx=0;aidx<cam->GetAverageAreas();aidx++)
608 {
609 cam->GetAverageHiGainArea(aidx).DrawClone("fourierevents");
610 cam->GetAverageLoGainArea(aidx).DrawClone("fourierevents");
611 }
612
613 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
614 {
615 cam->GetAverageHiGainSector(sector).DrawClone("fourierevents");
616 cam->GetAverageLoGainSector(sector).DrawClone("fourierevents");
617 }
618
619 }
620
621 return;
622 }
623}
624
625
626
627// --------------------------------------------------------------------------
628//
629// Find the colour of the pulsing LED:
630// - If the run number is smaller than gkIFAEBoxInaugurationRun, take MCalibrationCam::kCT1
631// - Otherwise find the colour out of the run name
632// - If no colour is found, return kFALSE
633//
634Bool_t MJCalibration::FindColor()
635{
636
637 const UInt_t nruns = fRuns->GetNumRuns();
638
639 if (nruns == 0)
640 return kFALSE;
641
642 TArrayI arr = fRuns->GetRuns();
643
644 if (arr[nruns-1] < gkIFAEBoxInaugurationRun)
645 {
646 *fLog << "Found colour kCT1 in runs: " << fRuns->GetRunsAsString() << endl;
647 fColor = MCalibrationCam::kCT1;
648 return kTRUE;
649 }
650
651 TString filenames;
652
653 while (!(filenames=((MDirIter*)fRuns)->Next()).IsNull())
654 {
655
656 filenames.ToLower();
657
658 if (filenames.Contains("green"))
659 if (fColor == MCalibrationCam::kNONE)
660 {
661 *fLog << "Found colour: kGREEN in " << filenames << endl;
662 fColor = MCalibrationCam::kGREEN;
663 }
664 else if (fColor != MCalibrationCam::kGREEN)
665 {
666 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
667 return kFALSE;
668 }
669
670 if (filenames.Contains("blue"))
671 if (fColor == MCalibrationCam::kNONE)
672 {
673 *fLog << "Found colour: kBLUE in " << filenames << endl;
674 fColor = MCalibrationCam::kBLUE;
675 }
676 else if (fColor != MCalibrationCam::kBLUE)
677 {
678 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
679 return kFALSE;
680 }
681
682 if (filenames.Contains("uv"))
683 if (fColor == MCalibrationCam::kNONE)
684 {
685 *fLog << "Found colour: kUV in " << filenames << endl;
686 fColor = MCalibrationCam::kUV;
687 }
688 else if (fColor != MCalibrationCam::kUV)
689 {
690 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
691 return kFALSE;
692 }
693
694 if (filenames.Contains("ct1"))
695 if (fColor == MCalibrationCam::kNONE)
696 {
697 *fLog << "Found colour: kCT1 in " << filenames << endl;
698 fColor = MCalibrationCam::kCT1;
699 }
700 else if (fColor != MCalibrationCam::kCT1)
701 {
702 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
703 return kFALSE;
704 }
705
706 }
707
708 if (fColor == MCalibrationCam::kNONE)
709 {
710 *fLog << "No colour found in filenames of runs: " << fRuns->GetRunsAsString()
711 << "... abort" << endl;
712 return kFALSE;
713 }
714
715 return kTRUE;
716}
717
718
719
720
721// --------------------------------------------------------------------------
722//
723// Retrieve the output file written by WriteResult()
724//
725TString MJCalibration::GetOutputFile() const
726{
727 if (!fRuns)
728 return "";
729
730 return Form("%s/%s-F1.root", (const char*)fOutputPath, (const char*)fRuns->GetRunsAsFileName());
731}
732
733
734// --------------------------------------------------------------------------
735//
736// Call the ProcessFile(MPedestalCam)
737//
738Bool_t MJCalibration::Process(MPedestalCam &pedcam)
739{
740 if (!ReadCalibrationCam())
741 return ProcessFile(pedcam);
742
743 return kTRUE;
744}
745
746// --------------------------------------------------------------------------
747//
748// Execute the task list and the eventloop:
749//
750// - Check if there are fRuns, otherwise return
751// - Check the colour of the files in fRuns (FindColor()), otherwise return
752// - Check for consistency between run numbers and number of files
753// - Add fRuns to MReadMarsFile
754// - Put into MParList:
755// 1) MPedestalCam (pedcam)
756// 2) MCalibrationQECam (fQECam)
757// 3) MCalibrationChargeCam (fCalibrationCam)
758// 4) MCalibrationRelTimeCam (fRelTimeCam) (only if flag fRelTimes is chosen)
759// 5) MBadPixelsCam (fBadPixels)
760// 6) MCalibrationChargePINDiode
761// 7) MCalibrationChargeBlindPix
762// - Put into the MTaskList:
763// 1) MReadMarsFile
764// 2) MBadPixelsMerge
765// 3) MGeomApply
766// 4) MExtractor
767// 5) MExtractPINDiode
768// 6) MExtractBlindPixel
769// 7) MTimeExtractor (only if flag fRelTimes is chosen)
770// 8) MContinue(MFCosmics)
771// 9) MFillH("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode")
772// 10) MFillH("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel")
773// 11) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam")
774// 12) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam")
775// 13) MCalibrationChargeCalc
776// 14) MFillH("MHCalibrationRelTimeCam", "MArrivalTimeCam") (only if flag fRelTimes is chosen)
777// - Execute MEvtLoop
778// - DisplayResult()
779// - WriteResult()
780//
781Bool_t MJCalibration::ProcessFile(MPedestalCam &pedcam)
782{
783 if (!fRuns)
784 {
785 *fLog << err << "No Runs choosen... abort." << endl;
786 return kFALSE;
787 }
788
789 if (fRuns->GetNumRuns() != fRuns->GetNumEntries())
790 {
791 *fLog << err << "Number of files found doesn't match number of runs... abort."
792 << fRuns->GetNumRuns() << " vs. " << fRuns->GetNumEntries() << endl;
793 return kFALSE;
794 }
795
796 *fLog << inf;
797 fLog->Separator(GetDescriptor());
798
799 if (!FindColor())
800 return kFALSE;
801
802 *fLog << "Calculate MCalibrationCam from Runs " << fRuns->GetRunsAsString() << endl;
803 *fLog << endl;
804
805 // Setup Tasklist
806 MParList plist;
807 MTaskList tlist;
808 plist.AddToList(&tlist);
809
810 MReadMarsFile read("Events");
811 MRawFileRead rawread("");
812
813 if (fDataCheck)
814 {
815// rawread.AddFiles(*fRuns);
816 tlist.AddToList(&rawread);
817 }
818 else
819 {
820 read.DisableAutoScheme();
821 static_cast<MRead&>(read).AddFiles(*fRuns);
822 tlist.AddToList(&read);
823 }
824
825
826 MCalibrationChargePINDiode pindiode;
827 MCalibrationChargeBlindPix blindpix;
828
829 plist.AddToList(&pedcam);
830 plist.AddToList(&fBadPixels);
831 plist.AddToList(&fQECam);
832 plist.AddToList(&fCalibrationCam);
833 plist.AddToList(&pindiode);
834 plist.AddToList(&blindpix);
835
836 MGeomApply apply;
837 // MBadPixelsMerge merge(&fBadPixels);
838 MExtractPINDiode pinext;
839 MExtractBlindPixel blindext;
840 MExtractSlidingWindow extract2;
841 MExtractTimeSpline timespline;
842 MCalibrationChargeCalc calcalc;
843
844 //
845 // As long as there are no DM's, have to colour by hand
846 //
847 calcalc.SetPulserColor(fColor);
848
849 MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode");
850 MFillH fillbnd("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel");
851 MFillH fillcam("MHCalibrationChargeCam", "MExtractedSignalCam");
852 MFillH filltme("MHCalibrationRelTimeCam", "MArrivalTimeCam");
853 fillpin.SetNameTab("PINDiode");
854 fillbnd.SetNameTab("BlindPix");
855 fillcam.SetNameTab("Charge");
856 filltme.SetNameTab("RelTimes");
857
858
859 //
860 // Apply a filter against cosmics
861 // (will have to be needed in the future
862 // when the calibration hardware-trigger is working)
863 //
864 MFCosmics cosmics;
865 MContinue cont(&cosmics);
866
867 // tlist.AddToList(&merge);
868 tlist.AddToList(&apply);
869
870 if (fExtractor)
871 tlist.AddToList(fExtractor);
872 else
873 {
874 *fLog << warn << GetDescriptor()
875 << ": No extractor has been chosen, take default MExtractSlidingWindow " << endl;
876 tlist.AddToList(&extract2);
877 }
878
879
880 tlist.AddToList(&pinext);
881 tlist.AddToList(&blindext);
882
883 if (fRelTimes)
884 {
885 plist.AddToList(&fRelTimeCam);
886
887 if (fTimeExtractor)
888 tlist.AddToList(fTimeExtractor);
889 else
890 {
891 *fLog << warn << GetDescriptor()
892 << ": No extractor has been chosen, take default MTimeExtractSpine " << endl;
893 tlist.AddToList(&timespline);
894 }
895 }
896
897 tlist.AddToList(&cont);
898 tlist.AddToList(&fillcam);
899 tlist.AddToList(&fillpin);
900 tlist.AddToList(&fillbnd);
901 tlist.AddToList(&calcalc);
902
903 if (fRelTimes)
904 tlist.AddToList(&filltme);
905
906 // Create and setup the eventloop
907 MEvtLoop evtloop(fName);
908 evtloop.SetParList(&plist);
909 evtloop.SetDisplay(fDisplay);
910 evtloop.SetLogStream(fLog);
911
912 // Execute first analysis
913 if (!evtloop.Eventloop())
914 {
915 *fLog << err << GetDescriptor() << ": Failed." << endl;
916 return kFALSE;
917 }
918
919 tlist.PrintStatistics();
920
921 DisplayResult(plist);
922
923 if (!WriteResult())
924 return kFALSE;
925
926 *fLog << inf << GetDescriptor() << ": Done." << endl;
927
928 return kTRUE;
929}
930
931// --------------------------------------------------------------------------
932//
933// Read the following containers from GetOutputFile()
934// - MCalibrationChargeCam
935// - MCalibrationQECam
936// - MBadPixelsCam
937//
938Bool_t MJCalibration::ReadCalibrationCam()
939{
940 const TString fname = GetOutputFile();
941
942 if (gSystem->AccessPathName(fname, kFileExists))
943 {
944 *fLog << err << "Input file " << fname << " doesn't exist." << endl;
945 return kFALSE;
946 }
947
948 *fLog << inf << "Reading from file: " << fname << endl;
949
950 TFile file(fname, "READ");
951 if (fCalibrationCam.Read()<=0)
952 {
953 *fLog << err << "Unable to read MCalibrationChargeCam from " << fname << endl;
954 return kFALSE;
955 }
956
957 if (fQECam.Read()<=0)
958 {
959 *fLog << err << "Unable to read MCalibrationQECam from " << fname << endl;
960 return kFALSE;
961 }
962
963 if (file.FindKey("MBadPixelsCam"))
964 {
965 MBadPixelsCam bad;
966 if (bad.Read()<=0)
967 {
968 *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
969 return kFALSE;
970 }
971 fBadPixels.Merge(bad);
972 }
973
974 if (fDisplay /*&& !fDisplay->GetCanvas("Pedestals")*/) // FIXME!
975 fDisplay->Read();
976
977 return kTRUE;
978}
979
980
981// --------------------------------------------------------------------------
982//
983// Set the path for output files, written by WriteResult()
984//
985void MJCalibration::SetOutputPath(const char *path)
986{
987 fOutputPath = path;
988 if (fOutputPath.EndsWith("/"))
989 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
990}
991
992
993// --------------------------------------------------------------------------
994//
995// Write the result into the output file GetOutputFile(), if fOutputPath exists.
996//
997// The following containers are written:
998// - MStatusDisplay
999// - MCalibrationChargeCam
1000// - MCalibrationQECam
1001// - MBadPixelsCam
1002//
1003Bool_t MJCalibration::WriteResult()
1004{
1005 if (fOutputPath.IsNull())
1006 return kTRUE;
1007
1008 const TString oname(GetOutputFile());
1009
1010 *fLog << inf << "Writing to file: " << oname << endl;
1011
1012 TFile file(oname, "UPDATE");
1013
1014 if (fDisplay && fDisplay->Write()<=0)
1015 {
1016 *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
1017 return kFALSE;
1018 }
1019
1020 if (fCalibrationCam.Write()<=0)
1021 {
1022 *fLog << err << "Unable to write MCalibrationChargeCam to " << oname << endl;
1023 return kFALSE;
1024 }
1025
1026 if (fQECam.Write()<=0)
1027 {
1028 *fLog << err << "Unable to write MCalibrationQECam to " << oname << endl;
1029 return kFALSE;
1030 }
1031
1032 if (fBadPixels.Write()<=0)
1033 {
1034 *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl;
1035 return kFALSE;
1036 }
1037
1038 return kTRUE;
1039
1040}
1041
Note: See TracBrowser for help on using the repository browser.