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

Last change on this file since 3976 was 3975, checked in by gaug, 21 years ago
*** empty log message ***
File size: 35.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 MExtractTime 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 "MHCalibrationChargeBlindPix.h"
114#include "MHCalibrationRelTimeCam.h"
115#include "MCalibrationRelTimeCam.h"
116#include "MCalibrationRelTimeCalc.h"
117
118#include "MReadMarsFile.h"
119#include "MRawFileRead.h"
120#include "MGeomApply.h"
121#include "MBadPixelsMerge.h"
122#include "MBadPixelsCam.h"
123#include "MExtractTime.h"
124#include "MExtractor.h"
125#include "MExtractPINDiode.h"
126#include "MExtractBlindPixel.h"
127#include "MExtractSlidingWindow.h"
128#include "MExtractTimeSpline.h"
129#include "MFCosmics.h"
130#include "MContinue.h"
131#include "MFillH.h"
132
133#include "MArrivalTimeCam.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, Form("%s%s","Charge",(fRuns->GetRunsAsFileName()).Data()),
222 "Fitted Mean Charges");
223 MHCamera disp2 (geomcam, Form("%s%s","SigmaCharge",(fRuns->GetRunsAsFileName()).Data()),
224 "Sigma of Fitted Charges");
225 MHCamera disp3 (geomcam, Form("%s%s","RSigma",(fRuns->GetRunsAsFileName()).Data()),
226 "Reduced Sigmas");
227 MHCamera disp4 (geomcam, Form("%s%s","RSigmaPerCharge",(fRuns->GetRunsAsFileName()).Data()),
228 "Reduced Sigma per Charge");
229 MHCamera disp5 (geomcam, Form("%s%s","NumPhes",(fRuns->GetRunsAsFileName()).Data()),
230 "Nr. of Phe's (F-Factor Method)");
231 MHCamera disp6 (geomcam, Form("%s%s","ConvFADC2Phes",(fRuns->GetRunsAsFileName()).Data()),
232 "Conversion Factor (F-Factor Method)");
233 MHCamera disp7 (geomcam, Form("%s%s","TotalFFactor",(fRuns->GetRunsAsFileName()).Data()),
234 "Total F-Factor (F-Factor Method)");
235 MHCamera disp8 (geomcam, Form("%s%s","CascadesQEFFactor",(fRuns->GetRunsAsFileName()).Data()),
236 "Cascades QE (F-Factor Method)");
237 MHCamera disp9 (geomcam, Form("%s%s","CascadesQEBlindPix",(fRuns->GetRunsAsFileName()).Data()),
238 "Cascades QE (Blind Pixel Method)");
239 MHCamera disp10(geomcam, Form("%s%s","CascadesQEPINDiode",(fRuns->GetRunsAsFileName()).Data()),
240 "Cascades QE (PIN Diode Method)");
241 MHCamera disp11(geomcam, Form("%s%s","CascadesQECombined",(fRuns->GetRunsAsFileName()).Data()),
242 "Cascades QE (Combined Method)");
243 MHCamera disp12(geomcam, Form("%s%s","FFactorValid",(fRuns->GetRunsAsFileName()).Data()),
244 "Pixels with valid F-Factor calibration");
245 MHCamera disp13(geomcam, Form("%s%s","BlindPixelValid",(fRuns->GetRunsAsFileName()).Data()),
246 "Pixels with valid BlindPixel calibration");
247 MHCamera disp14(geomcam, Form("%s%s","PINdiodeValid",(fRuns->GetRunsAsFileName()).Data()),
248 "Pixels with valid PINDiode calibration");
249 MHCamera disp15(geomcam, Form("%s%s","CombinedValid",(fRuns->GetRunsAsFileName()).Data()),
250 "Pixels with valid Combined calibration");
251 MHCamera disp16(geomcam, Form("%s%s","Saturation",(fRuns->GetRunsAsFileName()).Data()),
252 "Pixels with saturated Hi Gain");
253 MHCamera disp17(geomcam, Form("%s%s","ConversionMeans",(fRuns->GetRunsAsFileName()).Data()),
254 "Conversion HiGain.vs.LoGain Means");
255 MHCamera disp18(geomcam, Form("%s%s","ConversionSigmas",(fRuns->GetRunsAsFileName()).Data()),
256 "Conversion HiGain.vs.LoGain Sigmas");
257 MHCamera disp19(geomcam, Form("%s%s","HiGainPickup",(fRuns->GetRunsAsFileName()).Data()),
258 "Number Pickup events Hi Gain");
259 MHCamera disp20(geomcam, Form("%s%s","LoGainPickup",(fRuns->GetRunsAsFileName()).Data()),
260 "Number Pickup events Lo Gain");
261 MHCamera disp21(geomcam, Form("%s%s","HiGainBlackout",(fRuns->GetRunsAsFileName()).Data()),
262 "Number Blackout events Hi Gain");
263 MHCamera disp22(geomcam, Form("%s%s","LoGainBlackout",(fRuns->GetRunsAsFileName()).Data()),
264 "Number Blackout events Lo Gain");
265 MHCamera disp23(geomcam, Form("%s%s","Excluded",(fRuns->GetRunsAsFileName()).Data()),
266 "Pixels previously excluded");
267 MHCamera disp24(geomcam, Form("%s%s","UnSuitable",(fRuns->GetRunsAsFileName()).Data()),
268 "Pixels not suited for further analysis");
269 MHCamera disp25(geomcam, Form("%s%s","UnReliable",(fRuns->GetRunsAsFileName()).Data()),
270 "Pixels not reliable for further analysis");
271 MHCamera disp26(geomcam, Form("%s%s","HiGainOscillating",(fRuns->GetRunsAsFileName()).Data()),
272 "Oscillating Pixels High Gain");
273 MHCamera disp27(geomcam, Form("%s%s","LoGainOscillating",(fRuns->GetRunsAsFileName()).Data()),
274 "Oscillating Pixels Low Gain");
275 MHCamera disp28(geomcam, Form("%s%s","AbsTimeMean",(fRuns->GetRunsAsFileName()).Data()),
276 "Abs. Arrival Times");
277 MHCamera disp29(geomcam, Form("%s%s","AbsTimeRms",(fRuns->GetRunsAsFileName()).Data()),
278 "RMS of Arrival Times");
279 MHCamera disp30(geomcam, Form("%s%s","MeanTime",(fRuns->GetRunsAsFileName()).Data()),
280 "Mean Rel. Arrival Times");
281 MHCamera disp31(geomcam, Form("%s%s","SigmaTime",(fRuns->GetRunsAsFileName()).Data()),
282 "Sigma Rel. Arrival Times");
283 MHCamera disp32(geomcam, Form("%s%s","TimeProb",(fRuns->GetRunsAsFileName()).Data()),
284 "Probability of Time Fit");
285 MHCamera disp33(geomcam, Form("%s%s","TimeNotFitValid",(fRuns->GetRunsAsFileName()).Data()),
286 "Pixels with not valid fit results");
287 MHCamera disp34(geomcam, Form("%s%s","TimeOscillating",(fRuns->GetRunsAsFileName()).Data()),
288 "Oscillating Pixels");
289
290 // Fitted charge means and sigmas
291 disp1.SetCamContent(fCalibrationCam, 0);
292 disp1.SetCamError( fCalibrationCam, 1);
293 disp2.SetCamContent(fCalibrationCam, 2);
294 disp2.SetCamError( fCalibrationCam, 3);
295
296 // Reduced Sigmas and reduced sigmas per charge
297 disp3.SetCamContent(fCalibrationCam, 5);
298 disp3.SetCamError( fCalibrationCam, 6);
299 disp4.SetCamContent(fCalibrationCam, 7);
300 disp4.SetCamError( fCalibrationCam, 8);
301
302 // F-Factor Method
303 disp5.SetCamContent(fCalibrationCam, 9);
304 disp5.SetCamError( fCalibrationCam, 10);
305 disp6.SetCamContent(fCalibrationCam, 11);
306 disp6.SetCamError( fCalibrationCam, 12);
307 disp7.SetCamContent(fCalibrationCam, 13);
308 disp7.SetCamError( fCalibrationCam, 14);
309
310 // Quantum Efficiencies
311 disp8.SetCamContent (fQECam, 0 );
312 disp8.SetCamError (fQECam, 1 );
313 disp9.SetCamContent (fQECam, 2 );
314 disp9.SetCamError (fQECam, 3 );
315 disp10.SetCamContent(fQECam, 4 );
316 disp10.SetCamError (fQECam, 5 );
317 disp11.SetCamContent(fQECam, 6 );
318 disp11.SetCamError (fQECam, 7 );
319
320 // Valid flags
321 disp12.SetCamContent(fQECam, 8 );
322 disp13.SetCamContent(fQECam, 9 );
323 disp14.SetCamContent(fQECam, 10);
324 disp15.SetCamContent(fQECam, 11);
325
326 // Conversion Hi-Lo
327 disp16.SetCamContent(fCalibrationCam, 25);
328 disp17.SetCamContent(fCalibrationCam, 16);
329 disp17.SetCamError (fCalibrationCam, 17);
330 disp18.SetCamContent(fCalibrationCam, 18);
331 disp18.SetCamError (fCalibrationCam, 19);
332
333 // Pickup and Blackout
334 disp19.SetCamContent(fCalibrationCam, 21);
335 disp20.SetCamContent(fCalibrationCam, 22);
336 disp21.SetCamContent(fCalibrationCam, 23);
337 disp22.SetCamContent(fCalibrationCam, 24);
338
339 // Pixels with defects
340 disp23.SetCamContent(fCalibrationCam, 20);
341 disp24.SetCamContent(fBadPixels, 1);
342 disp25.SetCamContent(fBadPixels, 3);
343
344 // Oscillations
345 disp26.SetCamContent(fBadPixels, 10);
346 disp27.SetCamContent(fBadPixels, 11);
347
348 // Arrival Times
349 disp28.SetCamContent(fCalibrationCam, 26);
350 disp28.SetCamError( fCalibrationCam, 27);
351 disp29.SetCamContent(fCalibrationCam, 27);
352
353 disp1.SetYTitle("Q [FADC counts]");
354 disp2.SetYTitle("\\sigma_{Q} [FADC counts]");
355
356 disp3.SetYTitle("\\sqrt{\\sigma^{2}_{Q} - RMS^{2}_{Ped}} [FADC Counts]");
357 disp4.SetYTitle("Red.Sigma/<Q> [1]");
358
359 disp5.SetYTitle("Nr. Phe's [1]");
360 disp6.SetYTitle("Conv.Factor [PhE/FADC counts]");
361 disp7.SetYTitle("Total F-Factor [1]");
362
363 disp8.SetYTitle("QE [1]");
364 disp9.SetYTitle("QE [1]");
365 disp10.SetYTitle("QE [1]");
366 disp11.SetYTitle("QE [1]");
367
368 disp12.SetYTitle("[1]");
369 disp13.SetYTitle("[1]");
370 disp14.SetYTitle("[1]");
371 disp15.SetYTitle("[1]");
372 disp16.SetYTitle("[1]");
373
374 disp17.SetYTitle("<Q>(High)/<Q>(Low) [1]");
375 disp18.SetYTitle("\\sigma_{Q}(High)/\\sigma_{Q}(Low) [1]");
376
377 disp19.SetYTitle("[1]");
378 disp20.SetYTitle("[1]");
379 disp21.SetYTitle("[1]");
380 disp22.SetYTitle("[1]");
381 disp23.SetYTitle("[1]");
382 disp24.SetYTitle("[1]");
383 disp25.SetYTitle("[1]");
384 disp26.SetYTitle("[1]");
385 disp27.SetYTitle("[1]");
386
387 disp28.SetYTitle("Mean Abs. Time [FADC slice]");
388 disp29.SetYTitle("RMS Abs. Time [FADC slices]");
389
390 if (fRelTimes)
391 {
392
393 disp30.SetCamContent(fRelTimeCam,0);
394 disp30.SetCamError( fRelTimeCam,1);
395 disp31.SetCamContent(fRelTimeCam,2);
396 disp31.SetCamError( fRelTimeCam,3);
397 disp32.SetCamContent(fRelTimeCam,4);
398 disp33.SetCamContent(fBadPixels,20);
399 disp34.SetCamContent(fBadPixels,21);
400
401 disp30.SetYTitle("Time Offset [FADC units]");
402 disp31.SetYTitle("Timing resolution [FADC units]");
403 disp32.SetYTitle("P_{Time} [1]");
404 disp33.SetYTitle("[1]");
405 disp34.SetYTitle("[1]");
406 }
407
408 if (fDisplayType == kDataCheckDisplay)
409 {
410 TCanvas &c1 = fDisplay->AddTab("Fit.Charge");
411 c1.Divide(3, 3);
412
413 CamDraw(c1, 1, 3, disp1, 2);
414 CamDraw(c1, 2, 3, disp4, 2);
415 CamDraw(c1, 3, 3, disp28, 2);
416
417 // F-Factor
418 TCanvas &c2 = fDisplay->AddTab("Phe's");
419 c2.Divide(3,4);
420
421 CamDraw(c2, 1, 3, disp6, 2, 1);
422 CamDraw(c2, 2, 3, disp7, 2, 1);
423 CamDraw(c2, 3, 3, disp8, 2, 1);
424
425 // QE's
426 TCanvas &c3 = fDisplay->AddTab("QE's");
427 c3.Divide(3,4);
428
429 CamDraw(c3, 1, 3, disp8, 2, 1);
430 CamDraw(c3, 2, 3, disp9, 2, 1);
431 CamDraw(c3, 3, 3, disp10, 2, 1);
432
433 // Defects
434 TCanvas &c4 = fDisplay->AddTab("Defect");
435 // c4.Divide(3,2);
436 c4.Divide(2,2);
437
438 // CamDraw(c4, 1, 3, disp23, 0);
439 // CamDraw(c4, 2, 3, disp24, 0);
440 // CamDraw(c4, 3, 3, disp25, 0);
441 CamDraw(c4, 1, 2, disp24, 0);
442 CamDraw(c4, 2, 2, disp25, 0);
443
444 if (fRelTimes)
445 {
446 // Rel. Times
447 TCanvas &c5 = fDisplay->AddTab("Rel. Times");
448 c5.Divide(2,4);
449
450 CamDraw(c5, 1, 2, disp30, 2);
451 CamDraw(c5, 2, 2, disp31, 2);
452 }
453
454
455 return;
456 }
457
458 if (fDisplayType == kNormalDisplay)
459 {
460
461 // Charges
462 TCanvas &c11 = fDisplay->AddTab("Fit.Charge");
463 c11.Divide(2, 4);
464
465 CamDraw(c11, 1, 2, disp1, 5, 1);
466 CamDraw(c11, 2, 2, disp2, 5, 1);
467
468 // Reduced Sigmas
469 TCanvas &c12 = fDisplay->AddTab("Red.Sigma");
470 c12.Divide(2,4);
471
472 CamDraw(c12, 1, 2, disp3, 5, 1);
473 CamDraw(c12, 2, 2, disp4, 5, 1);
474
475 // F-Factor
476 TCanvas &c13 = fDisplay->AddTab("Phe's");
477 c13.Divide(3,4);
478
479 CamDraw(c13, 1, 3, disp5, 5, 1);
480 CamDraw(c13, 2, 3, disp6, 5, 1);
481 CamDraw(c13, 3, 3, disp7, 5, 1);
482
483 // QE's
484 TCanvas &c14 = fDisplay->AddTab("QE's");
485 c14.Divide(4,4);
486
487 CamDraw(c14, 1, 4, disp8, 5, 1);
488 CamDraw(c14, 2, 4, disp9, 5, 1);
489 CamDraw(c14, 3, 4, disp10, 5, 1);
490 CamDraw(c14, 4, 4, disp11, 5, 1);
491
492 // Defects
493 TCanvas &c15 = fDisplay->AddTab("Defect");
494 // c15.Divide(5,2);
495 c15.Divide(4,2);
496
497 /*
498 CamDraw(c15, 1, 5, disp23, 0);
499 CamDraw(c15, 2, 5, disp24, 0);
500 CamDraw(c15, 3, 5, disp25, 0);
501 CamDraw(c15, 4, 5, disp26, 0);
502 CamDraw(c15, 5, 5, disp27, 0);
503 */
504 CamDraw(c15, 1, 4, disp24, 0);
505 CamDraw(c15, 2, 4, disp25, 0);
506 CamDraw(c15, 3, 4, disp26, 0);
507 CamDraw(c15, 4, 4, disp27, 0);
508
509 // Abs. Times
510 TCanvas &c16 = fDisplay->AddTab("Abs. Times");
511 c16.Divide(2,3);
512
513 CamDraw(c16, 1, 2, disp28, 5);
514 CamDraw(c16, 2, 2, disp29, 5);
515
516 if (fRelTimes)
517 {
518 // Rel. Times
519 TCanvas &c17 = fDisplay->AddTab("Rel. Times");
520 c17.Divide(2,4);
521
522 CamDraw(c17, 1, 2, disp30, 5, 1);
523 CamDraw(c17, 2, 2, disp31, 5, 1);
524 }
525
526 return;
527 }
528
529 if (fDisplayType == kFullDisplay)
530 {
531
532 MHCalibrationChargeBlindPix *blind = (MHCalibrationChargeBlindPix*)plist.FindObject("MHCalibrationChargeBlindPix");
533 blind->DrawClone();
534 gPad->GetCanvas()->SaveAs(Form("%s%s%s",blind->GetName(),(fRuns->GetRunsAsString()).Data(),".root"));
535
536 MHCalibrationChargeCam *ccm = (MHCalibrationChargeCam*)plist.FindObject("MHCalibrationChargeCam");
537 ccm->DrawClone();
538 gPad->GetCanvas()->SaveAs(Form("%s%s%s",ccm->GetName(),(fRuns->GetRunsAsString()).Data(),".root"));
539
540 MHCalibrationRelTimeCam *rel = (MHCalibrationRelTimeCam*)plist.FindObject("MHCalibrationRelTimeCam");
541 rel->DrawClone();
542 gPad->GetCanvas()->SaveAs(Form("%s%s%s",rel->GetName(),(fRuns->GetRunsAsString()).Data(),".root"));
543
544 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationChargeCam");
545
546 for (Int_t aidx=0;aidx<cam->GetAverageAreas();aidx++)
547 {
548 cam->GetAverageHiGainArea(aidx).DrawClone("all");
549 cam->GetAverageLoGainArea(aidx).DrawClone("all");
550 }
551
552 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
553 {
554 cam->GetAverageHiGainSector(sector).DrawClone("all");
555 cam->GetAverageLoGainSector(sector).DrawClone("all");
556 }
557
558 // Charges
559 TCanvas &c21 = fDisplay->AddTab("Fit.Charge");
560 c21.Divide(2, 4);
561
562 CamDraw(c21, 1, 2, disp1, 2, 1);
563 CamDraw(c21, 2, 2, disp2, 2, 1);
564
565 // Reduced Sigmas
566 TCanvas &c23 = fDisplay->AddTab("Red.Sigma");
567 c23.Divide(2,4);
568
569 CamDraw(c23, 1, 2, disp3, 2, 1);
570 CamDraw(c23, 2, 2, disp4, 2, 1);
571
572 // F-Factor
573 TCanvas &c24 = fDisplay->AddTab("Phe's");
574 c24.Divide(3,4);
575
576 CamDraw(c24, 1, 3, disp5, 2, 1);
577 CamDraw(c24, 2, 3, disp6, 2, 1);
578 CamDraw(c24, 3, 3, disp7, 2, 1);
579
580 // QE's
581 TCanvas &c25 = fDisplay->AddTab("QE's");
582 c25.Divide(4,4);
583
584 CamDraw(c25, 1, 4, disp8, 2, 1);
585 CamDraw(c25, 2, 4, disp9, 2, 1);
586 CamDraw(c25, 3, 4, disp10, 2, 1);
587 CamDraw(c25, 4, 4, disp11, 2, 1);
588
589 // Validity
590 TCanvas &c26 = fDisplay->AddTab("Valid");
591 c26.Divide(4,2);
592
593 CamDraw(c26, 1, 4, disp12, 0);
594 CamDraw(c26, 2, 4, disp13, 0);
595 CamDraw(c26, 3, 4, disp14, 0);
596 CamDraw(c26, 4, 4, disp15, 0);
597
598 // Other info
599 TCanvas &c27 = fDisplay->AddTab("HiLoGain");
600 c27.Divide(3,3);
601
602 CamDraw(c27, 1, 3, disp16, 0);
603 CamDraw(c27, 2, 3, disp17, 1);
604 CamDraw(c27, 3, 3, disp18, 1);
605
606 // Pickup
607 TCanvas &c28 = fDisplay->AddTab("Pickup");
608 c28.Divide(4,2);
609
610 CamDraw(c28, 1, 4, disp19, 0);
611 CamDraw(c28, 2, 4, disp20, 0);
612 CamDraw(c28, 3, 4, disp21, 0);
613 CamDraw(c28, 4, 4, disp22, 0);
614
615 // Defects
616 TCanvas &c29 = fDisplay->AddTab("Defect");
617 // c29.Divide(5,2);
618 c29.Divide(4,2);
619
620 CamDraw(c29, 1, 4, disp24, 0);
621 CamDraw(c29, 2, 4, disp25, 0);
622 CamDraw(c29, 3, 4, disp26, 0);
623 CamDraw(c29, 4, 4, disp27, 0);
624
625 // Abs. Times
626 TCanvas &c30 = fDisplay->AddTab("Abs. Times");
627 c30.Divide(2,3);
628
629 CamDraw(c30, 1, 2, disp28, 2);
630 CamDraw(c30, 2, 2, disp29, 1);
631
632 if (fRelTimes)
633 {
634 // Rel. Times
635 TCanvas &c31 = fDisplay->AddTab("Rel. Times");
636 c31.Divide(3,4);
637
638 CamDraw(c31, 1, 3, disp30, 2, 1);
639 CamDraw(c31, 2, 3, disp31, 2, 1);
640 CamDraw(c31, 3, 3, disp32, 4, 1);
641
642 // Time Defects
643 TCanvas &c32 = fDisplay->AddTab("Time Def.");
644 c32.Divide(2,2);
645
646 CamDraw(c32, 1, 2, disp33,0);
647 CamDraw(c32, 2, 2, disp34,0);
648
649 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationRelTimeCam");
650
651 for (Int_t aidx=0;aidx<cam->GetAverageAreas();aidx++)
652 {
653 cam->GetAverageHiGainArea(aidx).DrawClone("fourierevents");
654 cam->GetAverageLoGainArea(aidx).DrawClone("fourierevents");
655 }
656
657 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
658 {
659 cam->GetAverageHiGainSector(sector).DrawClone("fourierevents");
660 cam->GetAverageLoGainSector(sector).DrawClone("fourierevents");
661 }
662
663 }
664
665 return;
666 }
667}
668
669
670
671// --------------------------------------------------------------------------
672//
673// Find the colour of the pulsing LED:
674// - If the run number is smaller than gkIFAEBoxInaugurationRun, take MCalibrationCam::kCT1
675// - Otherwise find the colour out of the run name
676// - If no colour is found, return kFALSE
677//
678Bool_t MJCalibration::FindColor()
679{
680
681 const UInt_t nruns = fRuns->GetNumRuns();
682
683 if (nruns == 0)
684 return kFALSE;
685
686 TArrayI arr = fRuns->GetRuns();
687
688 if (arr[nruns-1] < gkIFAEBoxInaugurationRun)
689 {
690 *fLog << "Found colour kCT1 in runs: " << fRuns->GetRunsAsString() << endl;
691 fColor = MCalibrationCam::kCT1;
692 return kTRUE;
693 }
694
695 TString filenames;
696
697 while (!(filenames=((MDirIter*)fRuns)->Next()).IsNull())
698 {
699
700 filenames.ToLower();
701
702 if (filenames.Contains("green"))
703 if (fColor == MCalibrationCam::kNONE)
704 {
705 *fLog << "Found colour: kGREEN in " << filenames << endl;
706 fColor = MCalibrationCam::kGREEN;
707 }
708 else if (fColor != MCalibrationCam::kGREEN)
709 {
710 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
711 return kFALSE;
712 }
713
714 if (filenames.Contains("blue"))
715 if (fColor == MCalibrationCam::kNONE)
716 {
717 *fLog << "Found colour: kBLUE in " << filenames << endl;
718 fColor = MCalibrationCam::kBLUE;
719 }
720 else if (fColor != MCalibrationCam::kBLUE)
721 {
722 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
723 return kFALSE;
724 }
725
726 if (filenames.Contains("uv"))
727 if (fColor == MCalibrationCam::kNONE)
728 {
729 *fLog << "Found colour: kUV in " << filenames << endl;
730 fColor = MCalibrationCam::kUV;
731 }
732 else if (fColor != MCalibrationCam::kUV)
733 {
734 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
735 return kFALSE;
736 }
737
738 if (filenames.Contains("ct1"))
739 if (fColor == MCalibrationCam::kNONE)
740 {
741 *fLog << "Found colour: kCT1 in " << filenames << endl;
742 fColor = MCalibrationCam::kCT1;
743 }
744 else if (fColor != MCalibrationCam::kCT1)
745 {
746 *fLog << err << "Different colour found in " << filenames << "... abort" << endl;
747 return kFALSE;
748 }
749
750 }
751
752 if (fColor == MCalibrationCam::kNONE)
753 {
754 *fLog << "No colour found in filenames of runs: " << fRuns->GetRunsAsString()
755 << "... abort" << endl;
756 return kFALSE;
757 }
758
759 return kTRUE;
760}
761
762
763
764
765// --------------------------------------------------------------------------
766//
767// Retrieve the output file written by WriteResult()
768//
769TString MJCalibration::GetOutputFile() const
770{
771 if (!fRuns)
772 return "";
773
774 return Form("%s/%s-F1.root", (const char*)fOutputPath, (const char*)fRuns->GetRunsAsFileName());
775}
776
777
778// --------------------------------------------------------------------------
779//
780// Call the ProcessFile(MPedestalCam)
781//
782Bool_t MJCalibration::Process(MPedestalCam &pedcam)
783{
784 if (!ReadCalibrationCam())
785 return ProcessFile(pedcam);
786
787 return kTRUE;
788}
789
790// --------------------------------------------------------------------------
791//
792// Execute the task list and the eventloop:
793//
794// - Check if there are fRuns, otherwise return
795// - Check the colour of the files in fRuns (FindColor()), otherwise return
796// - Check for consistency between run numbers and number of files
797// - Add fRuns to MReadMarsFile
798// - Put into MParList:
799// 1) MPedestalCam (pedcam)
800// 2) MCalibrationQECam (fQECam)
801// 3) MCalibrationChargeCam (fCalibrationCam)
802// 4) MCalibrationRelTimeCam (fRelTimeCam) (only if flag fRelTimes is chosen)
803// 5) MBadPixelsCam (fBadPixels)
804// 6) MCalibrationChargePINDiode
805// 7) MCalibrationChargeBlindPix
806// - Put into the MTaskList:
807// 1) MReadMarsFile
808// 2) MBadPixelsMerge
809// 3) MGeomApply
810// 4) MExtractor
811// 5) MExtractPINDiode
812// 6) MExtractBlindPixel
813// 7) MExtractTime (only if flag fRelTimes is chosen)
814// 8) MContinue(MFCosmics)
815// 9) MFillH("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode")
816// 10) MFillH("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel")
817// 11) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam")
818// 12) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam")
819// 13) MCalibrationChargeCalc
820// 14) MFillH("MHCalibrationRelTimeCam", "MArrivalTimeCam") (only if flag fRelTimes is chosen)
821// 15) MCalibrationRelTimeCalc
822// - Execute MEvtLoop
823// - DisplayResult()
824// - WriteResult()
825//
826Bool_t MJCalibration::ProcessFile(MPedestalCam &pedcam)
827{
828 if (!fRuns)
829 {
830 *fLog << err << "No Runs choosen... abort." << endl;
831 return kFALSE;
832 }
833
834 if (fRuns->GetNumRuns() != fRuns->GetNumEntries())
835 {
836 *fLog << err << "Number of files found doesn't match number of runs... abort."
837 << fRuns->GetNumRuns() << " vs. " << fRuns->GetNumEntries() << endl;
838 return kFALSE;
839 }
840
841 *fLog << inf;
842 fLog->Separator(GetDescriptor());
843
844 if (!FindColor())
845 return kFALSE;
846
847 *fLog << "Calculate MCalibrationCam from Runs " << fRuns->GetRunsAsString() << endl;
848 *fLog << endl;
849
850 // Setup Tasklist
851 MParList plist;
852 MTaskList tlist;
853 plist.AddToList(&tlist);
854
855 MReadMarsFile read("Events");
856 MRawFileRead rawread("");
857
858 if (fDataCheck)
859 {
860// rawread.AddFiles(*fRuns);
861 tlist.AddToList(&rawread);
862 }
863 else
864 {
865 read.DisableAutoScheme();
866 static_cast<MRead&>(read).AddFiles(*fRuns);
867 tlist.AddToList(&read);
868 }
869
870 MCalibrationChargePINDiode pindiode;
871 MCalibrationChargeBlindPix blindpix;
872
873 plist.AddToList(&pedcam);
874 plist.AddToList(&fBadPixels);
875 plist.AddToList(&fQECam);
876 plist.AddToList(&fCalibrationCam);
877 plist.AddToList(&pindiode);
878 plist.AddToList(&blindpix);
879
880 MGeomApply apply;
881 // MBadPixelsMerge merge(&fBadPixels);
882 MExtractPINDiode pinext;
883 MExtractBlindPixel blindext;
884 MExtractSlidingWindow extract2;
885 MExtractTimeSpline timespline;
886 MCalibrationChargeCalc calcalc;
887 MCalibrationRelTimeCalc timecalc;
888
889 //
890 // As long as there are no DM's, have to colour by hand
891 //
892 calcalc.SetPulserColor(fColor);
893
894 MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode");
895 MFillH fillbnd("MHCalibrationChargeBlindPix", "MExtractedSignalBlindPixel");
896 MFillH fillcam("MHCalibrationChargeCam", "MExtractedSignalCam");
897 MFillH filltme("MHCalibrationRelTimeCam", "MArrivalTimeCam");
898 fillpin.SetNameTab("PINDiode");
899 fillbnd.SetNameTab("BlindPix");
900 fillcam.SetNameTab("Charge");
901 filltme.SetNameTab("RelTimes");
902
903
904 //
905 // Apply a filter against cosmics
906 // (will have to be needed in the future
907 // when the calibration hardware-trigger is working)
908 //
909 MFCosmics cosmics;
910 MContinue cont(&cosmics);
911
912 // tlist.AddToList(&merge);
913 tlist.AddToList(&apply);
914
915 if (fExtractor)
916 tlist.AddToList(fExtractor);
917 else
918 {
919 *fLog << warn << GetDescriptor()
920 << ": No extractor has been chosen, take default MExtractSlidingWindow " << endl;
921 tlist.AddToList(&extract2);
922 }
923
924
925 tlist.AddToList(&pinext);
926 tlist.AddToList(&blindext);
927
928 if (fRelTimes)
929 {
930 plist.AddToList(&fRelTimeCam);
931
932 if (fTimeExtractor)
933 tlist.AddToList(fTimeExtractor);
934 else
935 {
936 *fLog << warn << GetDescriptor()
937 << ": No extractor has been chosen, take default MTimeExtractSpline " << endl;
938 tlist.AddToList(&timespline);
939 }
940 }
941
942 tlist.AddToList(&cont);
943 tlist.AddToList(&fillcam);
944 tlist.AddToList(&fillpin);
945 tlist.AddToList(&fillbnd);
946 tlist.AddToList(&calcalc);
947
948 if (fRelTimes)
949 {
950 tlist.AddToList(&filltme);
951 tlist.AddToList(&timecalc);
952 }
953
954
955 // Create and setup the eventloop
956 MEvtLoop evtloop(fName);
957 evtloop.SetParList(&plist);
958 evtloop.SetDisplay(fDisplay);
959 evtloop.SetLogStream(fLog);
960
961 // Execute first analysis
962 if (!evtloop.Eventloop())
963 {
964 *fLog << err << GetDescriptor() << ": Failed." << endl;
965 return kFALSE;
966 }
967
968 tlist.PrintStatistics();
969
970 DisplayResult(plist);
971
972 if (!WriteResult())
973 return kFALSE;
974
975 *fLog << inf << GetDescriptor() << ": Done." << endl;
976
977 return kTRUE;
978}
979
980// --------------------------------------------------------------------------
981//
982// Read the following containers from GetOutputFile()
983// - MCalibrationChargeCam
984// - MCalibrationQECam
985// - MBadPixelsCam
986//
987Bool_t MJCalibration::ReadCalibrationCam()
988{
989 const TString fname = GetOutputFile();
990
991 if (gSystem->AccessPathName(fname, kFileExists))
992 {
993 *fLog << err << "Input file " << fname << " doesn't exist." << endl;
994 return kFALSE;
995 }
996
997 *fLog << inf << "Reading from file: " << fname << endl;
998
999 TFile file(fname, "READ");
1000 if (fCalibrationCam.Read()<=0)
1001 {
1002 *fLog << err << "Unable to read MCalibrationChargeCam from " << fname << endl;
1003 return kFALSE;
1004 }
1005
1006 if (fQECam.Read()<=0)
1007 {
1008 *fLog << err << "Unable to read MCalibrationQECam from " << fname << endl;
1009 return kFALSE;
1010 }
1011
1012 if (file.FindKey("MBadPixelsCam"))
1013 {
1014 MBadPixelsCam bad;
1015 if (bad.Read()<=0)
1016 {
1017 *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
1018 return kFALSE;
1019 }
1020 fBadPixels.Merge(bad);
1021 }
1022
1023 if (fDisplay /*&& !fDisplay->GetCanvas("Pedestals")*/) // FIXME!
1024 fDisplay->Read();
1025
1026 return kTRUE;
1027}
1028
1029
1030// --------------------------------------------------------------------------
1031//
1032// Set the path for output files, written by WriteResult()
1033//
1034void MJCalibration::SetOutputPath(const char *path)
1035{
1036 fOutputPath = path;
1037 if (fOutputPath.EndsWith("/"))
1038 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
1039}
1040
1041
1042// --------------------------------------------------------------------------
1043//
1044// Write the result into the output file GetOutputFile(), if fOutputPath exists.
1045//
1046// The following containers are written:
1047// - MStatusDisplay
1048// - MCalibrationChargeCam
1049// - MCalibrationQECam
1050// - MBadPixelsCam
1051//
1052Bool_t MJCalibration::WriteResult()
1053{
1054 if (fOutputPath.IsNull())
1055 return kTRUE;
1056
1057 const TString oname(GetOutputFile());
1058
1059 *fLog << inf << "Writing to file: " << oname << endl;
1060
1061 TFile file(oname, "UPDATE");
1062
1063 if (fDisplay && fDisplay->Write()<=0)
1064 {
1065 *fLog << err << "Unable to write MStatusDisplay to " << oname << endl;
1066 return kFALSE;
1067 }
1068
1069 if (fCalibrationCam.Write()<=0)
1070 {
1071 *fLog << err << "Unable to write MCalibrationChargeCam to " << oname << endl;
1072 return kFALSE;
1073 }
1074
1075 if (fQECam.Write()<=0)
1076 {
1077 *fLog << err << "Unable to write MCalibrationQECam to " << oname << endl;
1078 return kFALSE;
1079 }
1080
1081 if (fBadPixels.Write()<=0)
1082 {
1083 *fLog << err << "Unable to write MBadPixelsCam to " << oname << endl;
1084 return kFALSE;
1085 }
1086
1087 return kTRUE;
1088
1089}
1090
Note: See TracBrowser for help on using the repository browser.