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

Last change on this file since 6828 was 6828, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 70.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! Author(s): 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 MExtractTimeHighestIntegral
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// The absolute light calibration devices Blind Pixel and PIN Diode can be switched on
54// and off with the commands:
55//
56// - SetUseBlindPixel(Bool_t )
57// - SetUsePINDiode(Bool_t )
58//
59// See also: MHCalibrationChargePix, MHCalibrationChargeCam, MHGausEvents
60// MHCalibrationChargeBlindPix, MHCalibrationChargePINDiode
61// MCalibrationChargePix, MCalibrationChargeCam, MCalibrationChargeCalc
62// MCalibrationBlindPix, MCalibrationChargePINDiode,
63// MCalibrationQECam, MBadPixelsPix, MBadPixelsCam
64//
65// If the flag RelTimeCalibration() is set, a calibration of the relative arrival
66// times is also performed. The following containers (rectangular) and
67// tasks (ellipses) are called to produce an MCalibrationRelTimeCam used by
68// MCalibrateTime to correct timing offset between pixels: (MCalibrateTime is not
69// called from this class)
70//
71//Begin_Html
72/*
73<img src="images/RelTimeClasses.gif">
74*/
75//End_Html
76//
77// Different arrival time extractors can be set directly with the command
78// SetTimeExtractor(MExtractor *)
79//
80// Resource file entries are case sensitive!
81//
82// See also: MHCalibrationRelTimePix, MHCalibrationRelTimeCam, MHGausEvents
83// MCalibrationRelTimePix, MCalibrationRelTimeCam
84// MBadPixelsPix, MBadPixelsCam
85//
86/////////////////////////////////////////////////////////////////////////////
87#include "MJCalibration.h"
88
89#include <TFile.h>
90#include <TF1.h>
91#include <TStyle.h>
92#include <TCanvas.h>
93#include <TSystem.h>
94#include <TLine.h>
95#include <TLatex.h>
96#include <TLegend.h>
97#include <TRegexp.h>
98#include <TPaveText.h>
99#include <TPaveStats.h>
100#include <TEnv.h>
101
102#include "MLog.h"
103#include "MLogManip.h"
104
105#include "MRunIter.h"
106#include "MSequence.h"
107#include "MParList.h"
108#include "MTaskList.h"
109#include "MEvtLoop.h"
110
111#include "MHCamera.h"
112#include "MGeomCam.h"
113
114#include "MPedestalCam.h"
115#include "MCalibColorSteer.h"
116
117#include "MCalibrationIntensityChargeCam.h"
118#include "MCalibrationIntensityBlindCam.h"
119#include "MCalibrationIntensityRelTimeCam.h"
120#include "MCalibrationIntensityQECam.h"
121
122#include "MCalibrationPatternDecode.h"
123#include "MCalibrationCam.h"
124#include "MCalibrationHiLoCam.h"
125#include "MCalibrationHiLoPix.h"
126#include "MCalibrationQECam.h"
127#include "MCalibrationQEPix.h"
128#include "MCalibrationChargeCam.h"
129#include "MCalibrationChargePix.h"
130#include "MCalibrationChargePINDiode.h"
131#include "MCalibrationBlindPix.h"
132#include "MCalibrationBlindCam.h"
133#include "MCalibrationBlindCamOneOldStyle.h"
134#include "MCalibrationBlindCamTwoNewStyle.h"
135#include "MCalibrationBlindCamThreeNewStyle.h"
136#include "MCalibrationChargeCalc.h"
137#include "MCalibColorSet.h"
138#include "MCalibrationRelTimeCam.h"
139#include "MCalibrationRelTimeCalc.h"
140
141#include "MHGausEvents.h"
142#include "MHCalibrationCam.h"
143#include "MHCalibrationChargeCam.h"
144#include "MHCalibrationChargeBlindCam.h"
145#include "MHCalibrationChargePINDiode.h"
146#include "MHCalibrationRelTimeCam.h"
147#include "MHCalibrationPix.h"
148
149#include "MReadMarsFile.h"
150#include "MPedCalcPedRun.h"
151#include "MRawFileRead.h"
152#include "MGeomApply.h"
153#include "MTaskEnv.h"
154#include "MBadPixelsMerge.h"
155#include "MBadPixelsCam.h"
156#include "MExtractTime.h"
157#include "MExtractor.h"
158#include "MExtractPINDiode.h"
159#include "MExtractBlindPixel.h"
160#include "MExtractSlidingWindow.h"
161#include "MExtractTimeHighestIntegral.h"
162#include "MFCosmics.h"
163#include "MContinue.h"
164#include "MFillH.h"
165
166#include "MArrivalTimeCam.h"
167
168#include "MStatusDisplay.h"
169
170ClassImp(MJCalibration);
171
172using namespace std;
173
174const Int_t MJCalibration::gkIFAEBoxInaugurationRun = 20113;
175const Int_t MJCalibration::gkSecondBlindPixelInstallation = 31693;
176const Int_t MJCalibration::gkSpecialPixelsContInstallation = 34057;
177const Int_t MJCalibration::gkThirdBlindPixelInstallation = 43308;
178const TString MJCalibration::fgReferenceFile = "mjobs/calibrationref.rc";
179const TString MJCalibration::fgHiLoCalibFile = "mjobs/hilocalib_df4.root";
180// --------------------------------------------------------------------------
181//
182// Default constructor.
183//
184// - Sets fRuns to 0, fExtractor to NULL, fTimeExtractor to NULL, fColor to kNONE,
185// fDisplay to kNormalDisplay, kRelTimes to kFALSE, kataCheck to kFALSE, kDebug to kFALSE,
186// kIntensity to kFALSE
187// - SetUseBlindPixel()
188// - SetUsePINDiode()
189//
190MJCalibration::MJCalibration(const char *name, const char *title)
191 : fExtractor(NULL), fTimeExtractor(NULL),
192 fColor(MCalibrationCam::kNONE), fDisplayType(kDataCheckDisplay),
193 fGeometry("MGeomCamMagic")
194{
195
196 fName = name ? name : "MJCalibration";
197 fTitle = title ? title : "Tool to create the calibration constants for one calibration run";
198
199 SetUseBlindPixel(kFALSE);
200 SetUsePINDiode(kFALSE);
201
202 SetHiLoCalibration();
203 SetRelTimeCalibration();
204 SetDebug(kFALSE);
205 SetIntensity(kFALSE);
206
207 SetReferenceFile();
208 SetHiLoCalibFile();
209
210 fConvFADC2PheMin = 0.;
211 fConvFADC2PheMax = 1.5;
212 fConvFADC2PhotMin = 0.;
213 fConvFADC2PhotMax = 10.;
214 fQEMin = 0.;
215 fQEMax = 0.3;
216 fArrivalTimeMin = 1.;
217 fArrivalTimeMax = 10.;
218 fTimeOffsetMin = -3.;
219 fTimeOffsetMax = 3.;
220 fTimeResolutionMin = 0.;
221 fTimeResolutionMax = 1.;
222
223 fRefConvFADC2PheInner = 0.14;
224 fRefConvFADC2PheOuter = 0.4;
225 fRefConvFADC2PhotInner = 0.8;
226 fRefConvFADC2PhotOuter = 3.8;
227 fRefQEInner = 0.18;
228 fRefQEOuter = 0.12;
229 fRefArrivalTimeInner = 4.5;
230 fRefArrivalTimeOuter = 5.0;
231 fRefArrivalTimeRmsInner = 0.5;
232 fRefArrivalTimeRmsOuter = 0.5;
233 fRefTimeOffsetInner = -0.23;
234 fRefTimeOffsetOuter = 0.39;
235 fRefTimeResolutionInner = 0.153;
236 fRefTimeResolutionOuter = 0.128;
237
238}
239
240void MJCalibration::DrawTab(MParList &plist, const char *cont, const char *name, Option_t *opt)
241{
242 TObject *obj = plist.FindObject(cont);
243 if (!obj)
244 return;
245
246 fDisplay->AddTab(name);
247 obj->DrawClone(opt);
248}
249
250// --------------------------------------------------------------------------
251//
252// Display the results in MStatusDisplay:
253//
254// - Add "Calibration" to the MStatusDisplay title
255// - Retrieve the MGeomCam from MParList
256// - Initialize the following MHCamera's:
257// 1) MCalibrationPix::GetMean()
258// 2) MCalibrationPix::Sigma()
259// 3) MCalibrationChargePix::GetRSigma()
260// 4) MCalibrationChargePix::GetRSigmaPerCharge()
261// 5) MCalibrationChargePix::GetPheFFactorMethod()
262// 6) MCalibrationChargePix::GetMeanConvFADC2Phe()
263// 7) MCalibrationChargePix::GetMeanFFactorFADC2Phot()
264// 8) MCalibrationQEPix::GetQECascadesFFactor()
265// 9) MCalibrationQEPix::GetQECascadesBlindPixel()
266// 10) MCalibrationQEPix::GetQECascadesPINDiode()
267// 11) MCalibrationQEPix::GetQECascadesCombined()
268// 12) MCalibrationQEPix::IsAverageQEFFactorAvailable()
269// 13) MCalibrationQEPix::IsAverageQEBlindPixelAvailable()
270// 14) MCalibrationQEPix::IsAverageQEPINDiodeAvailable()
271// 15) MCalibrationQEPix::IsAverageQECombinedAvailable()
272// 16) MCalibrationChargePix::IsHiGainSaturation()
273// 17) MCalibrationPix::GetHiLoMeansDivided()
274// 18) MCalibrationPix::GetHiLoSigmasDivided()
275// 19) MCalibrationChargePix::GetHiGainPickup()
276// 20) MCalibrationChargePix::GetLoGainPickup()
277// 21) MCalibrationChargePix::GetHiGainBlackout()
278// 22) MCalibrationChargePix::GetLoGainBlackout()
279// 23) MCalibrationPix::IsExcluded()
280// 24) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnsuitableRun)
281// 25) MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnreliableRun)
282// 26) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kHiGainOscillating)
283// 27) MBadPixelsPix::IsUncalibrated(MBadPixelsPix::kLoGainOscillating)
284// 28) MCalibrationChargePix::GetAbsTimeMean()
285// 29) MCalibrationChargePix::GetAbsTimeRms()
286//
287// If the flag SetFullDisplay() is set, all MHCameras will be displayed.
288// if the flag SetDataCheckDisplay() is set, only the most important ones are displayed
289// and otherwise, (default: SetNormalDisplay()), a good selection of plots is given
290//
291void MJCalibration::DisplayResult(MParList &plist)
292{
293 if (!fDisplay)
294 return;
295
296 TString drawoption = "nonew ";
297 if (fDisplayType == kDataCheckDisplay)
298 drawoption += "datacheck";
299 if (fDisplayType == kFullDisplay)
300 drawoption += "all";
301
302 if (IsUsePINDiode())
303 DrawTab(plist, "MHCalibrationChargePINDiode", "PINDiode", drawoption);
304 if (IsUseBlindPixel())
305 DrawTab(plist, "MHCalibrationChargeBlindCam", "BlindPix", drawoption);
306 if (IsRelTimes())
307 DrawTab(plist, "MHCalibrationRelTimeCam", "Time", drawoption);
308 DrawTab(plist, "MHCalibrationChargeCam", "Charge", drawoption);
309
310 //
311 // Update display
312 //
313 TString title = fDisplay->GetTitle();
314 title += "-- Calibration ";
315 title += fSequence.IsValid() ? Form("calib%08d", fSequence.GetSequence()) : (const char*)fRuns->GetRunsAsString();
316 title += " --";
317 fDisplay->SetTitle(title);
318
319 //
320 // Get container from list
321 //
322 MGeomCam &geomcam = *(MGeomCam*)plist.FindObject("MGeomCam");
323
324 // Create histograms to display
325 MHCamera disp1 (geomcam, "Charge", "Fitted Mean Signal (Charges)");
326 MHCamera disp2 (geomcam, "SigmaCharge", "Sigma of Fitted Signal");
327 MHCamera disp3 (geomcam, "RSigma", "Reduced Sigmas");
328 MHCamera disp4 (geomcam, "RSigmaPerCharge", "Reduced Sigma per Charge");
329 MHCamera disp5 (geomcam, "NumPhes", "Number Photo-electrons");
330 MHCamera disp6 (geomcam, "ConvFADC2Phes", "Conversion Factor to Phes");
331 MHCamera disp7 (geomcam, "TotalFFactor", "Total F-Factor(F-Factor Method)");
332 MHCamera disp8 (geomcam, "CascadesQEFFactor", "Cascades QE (F-Factor Method)");
333 MHCamera disp9 (geomcam, "CascadesQEBlindPix","Cascades QE (Blind Pixel Method)");
334 MHCamera disp10(geomcam, "CascadesQEPINDiode","Cascades QE (PIN Diode Method)");
335 MHCamera disp11(geomcam, "CascadesQECombined","Cascades QE (Combined Method)");
336 MHCamera disp12(geomcam, "FFactorValid", "Pixels with valid F-Factor calibration");
337 MHCamera disp13(geomcam, "BlindPixelValid", "Pixels with valid BlindPixel calibration");
338 MHCamera disp14(geomcam, "PINdiodeValid", "Pixels with valid PINDiode calibration");
339 MHCamera disp15(geomcam, "CombinedValid", "Pixels with valid Combined calibration");
340 MHCamera disp16(geomcam, "Saturation", "Pixels with saturated Hi Gain");
341 MHCamera disp17(geomcam, "ConversionMeans", "Conversion HiGain.vs.LoGain Means");
342 MHCamera disp18(geomcam, "ConversionSigmas", "Conversion HiGain.vs.LoGain Sigmas");
343 MHCamera disp19(geomcam, "HiGainPickup", "Number Pickup events Hi Gain");
344 MHCamera disp20(geomcam, "LoGainPickup", "Number Pickup events Lo Gain");
345 MHCamera disp21(geomcam, "HiGainBlackout", "Number Blackout events Hi Gain");
346 MHCamera disp22(geomcam, "LoGainBlackout", "Number Blackout events Lo Gain");
347 MHCamera disp23(geomcam, "Excluded", "Pixels previously excluded");
348 MHCamera disp24(geomcam, "UnSuitable", "Pixels not suited for further analysis");
349 MHCamera disp25(geomcam, "UnReliable", "Pixels suitable, but not reliable for further analysis");
350 MHCamera disp26(geomcam, "HiGainOscillating", "Oscillating Pixels High Gain");
351 MHCamera disp27(geomcam, "LoGainOscillating", "Oscillating Pixels Low Gain");
352 MHCamera disp28(geomcam, "AbsTimeMean", "Abs. Arrival Times");
353 MHCamera disp29(geomcam, "AbsTimeRms", "RMS of Arrival Times");
354 MHCamera disp30(geomcam, "MeanTime", "Mean Rel. Arrival Times");
355 MHCamera disp31(geomcam, "SigmaTime", "Sigma Rel. Arrival Times");
356 MHCamera disp32(geomcam, "TimeProb", "Probability of Time Fit");
357 MHCamera disp33(geomcam, "TimeNotFitValid", "Pixels with not valid fit results");
358 MHCamera disp34(geomcam, "TimeOscillating", "Oscillating Pixels");
359 MHCamera disp35(geomcam, "TotalConv", "Conversion Factor to photons");
360 MHCamera disp36(geomcam, "RMSperMean", "Charge histogram RMS per Mean");
361
362 MCalibrationChargeCam *cam = NULL;
363 MCalibrationQECam *qecam = NULL;
364 MCalibrationRelTimeCam *relcam = NULL;
365 MBadPixelsCam *badcam = NULL;
366
367 if (IsIntensity())
368 {
369 cam = (MCalibrationChargeCam*) fIntensCalibCam.GetCam();
370 qecam = (MCalibrationQECam*) fIntensQECam.GetCam();
371 relcam = (MCalibrationRelTimeCam*)fIntensRelTimeCam.GetCam();
372 badcam = (MBadPixelsCam*) fIntensBadCam.GetCam();
373 }
374 else
375 {
376 cam = &fCalibrationCam;
377 qecam = &fQECam;
378 relcam = &fRelTimeCam;
379 badcam = &fBadPixels;
380 }
381
382 // Fitted charge means and sigmas
383 disp1.SetCamContent(*cam, 0);
384 disp1.SetCamError( *cam, 1);
385 disp2.SetCamContent(*cam, 2);
386 disp2.SetCamError( *cam, 3);
387
388 // Reduced Sigmas and reduced sigmas per charge
389 disp3.SetCamContent(*cam, 5);
390 disp3.SetCamError( *cam, 6);
391 disp4.SetCamContent(*cam, 7);
392 disp4.SetCamError( *cam, 8);
393
394 // F-Factor Method
395 disp5.SetCamContent(*cam, 9);
396 disp5.SetCamError( *cam, 10);
397 disp6.SetCamContent(*cam, 11);
398 disp6.SetCamError( *cam, 12);
399 disp7.SetCamContent(*cam, 13);
400 disp7.SetCamError( *cam, 14);
401
402 // Quantum Efficiencies
403 disp8.SetCamContent (*qecam, 0 );
404 disp8.SetCamError (*qecam, 1 );
405 disp9.SetCamContent (*qecam, 2 );
406 disp9.SetCamError (*qecam, 3 );
407 disp10.SetCamContent(*qecam, 4 );
408 disp10.SetCamError (*qecam, 5 );
409 disp11.SetCamContent(*qecam, 6 );
410 disp11.SetCamError (*qecam, 7 );
411
412 // Valid flags
413 disp12.SetCamContent(*qecam, 8 );
414 disp13.SetCamContent(*qecam, 9 );
415 disp14.SetCamContent(*qecam, 10);
416 disp15.SetCamContent(*qecam, 11);
417
418 // Conversion Hi-Lo
419 disp16.SetCamContent(*cam, 25);
420 disp17.SetCamContent(*cam, 16);
421 disp17.SetCamError (*cam, 17);
422 disp18.SetCamContent(*cam, 18);
423 disp18.SetCamError (*cam, 19);
424
425 // Pickup and Blackout
426 disp19.SetCamContent(*cam, 21);
427 disp20.SetCamContent(*cam, 22);
428 disp21.SetCamContent(*cam, 23);
429 disp22.SetCamContent(*cam, 24);
430
431 // Pixels with defects
432 disp23.SetCamContent(*cam, 20);
433 disp24.SetCamContent(*badcam, 6);
434 disp25.SetCamContent(*badcam, 7);
435
436 // Oscillations
437 disp26.SetCamContent(*badcam, 10);
438 disp27.SetCamContent(*badcam, 11);
439
440 // Arrival Times
441 disp28.SetCamContent(*cam, 26);
442 disp28.SetCamError( *cam, 27);
443 disp29.SetCamContent(*cam, 27);
444
445 // RMS and Mean
446 disp36.SetCamContent(*cam,32);
447
448 disp1.SetYTitle("Q [FADC cnts]");
449 disp2.SetYTitle("\\sigma_{Q} [FADC cnts]");
450
451 disp3.SetYTitle("\\sqrt{\\sigma^{2}_{Q} - RMS^{2}_{Ped}} [FADC cnts]");
452 disp4.SetYTitle("Red.Sigma/<Q> [1]");
453
454 disp5.SetYTitle("Photo-electons [1]");
455 disp6.SetYTitle("Conv.Factor [PhE/FADC cnts]");
456 disp7.SetYTitle("Total F-Factor [1]");
457
458 disp8.SetYTitle("QE [1]");
459 disp9.SetYTitle("QE [1]");
460 disp10.SetYTitle("QE [1]");
461 disp11.SetYTitle("QE [1]");
462
463 disp12.SetYTitle("[1]");
464 disp13.SetYTitle("[1]");
465 disp14.SetYTitle("[1]");
466 disp15.SetYTitle("[1]");
467 disp16.SetYTitle("[1]");
468
469 disp17.SetYTitle("<Q>(High)/<Q>(Low) [1]");
470 disp18.SetYTitle("\\sigma_{Q}(High)/\\sigma_{Q}(Low) [1]");
471
472 disp19.SetYTitle("[1]");
473 disp20.SetYTitle("[1]");
474 disp21.SetYTitle("[1]");
475 disp22.SetYTitle("[1]");
476 // disp23.SetYTitle("[1]");
477 // disp24.SetYTitle("[1]");
478 // disp25.SetYTitle("[1]");
479 disp26.SetYTitle("[1]");
480 disp27.SetYTitle("[1]");
481
482 disp28.SetYTitle("Mean Abs. Time [FADC sl.]");
483 disp29.SetYTitle("RMS Abs. Time [FADC sl.]");
484
485 disp35.SetYTitle("Conv.Factor [Ph/FADC cnts]");
486
487 disp36.SetYTitle("Charge RMS/<Q> [1]");
488
489 for (UInt_t i=0;i<geomcam.GetNumPixels();i++)
490 {
491
492 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*cam) [i];
493 MCalibrationQEPix &qe = (MCalibrationQEPix&) (*qecam)[i];
494
495 if (!pix.IsFFactorMethodValid())
496 continue;
497
498 const Float_t convphe = pix.GetMeanConvFADC2Phe();
499 const Float_t quaeff = qe.GetQECascadesFFactor(0.);
500
501 disp35.Fill(i,convphe/quaeff);
502 disp35.SetUsed(i);
503 }
504
505
506 if (IsRelTimes())
507 {
508 disp30.SetCamContent(*relcam,0);
509 disp30.SetCamError( *relcam,1);
510 disp31.SetCamContent(*relcam,2);
511 disp31.SetCamError( *relcam,3);
512 disp32.SetCamContent(*relcam,4);
513 disp33.SetCamContent(fBadPixels,20);
514 disp34.SetCamContent(fBadPixels,21);
515
516 disp30.SetYTitle("Time Offset [FADC units]");
517 disp31.SetYTitle("Timing resolution [FADC units]");
518 disp32.SetYTitle("P_{Time} [1]");
519 disp33.SetYTitle("[1]");
520 disp34.SetYTitle("[1]");
521 }
522
523 if (fDisplayType == kDataCheckDisplay)
524 {
525
526 TCanvas &c1 = fDisplay->AddTab("FitCharge");
527 c1.Divide(3, 3);
528
529 //
530 // MEAN CHARGES
531 //
532
533 c1.cd(1);
534 gPad->SetBorderMode(0);
535 gPad->SetTicks();
536 MHCamera *obj1=(MHCamera*)disp1.DrawCopy("hist");
537 //
538 // for the datacheck, fix the ranges!!
539 //
540 // obj1->SetMinimum(fChargeMin);
541 // obj1->SetMaximum(fChargeMax);
542 //
543 // Set the datacheck sizes:
544 //
545 FixDataCheckHist((TH1D*)obj1);
546 obj1->SetStats(kFALSE);
547 //
548 // set reference lines
549 //
550 // DisplayReferenceLines(obj1,0);
551
552 c1.cd(4);
553 gPad->SetBorderMode(0);
554 obj1->SetPrettyPalette();
555 obj1->Draw();
556
557 c1.cd(7);
558 gPad->SetBorderMode(0);
559 gPad->SetTicks();
560 TH1D *obj2 = (TH1D*)obj1->Projection(obj1->GetName());
561 obj2->Draw();
562 obj2->SetBit(kCanDelete);
563 obj2->Fit("gaus","Q");
564 TF1 *fun2 = obj2->GetFunction("gaus");
565 fun2->SetLineColor(kYellow);
566 gPad->Modified();
567 gPad->Update();
568 TPaveStats *st = (TPaveStats*)obj2->GetListOfFunctions()->FindObject("stats");
569 st->SetY1NDC(0.55);
570 st->SetY2NDC(0.89);
571 st->SetX1NDC(0.65);
572 st->SetX2NDC(0.99);
573 gPad->Modified();
574 gPad->Update();
575 //
576 // Set the datacheck sizes:
577 //
578 FixDataCheckHist(obj2);
579 obj2->SetStats(1);
580
581 //
582 // Display the outliers as dead and noisy pixels
583 //
584 DisplayOutliers(obj2,"low-ampl.","high-ampl.");
585 TLatex flattex;
586 flattex.SetTextSize(0.07);
587 const Double_t minl = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
588 const Double_t maxl = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
589 flattex.DrawLatex(minl+0.015*(maxl-minl),obj2->GetBinContent(obj2->GetMaximumBin())/1.35,
590 Form("Flatfield precision: %4.2f%%",
591 fun2->GetParameter(2)/fun2->GetParameter(1)*100.));
592
593 //
594 // RMS per Charge
595 //
596
597 c1.cd(2);
598 gPad->SetBorderMode(0);
599 gPad->SetTicks();
600 MHCamera *obj3=(MHCamera*)disp36.DrawCopy("hist");
601 //
602 // for the datacheck, fix the ranges!!
603 //
604 // obj3->SetMinimum(0.);
605 // obj3->SetMaximum(fChargeMax);
606 //
607 // Set the datacheck sizes:
608 //
609 FixDataCheckHist((TH1D*)obj3);
610 obj3->SetStats(kFALSE);
611 //
612 // set reference lines
613 //
614 // DisplayReferenceLines(obj3,0);
615
616 c1.cd(5);
617 gPad->SetBorderMode(0);
618 obj3->SetPrettyPalette();
619 obj3->Draw();
620
621 c1.cd(8);
622 gPad->SetBorderMode(0);
623 if (geomcam.InheritsFrom("MGeomCamMagic"))
624 DisplayDoubleProject(&disp36, "dead", "noisy");
625
626 //
627 // PHOTO ELECTRONS
628 //
629
630 c1.cd(3);
631 gPad->SetBorderMode(0);
632 gPad->SetTicks();
633 MHCamera *obj4=(MHCamera*)disp5.DrawCopy("hist");
634 //
635 // for the datacheck, fix the ranges!!
636 //
637 // obj3->SetMinimum(fChargeMin);
638 // obj3->SetMaximum(fChargeMax);
639 //
640 // Set the datacheck sizes:
641 //
642 FixDataCheckHist((TH1D*)obj4);
643 obj4->SetStats(kFALSE);
644 //
645 // set reference lines
646 //
647 // DisplayReferenceLines(obj3,0);
648
649 c1.cd(6);
650 gPad->SetBorderMode(0);
651 obj4->SetPrettyPalette();
652 obj4->Draw();
653
654 c1.cd(9);
655 gPad->SetBorderMode(0);
656 if (geomcam.InheritsFrom("MGeomCamMagic"))
657 DisplayDoubleProject(&disp5, "dead", "noisy");
658
659 //
660 // CONVERSION FACTORS
661 //
662 TCanvas &c2 = fDisplay->AddTab("Conversion");
663 c2.Divide(3,3);
664
665 c2.cd(1);
666 gPad->SetBorderMode(0);
667 gPad->SetTicks();
668 MHCamera *obj5=(MHCamera*)disp6.DrawCopy("hist");
669 //
670 // for the datacheck, fix the ranges!!
671 //
672 obj5->SetMinimum(fConvFADC2PheMin);
673 obj5->SetMaximum(fConvFADC2PheMax);
674 //
675 // Set the datacheck sizes:
676 //
677 FixDataCheckHist((TH1D*)obj5);
678 obj5->SetStats(kFALSE);
679 //
680 // set reference lines
681 //
682 DisplayReferenceLines(obj5,2);
683
684 c2.cd(4);
685 gPad->SetBorderMode(0);
686 obj5->SetPrettyPalette();
687 obj5->Draw();
688
689 c2.cd(7);
690 gPad->SetBorderMode(0);
691 if (geomcam.InheritsFrom("MGeomCamMagic"))
692 DisplayDoubleProject(&disp6, "noisy", "dead");
693
694 //
695 // QUANTUM EFFICIENCY
696 //
697 c2.cd(2);
698 gPad->SetBorderMode(0);
699 gPad->SetTicks();
700 MHCamera *obj6=(MHCamera*)disp8.DrawCopy("hist");
701 //
702 // for the datacheck, fix the ranges!!
703 //
704 obj6->SetMinimum(fQEMin);
705 obj6->SetMaximum(fQEMax);
706 //
707 // Set the datacheck sizes:
708 //
709 FixDataCheckHist((TH1D*)obj6);
710 obj6->SetStats(kFALSE);
711 //
712 // set reference lines
713 //
714 DisplayReferenceLines(obj6,0);
715
716 c2.cd(5);
717 gPad->SetBorderMode(0);
718 obj6->SetPrettyPalette();
719 obj6->Draw();
720
721 c2.cd(8);
722 gPad->SetBorderMode(0);
723 if (geomcam.InheritsFrom("MGeomCamMagic"))
724 DisplayDoubleProject(&disp8, "noisy", "dead");
725
726 //
727 // CONVERSION FADC TO PHOTONS
728 //
729
730 c2.cd(3);
731 gPad->SetBorderMode(0);
732 gPad->SetTicks();
733 MHCamera *obj7=(MHCamera*)disp35.DrawCopy("hist");
734 //
735 // for the datacheck, fix the ranges!!
736 //
737 obj7->SetMinimum(fConvFADC2PhotMin);
738 obj7->SetMaximum(fConvFADC2PhotMax);
739 //
740 // Set the datacheck sizes:
741 //
742 FixDataCheckHist((TH1D*)obj7);
743 obj7->SetStats(kFALSE);
744 //
745 // set reference lines
746 //
747 DisplayReferenceLines(obj7,1);
748
749 c2.cd(6);
750 gPad->SetBorderMode(0);
751 obj7->SetPrettyPalette();
752 obj7->Draw();
753 c2.cd(9);
754 gPad->SetBorderMode(0);
755 if (geomcam.InheritsFrom("MGeomCamMagic"))
756 DisplayDoubleProject(&disp35, "noisy", "dead");
757
758 //
759 // ARRIVAL TIMES
760 //
761 TCanvas &c3 = fDisplay->AddTab("AbsTimes");
762 c3.Divide(2,3);
763
764 c3.cd(1);
765 gPad->SetBorderMode(0);
766 gPad->SetTicks();
767 MHCamera *obj10=(MHCamera*)disp28.DrawCopy("hist");
768 //
769 // for the datacheck, fix the ranges!!
770 //
771 obj10->SetMinimum(fArrivalTimeMin);
772 obj10->SetMaximum(fArrivalTimeMax);
773 //
774 // Set the datacheck sizes:
775 //
776 FixDataCheckHist((TH1D*)obj10);
777 obj10->SetStats(kFALSE);
778 //
779 // set reference lines
780 //
781 DisplayReferenceLines(obj10,3);
782
783 c3.cd(3);
784 gPad->SetBorderMode(0);
785 obj10->SetPrettyPalette();
786 obj10->Draw();
787
788 c3.cd(5);
789 gPad->SetBorderMode(0);
790 if (geomcam.InheritsFrom("MGeomCamMagic"))
791 DisplayDoubleProject(&disp28, "early", "late");
792
793 //
794 // ARRIVAL TIMES JITTER
795 //
796 c3.cd(2);
797 gPad->SetBorderMode(0);
798 gPad->SetTicks();
799 MHCamera *obj11=(MHCamera*)disp29.DrawCopy("hist");
800 //
801 // for the datacheck, fix the ranges!!
802 //
803 // obj11->SetMinimum(fArrivalTimeMin);
804 // obj11->SetMaximum(fArrivalTimeMax);
805 //
806 // Set the datacheck sizes:
807 //
808 FixDataCheckHist((TH1D*)obj11);
809 obj11->SetStats(kFALSE);
810 //
811 // set reference lines
812 //
813 DisplayReferenceLines(obj11,4);
814
815 c3.cd(4);
816 gPad->SetBorderMode(0);
817 obj11->SetPrettyPalette();
818 obj11->Draw();
819
820 c3.cd(6);
821 gPad->SetBorderMode(0);
822 if (geomcam.InheritsFrom("MGeomCamMagic"))
823 DisplayDoubleProject(&disp29, "", "jittering");
824
825 //
826 // UNSUITABLE PIXELS
827 //
828 TCanvas &c4 = fDisplay->AddTab("Defect");
829 c4.Divide(2,2);
830
831 c4.cd(1);
832 gPad->SetBorderMode(0);
833 gPad->SetTicks();
834 MHCamera *obj8=(MHCamera*)disp24.DrawCopy("hist");
835 //
836 // for the datacheck, fix the ranges!!
837 //
838 const Double_t max = 11.;
839 obj8->SetMinimum(0.);
840 obj8->SetMaximum(11.);
841 //
842 // Set the datacheck sizes:
843 //
844 FixDataCheckHist((TH1D*)obj8);
845 obj8->SetStats(kFALSE);
846
847 gStyle->SetPalette(1);
848 const Int_t numcol = gStyle->GetNumberOfColors()-3;
849
850 TPaveText *pave = new TPaveText(0.0,0.0,0.99,0.99);
851 pave->SetBit(kCanDelete);
852 pave->ConvertNDCtoPad();
853 pave->SetTextSize(0.05);
854 pave->AddText(" ");
855 TText *t1 = pave->AddText(Form("Signal smaller 3 Pedestal RMS: %3i pixels",
856 CountBadPixels(&disp24,1)));
857 t1->SetTextColor(gStyle->GetColorPalette(Int_t(1./max*numcol + 1.)));
858 t1->SetTextAlign(12);
859 TText *t2 = pave->AddText(Form("%s%3i%s","Signal Rel. error too large: ",
860 CountBadPixels(&disp24,2)," pixels"));
861 t2->SetTextColor(gStyle->GetColorPalette(Int_t(2./max*numcol + 1.)));
862 t2->SetTextAlign(12);
863 TText *t4 = pave->AddText(Form("Low Gain Saturation: %3i pixels",
864 CountBadPixels(&disp24,3)));
865 t4->SetTextColor(gStyle->GetColorPalette(Int_t(3./max*numcol + 1.)));
866 t4->SetTextAlign(12);
867 TText *t5 = pave->AddText(Form("Mean Arr. Time In First Extraction Bin: %3i pixels",
868 CountBadPixels(&disp24,4)));
869 t5->SetTextColor(gStyle->GetColorPalette(Int_t(4./max*numcol + 1.)));
870 t5->SetTextAlign(12);
871 TText *t6 = pave->AddText(Form("Mean Arr. Time In Last 2 Extraction Bins: %3i pixels",
872 CountBadPixels(&disp24,5)));
873 t6->SetTextColor(gStyle->GetColorPalette(Int_t(5./max*numcol + 1.)));
874 t6->SetTextAlign(12);
875 TText *t9 = pave->AddText(Form("Deviating Number of Photons: %3i pixels",
876 CountBadPixels(&disp24,6)));
877 t9->SetTextColor(gStyle->GetColorPalette(Int_t(6./max*numcol + 1.)));
878 t9->SetTextAlign(12);
879 TText *t10= pave->AddText(Form("High-Gain Histogram Overflow: %3i pixels",
880 CountBadPixels(&disp24,7 )));
881 t10->SetTextColor(gStyle->GetColorPalette(Int_t(7./max*numcol + 1.)));
882 t10->SetTextAlign(12);
883 TText *t11= pave->AddText(Form("Low-Gain Histogram Overflow: %3i pixels",
884 CountBadPixels(&disp24,8 )));
885 t11->SetTextColor(gStyle->GetColorPalette(Int_t(8./max*numcol + 1.)));
886 t11->SetTextAlign(12);
887 TText *t12= pave->AddText(Form("Previously Excluded: %3i pixels",
888 CountBadPixels(&disp24,9)));
889 t12->SetTextColor(gStyle->GetColorPalette(Int_t(9./max*numcol + 1.)));
890 t12->SetTextAlign(12);
891 pave->Draw();
892
893 c4.cd(3);
894 gPad->SetBorderMode(0);
895 obj8->Draw();
896 obj8->SetPrettyPalette();
897
898 //
899 // UNRELIABLE PIXELS
900 //
901
902 c4.cd(2);
903 gPad->SetBorderMode(0);
904 gPad->SetTicks();
905 MHCamera *obj9=(MHCamera*)disp25.DrawCopy("hist");
906 //
907 // for the datacheck, fix the ranges!!
908 //
909 const Double_t max2 = 8.;
910 obj9->SetMinimum(0.);
911 obj9->SetMaximum(max2);
912 //
913 // Set the datacheck sizes:
914 //
915 FixDataCheckHist((TH1D*)obj9);
916 obj9->SetStats(kFALSE);
917
918 gStyle->SetPalette(1);
919
920 TPaveText *pave2 = new TPaveText(0.0,0.0,0.99,0.99);
921 pave2->SetBit(kCanDelete);
922 pave2->ConvertNDCtoPad();
923 pave2->SetTextSize(0.05);
924 pave2->AddText(" ");
925 TText *t3 = pave2->AddText(Form("Signal Sigma smaller Pedestal RMS: %3i pixels",
926 CountBadPixels(&disp25,1)));
927 t3->SetTextColor(gStyle->GetColorPalette(Int_t(1./max*numcol + 1.)));
928 t3->SetTextAlign(12);
929
930 TText *t7 = pave2->AddText(Form("Deviating Number of Photo-electrons: %3i pixels",
931 CountBadPixels(&disp25,2)));
932 t7->SetTextColor(gStyle->GetColorPalette(Int_t(2./max*numcol + 1.)));
933 t7->SetTextAlign(12);
934
935 TText *tt1 = pave2->AddText(Form("High Gain Signals could not be fitted: %3i pixels",
936 CountBadPixels(&disp25,3)));
937 tt1->SetTextColor(gStyle->GetColorPalette(Int_t(3./max2*numcol + 1.)));
938 tt1->SetTextAlign(12);
939 TText *tt2 = pave2->AddText(Form("Low Gain Signals could not be fitted: %3i pixels",
940 CountBadPixels(&disp25,4)));
941 tt2->SetTextColor(gStyle->GetColorPalette(Int_t(4./max2*numcol + 1.)));
942 tt2->SetTextAlign(12);
943 TText *tt3 = pave2->AddText(Form("Relative Arr. Times could not be fitted: %3i pixels",
944 CountBadPixels(&disp25,5)));
945 tt3->SetTextColor(gStyle->GetColorPalette(Int_t(5./max2*numcol + 1.)));
946 tt3->SetTextAlign(12);
947 TText *tt4 = pave2->AddText(Form("High Gain Signals Oscillation: %3i pixels",
948 CountBadPixels(&disp25,6)));
949 tt4->SetTextColor(gStyle->GetColorPalette(Int_t(6./max2*numcol + 1.)));
950 tt4->SetTextAlign(12);
951 TText *tt5 = pave2->AddText(Form("Low Gain Signals Oscillation: %3i pixels",
952 CountBadPixels(&disp25,7)));
953 tt5->SetTextColor(gStyle->GetColorPalette(Int_t(7./max2*numcol + 1.)));
954 tt5->SetTextAlign(12);
955 TText *tt6 = pave2->AddText(Form("Relative Arr. Times Oscillation: %3i pixels",
956 CountBadPixels(&disp25,8)));
957 tt6->SetTextColor(gStyle->GetColorPalette(Int_t(8./max2*numcol + 1.)));
958 tt6->SetTextAlign(12);
959 TText *tt8 = pave2->AddText(Form("Deviating global F-Factor: %3i pixels",
960 CountBadPixels(&disp25,9)));
961 tt8->SetTextColor(gStyle->GetColorPalette(Int_t(9./max2*numcol + 1.)));
962 tt8->SetTextAlign(12);
963 pave2->Draw();
964
965 c4.cd(4);
966 gPad->SetBorderMode(0);
967 obj9->SetPrettyPalette();
968 obj9->Draw();
969
970 if (IsRelTimes())
971 {
972 TCanvas &c5 = fDisplay->AddTab("RelTimes");
973 c5.Divide(2,3);
974
975 //
976 // MEAN REL. ARR. TIMES
977 //
978 c5.cd(1);
979 gPad->SetBorderMode(0);
980 gPad->SetTicks();
981 MHCamera *obj10=(MHCamera*)disp30.DrawCopy("hist");
982 //
983 // for the datacheck, fix the ranges!!
984 //
985 obj10->SetMinimum(fTimeOffsetMin);
986 obj10->SetMaximum(fTimeOffsetMax);
987 //
988 // Set the datacheck sizes:
989 //
990 FixDataCheckHist((TH1D*)obj10);
991 obj10->SetStats(kFALSE);
992 //
993 // set reference lines
994 //
995 DisplayReferenceLines(obj10,5);
996
997 c5.cd(3);
998 gPad->SetBorderMode(0);
999 obj10->SetPrettyPalette();
1000 obj10->Draw();
1001
1002 c5.cd(5);
1003 gPad->SetBorderMode(0);
1004 if (geomcam.InheritsFrom("MGeomCamMagic"))
1005 DisplayDoubleProject(&disp30, "early", "late");
1006
1007 //
1008 // JITTER Rel. Arr. Times
1009 //
1010 c5.cd(2);
1011 gPad->SetBorderMode(0);
1012 gPad->SetTicks();
1013 MHCamera *obj11=(MHCamera*)disp31.DrawCopy("hist");
1014 //
1015 // for the datacheck, fix the ranges!!
1016 //
1017 obj11->SetMinimum(fTimeResolutionMin);
1018 obj11->SetMaximum(fTimeResolutionMax);
1019 //
1020 // Set the datacheck sizes:
1021 //
1022 FixDataCheckHist((TH1D*)obj11);
1023 obj11->SetStats(kFALSE);
1024 //
1025 // set reference lines
1026 //
1027 DisplayReferenceLines(obj11,6);
1028
1029 c5.cd(4);
1030 gPad->SetBorderMode(0);
1031 obj11->SetPrettyPalette();
1032 obj11->Draw();
1033
1034 c5.cd(6);
1035 gPad->SetBorderMode(0);
1036 if (geomcam.InheritsFrom("MGeomCamMagic"))
1037 DisplayDoubleProject(&disp31, "too stable", "jittering");
1038
1039 }
1040 return;
1041 }
1042
1043 if (fDisplayType == kNormalDisplay)
1044 {
1045
1046 // Charges
1047 TCanvas &c11 = fDisplay->AddTab("FitCharge");
1048 c11.Divide(2, 4);
1049
1050 disp1.CamDraw(c11, 1, 2, 5, 1);
1051 disp2.CamDraw(c11, 2, 2, 5, 1);
1052
1053 // Reduced Sigmas
1054 TCanvas &c12 = fDisplay->AddTab("RedSigma");
1055 c12.Divide(2,4);
1056
1057 disp3.CamDraw(c12, 1, 2, 5, 1);
1058 disp4.CamDraw(c12, 2, 2, 5, 1);
1059
1060 // F-Factor
1061 TCanvas &c13 = fDisplay->AddTab("Phe's");
1062 c13.Divide(3,4);
1063
1064 disp5.CamDraw(c13, 1, 3, 5, 1);
1065 disp6.CamDraw(c13, 2, 3, 5, 1);
1066 disp7.CamDraw(c13, 3, 3, 5, 1);
1067
1068 // QE's
1069 TCanvas &c14 = fDisplay->AddTab("QE's");
1070 c14.Divide(4,4);
1071
1072 disp8.CamDraw(c14, 1, 4, 5, 1);
1073 disp9.CamDraw(c14, 2, 4, 5, 1);
1074 disp10.CamDraw(c14, 3, 4, 5, 1);
1075 disp11.CamDraw(c14, 4, 4, 5, 1);
1076
1077 // Defects
1078 TCanvas &c15 = fDisplay->AddTab("Defect");
1079 // c15.Divide(5,2);
1080 c15.Divide(4,2);
1081
1082 /*
1083 disp23.CamDraw(c15, 1, 5, 0);
1084 disp24.CamDraw(c15, 2, 5, 0);
1085 disp25.CamDraw(c15, 3, 5, 0);
1086 disp26.CamDraw(c15, 4, 5, 0);
1087 disp27.CamDraw(c15, 5, 5, 0);
1088 */
1089 disp24.CamDraw(c15, 1, 4, 0);
1090 disp25.CamDraw(c15, 2, 4, 0);
1091 disp26.CamDraw(c15, 3, 4, 0);
1092 disp27.CamDraw(c15, 4, 4, 0);
1093
1094 // Abs. Times
1095 TCanvas &c16 = fDisplay->AddTab("AbsTimes");
1096 c16.Divide(2,3);
1097
1098 disp28.CamDraw(c16, 1, 2, 5);
1099 disp29.CamDraw(c16, 2, 2, 5);
1100
1101 if (IsRelTimes())
1102 {
1103 // Rel. Times
1104 TCanvas &c17 = fDisplay->AddTab("RelTimes");
1105 c17.Divide(2,4);
1106
1107 disp30.CamDraw(c17, 1, 2, 5, 1);
1108 disp31.CamDraw(c17, 2, 2, 5, 1);
1109 }
1110
1111 return;
1112 }
1113
1114 if (fDisplayType == kFullDisplay)
1115 {
1116 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationChargeCam");
1117
1118 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
1119 {
1120 cam->GetAverageHiGainSector(sector).DrawClone("all");
1121 cam->GetAverageLoGainSector(sector).DrawClone("all");
1122 }
1123
1124 // Charges
1125 TCanvas &c21 = fDisplay->AddTab("FitCharge");
1126 c21.Divide(2, 4);
1127
1128 disp1.CamDraw(c21, 1, 2, 2, 1);
1129 disp2.CamDraw(c21, 2, 2, 2, 1);
1130
1131 // Reduced Sigmas
1132 TCanvas &c23 = fDisplay->AddTab("RedSigma");
1133 c23.Divide(2,4);
1134
1135 disp3.CamDraw(c23, 1, 2, 2, 1);
1136 disp4.CamDraw(c23, 2, 2, 2, 1);
1137
1138 // F-Factor
1139 TCanvas &c24 = fDisplay->AddTab("Phe's");
1140 c24.Divide(3,5);
1141
1142 disp5.CamDraw(c24, 1, 3, 2, 1, 1);
1143 disp6.CamDraw(c24, 2, 3, 2, 1, 1);
1144 disp7.CamDraw(c24, 3, 3, 2, 1, 1);
1145
1146 // QE's
1147 TCanvas &c25 = fDisplay->AddTab("QE's");
1148 c25.Divide(4,5);
1149
1150 disp8.CamDraw(c25, 1, 4, 2, 1, 1);
1151 disp9.CamDraw(c25, 2, 4, 2, 1, 1);
1152 disp10.CamDraw(c25, 3, 4, 2, 1, 1);
1153 disp11.CamDraw(c25, 4, 4, 2, 1, 1);
1154
1155 // Validity
1156 TCanvas &c26 = fDisplay->AddTab("Valid");
1157 c26.Divide(4,2);
1158
1159 disp12.CamDraw(c26, 1, 4, 0);
1160 disp13.CamDraw(c26, 2, 4, 0);
1161 disp14.CamDraw(c26, 3, 4, 0);
1162 disp15.CamDraw(c26, 4, 4, 0);
1163
1164 // Other info
1165 TCanvas &c27 = fDisplay->AddTab("HiLoGain");
1166 c27.Divide(3,3);
1167
1168 disp16.CamDraw(c27, 1, 3, 0);
1169 disp17.CamDraw(c27, 2, 3, 1);
1170 disp18.CamDraw(c27, 3, 3, 1);
1171
1172 // Pickup
1173 TCanvas &c28 = fDisplay->AddTab("Pickup");
1174 c28.Divide(4,2);
1175
1176 disp19.CamDraw(c28, 1, 4, 0);
1177 disp20.CamDraw(c28, 2, 4, 0);
1178 disp21.CamDraw(c28, 3, 4, 0);
1179 disp22.CamDraw(c28, 4, 4, 0);
1180
1181 // Defects
1182 TCanvas &c29 = fDisplay->AddTab("Defect");
1183 // c29.Divide(5,2);
1184 c29.Divide(4,2);
1185
1186 disp24.CamDraw(c29, 1, 4, 0);
1187 disp25.CamDraw(c29, 2, 4, 0);
1188 disp26.CamDraw(c29, 3, 4, 0);
1189 disp27.CamDraw(c29, 4, 4, 0);
1190
1191 // Abs. Times
1192 TCanvas &c30 = fDisplay->AddTab("AbsTimes");
1193 c30.Divide(2,3);
1194
1195 disp28.CamDraw(c30, 1, 2, 2);
1196 disp29.CamDraw(c30, 2, 2, 1);
1197
1198 if (IsRelTimes())
1199 {
1200 // Rel. Times
1201 TCanvas &c31 = fDisplay->AddTab("RelTimes");
1202 c31.Divide(3,5);
1203
1204 disp30.CamDraw(c31, 1, 3, 2, 1, 1);
1205 disp31.CamDraw(c31, 2, 3, 2, 1, 1);
1206 disp32.CamDraw(c31, 3, 3, 4, 1, 1);
1207
1208 // Time Defects
1209 TCanvas &c32 = fDisplay->AddTab("DefTime");
1210 c32.Divide(2,2);
1211
1212 disp33.CamDraw(c32, 1, 2, 0);
1213 disp34.CamDraw(c32, 2, 2, 0);
1214
1215 MHCalibrationCam *cam = (MHCalibrationCam*)plist.FindObject("MHCalibrationRelTimeCam");
1216
1217 for (Int_t sector=1;sector<cam->GetAverageSectors();sector++)
1218 {
1219 cam->GetAverageHiGainSector(sector).DrawClone("fourierevents");
1220 cam->GetAverageLoGainSector(sector).DrawClone("fourierevents");
1221 }
1222
1223 }
1224
1225 return;
1226 }
1227}
1228
1229
1230void MJCalibration::DisplayReferenceLines(MHCamera *cam, const Int_t what) const
1231{
1232
1233 const MGeomCam *geom = cam->GetGeometry();
1234
1235 Double_t x = geom->InheritsFrom("MGeomCamMagic") ? 397 : cam->GetNbinsX() ;
1236
1237 TLine line;
1238 line.SetLineStyle(kDashed);
1239 line.SetLineWidth(3);
1240 line.SetLineColor(kBlue);
1241
1242 TLine *l1 = NULL;
1243
1244 switch (what)
1245 {
1246 case 0:
1247 l1 = line.DrawLine(0, fRefQEInner, x, fRefQEInner);
1248 break;
1249 case 1:
1250 l1 = line.DrawLine(0, fRefConvFADC2PhotInner, x, fRefConvFADC2PhotInner);
1251 break;
1252 case 2:
1253 l1 = line.DrawLine(0, fRefConvFADC2PheInner, x, fRefConvFADC2PheInner );
1254 break;
1255 case 3:
1256 l1 = line.DrawLine(0, fRefArrivalTimeInner, x, fRefArrivalTimeInner );
1257 break;
1258 case 4:
1259 l1 = line.DrawLine(0, fRefArrivalTimeRmsInner, x, fRefArrivalTimeRmsInner );
1260 break;
1261 case 5:
1262 l1 = line.DrawLine(0, fRefTimeOffsetInner, x, fRefTimeOffsetInner );
1263 break;
1264 case 6:
1265 l1 = line.DrawLine(0, fRefTimeResolutionInner, x, fRefTimeResolutionInner );
1266 break;
1267 default:
1268 break;
1269 }
1270
1271 if (geom->InheritsFrom("MGeomCamMagic"))
1272 {
1273 const Double_t x2 = cam->GetNbinsX();
1274
1275 switch (what)
1276 {
1277 case 0:
1278 line.DrawLine(x2, fRefQEOuter, 398, fRefQEOuter);
1279 break;
1280 case 1:
1281 line.DrawLine(x2, fRefConvFADC2PhotOuter, 398, fRefConvFADC2PhotOuter );
1282 break;
1283 case 2:
1284 line.DrawLine(x2, fRefConvFADC2PheOuter, 398, fRefConvFADC2PheOuter);
1285 break;
1286 case 3:
1287 line.DrawLine(x2, fRefArrivalTimeOuter, 398, fRefArrivalTimeOuter);
1288 break;
1289 case 4:
1290 line.DrawLine(x2, fRefArrivalTimeRmsOuter, 398, fRefArrivalTimeRmsOuter);
1291 break;
1292 case 5:
1293 line.DrawLine(x2, fRefTimeOffsetOuter, 398, fRefTimeOffsetOuter);
1294 break;
1295 case 6:
1296 line.DrawLine(x2, fRefTimeResolutionOuter, 398, fRefTimeResolutionOuter);
1297 break;
1298 default:
1299 break;
1300 }
1301 }
1302
1303 TLegend *leg = new TLegend(0.6,0.85,0.9 ,0.95);
1304 leg->SetBit(kCanDelete);
1305 leg->AddEntry(l1, "Reference","l");
1306 leg->Draw();
1307}
1308
1309void MJCalibration::DisplayOutliers(TH1D *hist, const char* whatsmall, const char* whatbig) const
1310{
1311
1312 const Int_t kNotDraw = 1<<9;
1313 TF1 *f = hist->GetFunction("gaus");
1314 f->ResetBit(kNotDraw);
1315
1316 const Float_t mean = f->GetParameter(1);
1317 const Float_t lolim = mean - 4.0*f->GetParameter(2);
1318 const Float_t uplim = mean + 4.0*f->GetParameter(2);
1319 const Stat_t dead = hist->Integral(0,hist->FindBin(lolim)-1);
1320 const Stat_t noisy = hist->Integral(hist->FindBin(uplim)+1,hist->GetNbinsX()+1);
1321
1322 const Double_t max = hist->GetBinContent(hist->GetMaximumBin());
1323
1324 const Double_t minl = hist->GetBinCenter(hist->GetXaxis()->GetFirst());
1325 const Double_t maxl = hist->GetBinCenter(hist->GetXaxis()->GetLast());
1326
1327 TLatex deadtex;
1328 deadtex.SetTextSize(0.07);
1329 deadtex.DrawLatex(minl+0.015*(maxl-minl),max/1.1,
1330 Form("%3i %s pixels",(Int_t)dead,whatsmall));
1331
1332 TLatex noisytex;
1333 noisytex.SetTextSize(0.07);
1334 noisytex.DrawLatex(minl+0.015*(maxl-minl),max/1.2,
1335 Form("%3i %s pixels",(Int_t)noisy,whatbig));
1336
1337}
1338
1339void MJCalibration::FixDataCheckHist(TH1D *hist) const
1340{
1341
1342 hist->SetDirectory(NULL);
1343
1344 //
1345 // set the labels bigger
1346 //
1347 TAxis *xaxe = hist->GetXaxis();
1348 TAxis *yaxe = hist->GetYaxis();
1349
1350 xaxe->CenterTitle();
1351 yaxe->CenterTitle();
1352 xaxe->SetTitleSize(0.06);
1353 yaxe->SetTitleSize(0.06);
1354 xaxe->SetTitleOffset(0.8);
1355 yaxe->SetTitleOffset(0.85);
1356 xaxe->SetLabelSize(0.05);
1357 yaxe->SetLabelSize(0.05);
1358
1359}
1360
1361const Int_t MJCalibration::CountBadPixels ( MHCamera *cam , const Int_t what ) const
1362{
1363
1364 Int_t cnt = 0;
1365
1366 for (UInt_t pix=0; pix<cam->GetNumPixels();pix++)
1367 if ((Int_t)cam->GetPixContent(pix) == what)
1368 cnt++;
1369
1370 return cnt;
1371}
1372
1373// --------------------------------------------------------------------------
1374//
1375// Retrieve the output file written by WriteResult()
1376//
1377const char* MJCalibration::GetOutputFile() const
1378{
1379 const TString name(GetOutputFileName());
1380 if (name.IsNull())
1381 return "";
1382
1383 return Form("%s/%s", fPathOut.Data(), name.Data());
1384}
1385
1386const char* MJCalibration::GetOutputFileName() const
1387{
1388
1389 if (fSequence.IsValid())
1390 return Form("calib%08d.root", fSequence.GetSequence());
1391 if (!fRuns)
1392 return "";
1393
1394 return Form("%s-F1.root", (const char*)fRuns->GetRunsAsFileName());
1395}
1396
1397// --------------------------------------------------------------------------
1398//
1399// Read the following values from resource file:
1400//
1401// ConvFADC2PheMin
1402// ConvFADC2PheMax
1403// ConvFADC2PhotMin
1404// ConvFADC2PhotMax
1405//
1406// QEMin
1407// QEMax
1408//
1409// ArrivalTimeMin
1410// ArrivalTimeMax
1411//
1412// TimeOffsetMin
1413// TimeOffsetMax
1414// TimeResolutionMin
1415// TimeResolutionMax
1416//
1417// RefConvFADC2PheInner
1418// RefConvFADC2PheOuter
1419// RefConvFADC2PhotInner
1420// RefConvFADC2PhotOuter
1421//
1422// RefQEInner
1423// RefQEOuter
1424//
1425// RefArrivalTimeInner
1426// RefArrivalTimeOuter
1427// RefArrivalTimeRmsInner
1428// RefArrivalTimeRmsOuter
1429//
1430// RefTimeOffsetInner
1431// RefTimeOffsetOuter
1432// RefTimeResolutionInner
1433// RefTimeResolutionOuter
1434//
1435void MJCalibration::ReadReferenceFile()
1436{
1437 TEnv refenv(fReferenceFile);
1438
1439 fConvFADC2PheMin = refenv.GetValue("ConvFADC2PheMin",fConvFADC2PheMin);
1440 fConvFADC2PheMax = refenv.GetValue("ConvFADC2PheMax",fConvFADC2PheMax);
1441 fConvFADC2PhotMin = refenv.GetValue("ConvFADC2PhotMin",fConvFADC2PhotMin);
1442 fConvFADC2PhotMax = refenv.GetValue("ConvFADC2PhotMax",fConvFADC2PhotMax);
1443 fQEMin = refenv.GetValue("QEMin",fQEMin);
1444 fQEMax = refenv.GetValue("QEMax",fQEMax);
1445 fArrivalTimeMin = refenv.GetValue("ArrivalTimeMin",fArrivalTimeMin);
1446 fArrivalTimeMax = refenv.GetValue("ArrivalTimeMax",fArrivalTimeMax);
1447 fTimeOffsetMin = refenv.GetValue("TimeOffsetMin",fTimeOffsetMin);
1448 fTimeOffsetMax = refenv.GetValue("TimeOffsetMax",fTimeOffsetMax);
1449 fTimeResolutionMin = refenv.GetValue("TimeResolutionMin",fTimeResolutionMin);
1450 fTimeResolutionMax = refenv.GetValue("TimeResolutionMax",fTimeResolutionMax);
1451
1452 fRefConvFADC2PheInner = refenv.GetValue("RefConvFADC2PheInner",fRefConvFADC2PheInner);
1453 fRefConvFADC2PheOuter = refenv.GetValue("RefConvFADC2PheOuter",fRefConvFADC2PheOuter);
1454 fRefConvFADC2PhotInner = refenv.GetValue("RefConvFADC2PhotInner",fRefConvFADC2PhotInner);
1455 fRefConvFADC2PhotOuter = refenv.GetValue("RefConvFADC2PhotOuter",fRefConvFADC2PhotOuter);
1456 fRefQEInner = refenv.GetValue("RefQEInner",fRefQEInner);
1457 fRefQEOuter = refenv.GetValue("RefQEOuter",fRefQEOuter);
1458 fRefArrivalTimeInner = refenv.GetValue("RefArrivalTimeInner",fRefArrivalTimeInner);
1459 fRefArrivalTimeOuter = refenv.GetValue("RefArrivalTimeOuter",fRefArrivalTimeOuter);
1460 fRefArrivalTimeRmsInner = refenv.GetValue("RefArrivalTimeRmsInner",fRefArrivalTimeRmsInner);
1461 fRefArrivalTimeRmsOuter = refenv.GetValue("RefArrivalTimeRmsOuter",fRefArrivalTimeRmsOuter);
1462 fRefTimeOffsetInner = refenv.GetValue("RefTimeOffsetInner",fRefTimeOffsetInner);
1463 fRefTimeOffsetOuter = refenv.GetValue("RefTimeOffsetOuter",fRefTimeOffsetOuter);
1464 fRefTimeResolutionInner = refenv.GetValue("RefTimeResolutionInner",fRefTimeResolutionInner);
1465 fRefTimeResolutionOuter = refenv.GetValue("RefTimeResolutionOuter",fRefTimeResolutionOuter);
1466}
1467
1468// --------------------------------------------------------------------------
1469//
1470// Read the following values from resource file:
1471//
1472// MCalibrationHiLoCam
1473//
1474Bool_t MJCalibration::ReadHiLoCalibFile()
1475{
1476
1477 if (!fIsHiLoCalibration)
1478 return kTRUE;
1479
1480 TFile file(fHiLoCalibFile,"READ");
1481
1482 MCalibrationHiLoCam hilocam;
1483
1484 *fLog << all << "Initializing High gain vs. Low gain intercalibration from " << fHiLoCalibFile << endl;
1485 *fLog << all << endl;
1486
1487 if (hilocam.Read()<=0)
1488 {
1489 *fLog << err << "Unable to read MCalibrationHiLoCam from " << fHiLoCalibFile << endl;
1490 return kFALSE;
1491 }
1492
1493 if (hilocam.GetSize() < 1)
1494 {
1495 *fLog << err << "MCalibationHiLoCam is un-initialized in file " << fHiLoCalibFile << endl;
1496 return kFALSE;
1497 }
1498
1499 if (fCalibrationCam.GetSize() < 1)
1500 fCalibrationCam.InitSize(hilocam.GetSize());
1501
1502 if (fBadPixels.GetSize() < 1)
1503 fBadPixels.InitSize(hilocam.GetSize());
1504
1505 if (fCalibrationCam.GetSize() != hilocam.GetSize())
1506 {
1507 *fLog << err << "Size mismatch MCalibationHiLoCam and MCalibrationChargeCam " << endl;
1508 return kFALSE;
1509 }
1510
1511 for (Int_t i=0;i<hilocam.GetSize();i++)
1512 {
1513 const MCalibrationHiLoPix &pix = (MCalibrationHiLoPix&)hilocam[i];
1514
1515 const Float_t ratio = pix.GetHiLoChargeRatio();
1516 const Float_t sigma = pix.GetHiLoChargeRatioSigma();
1517
1518 if (ratio < 0.)
1519 {
1520 fBadPixels[i].SetUncalibrated(MBadPixelsPix::kConversionHiLoNotValid);
1521 continue;
1522 }
1523
1524 MCalibrationChargePix &cpix = (MCalibrationChargePix&)fCalibrationCam[i];
1525
1526 cpix.SetConversionHiLo(ratio);
1527 cpix.SetConversionHiLoErr(sigma);
1528 }
1529
1530 return kTRUE;
1531}
1532
1533// --------------------------------------------------------------------------
1534//
1535// MJCalibration allows to setup several option by a resource file:
1536// MJCalibration.Display: full, datacheck, normal
1537// MJCalibration.RelTimeCalibration: yes,no
1538// MJCalibration.DataCheck: yes,no
1539// MJCalibration.Debug: yes,no
1540// MJCalibration.Intensity: yes,no
1541// MJCalibration.UseBlindPixel: yes,no
1542// MJCalibration.UsePINDiode: yes,no
1543// MJCalibration.Geometry: MGeomCamMagic, MGeomCamECO1000
1544//
1545// Name of a file containing reference values (see ReadReferenceFile)
1546// Prefix.ReferenceFile: filename
1547// (see ReadReferenceFile)
1548//
1549// For more details see the class description and the corresponding Getters
1550//
1551Bool_t MJCalibration::CheckEnvLocal()
1552{
1553
1554 TString dis = GetEnv("Display", "");
1555 if (dis.BeginsWith("Full", TString::kIgnoreCase))
1556 SetFullDisplay();
1557 if (dis.BeginsWith("DataCheck", TString::kIgnoreCase))
1558 SetDataCheckDisplay();
1559 if (dis.BeginsWith("Normal", TString::kIgnoreCase))
1560 SetNormalDisplay();
1561
1562 if (!MJCalib::CheckEnvLocal())
1563 return kFALSE;
1564
1565 SetRelTimeCalibration(GetEnv("RelTimeCalibration", IsRelTimes()));
1566 SetIntensity(GetEnv("IntensityCalibration", IsIntensity()));
1567 SetDebug(GetEnv("Debug", IsDebug()));
1568
1569 SetUseBlindPixel(GetEnv("UseBlindPixel", IsUseBlindPixel()));
1570 SetUsePINDiode(GetEnv("UsePINDiode", IsUsePINDiode()));
1571 SetGeometry(GetEnv("Geometry", fGeometry));
1572
1573 fReferenceFile = GetEnv("ReferenceFile",fReferenceFile.Data());
1574 ReadReferenceFile();
1575
1576 fHiLoCalibFile = GetEnv("HiLoCalibFile",fHiLoCalibFile.Data());
1577
1578 return ReadHiLoCalibFile();
1579}
1580
1581// --------------------------------------------------------------------------
1582//
1583// Call the ProcessFile(MPedestalCam)
1584//
1585Bool_t MJCalibration::Process(MPedestalCam &pedcam)
1586{
1587 if (!ReadCalibrationCam())
1588 return ProcessFile(pedcam);
1589
1590 return kTRUE;
1591}
1592
1593void MJCalibration::InitBlindPixel(MExtractBlindPixel &blindext,
1594 MHCalibrationChargeBlindCam &blindcam)
1595{
1596
1597 Int_t run = fSequence.IsValid() ? fSequence.GetLastRun() : fRuns->GetRuns()[fRuns->GetNumRuns()-1];
1598
1599 //
1600 // Initialize the blind pixel. Unfortunately, there is a hardware difference
1601 // in the first blind pixel until run "gkSecondBlindPixelInstallation" and the
1602 // later setup. The first needs to use a filter because of the length of
1603 // spurious NSB photon signals. The latter get better along extracting the amplitude
1604 // from a small window.
1605 //
1606 if (run < gkSecondBlindPixelInstallation)
1607 {
1608 MCalibrationBlindCamOneOldStyle blindresults;
1609 if (IsIntensity())
1610 blindresults.Copy(*fIntensBlindCam.GetCam());
1611 else
1612 blindresults.Copy(fCalibrationBlindCam);
1613
1614 blindext.SetExtractionType(MExtractBlindPixel::kIntegral);
1615 blindext.SetExtractionType(MExtractBlindPixel::kFilter);
1616 blindext.SetRange(10,19,0,6);
1617 blindext.SetNSBFilterLimit(70);
1618 }
1619 else if (run < gkThirdBlindPixelInstallation)
1620 {
1621
1622 MCalibrationBlindCamTwoNewStyle blindresults;
1623
1624 if (IsIntensity())
1625 blindresults.Copy(*fIntensBlindCam.GetCam());
1626 else
1627 blindresults.Copy(fCalibrationBlindCam);
1628
1629 blindext.SetNumBlindPixels(blindresults.GetSize());
1630 for (Int_t i=0;i<blindresults.GetSize();i++)
1631 blindext.SetBlindPixelIdx(blindresults[i].GetPixId(),i);
1632
1633 blindext.SetExtractionType(MExtractBlindPixel::kAmplitude);
1634 blindext.SetExtractionType(MExtractBlindPixel::kFilter);
1635 blindext.SetRange(5,8,0,2);
1636 blindext.SetNSBFilterLimit(38);
1637
1638 }
1639 else
1640 {
1641
1642 MCalibrationBlindCamThreeNewStyle blindresults;
1643
1644 if (IsIntensity())
1645 blindresults.Copy(*fIntensBlindCam.GetCam());
1646 else
1647 blindresults.Copy(fCalibrationBlindCam);
1648
1649 blindext.SetNumBlindPixels(blindresults.GetSize());
1650
1651 for (Int_t i=0;i<blindresults.GetSize();i++)
1652 blindext.SetBlindPixelIdx(blindresults[i].GetPixId(),i);
1653
1654 blindext.SetExtractionType(MExtractBlindPixel::kAmplitude);
1655 blindext.SetExtractionType(MExtractBlindPixel::kFilter);
1656 blindext.SetDataType(MExtractBlindPixel::kRawEvt2);
1657 blindext.SetRange(5,8,0,2);
1658 blindext.SetNSBFilterLimit(38);
1659
1660 }
1661
1662}
1663
1664// --------------------------------------------------------------------------
1665//
1666// Execute the task list and the eventloop:
1667//
1668// - Check if there are fRuns, otherwise return
1669// - Check the colour of the files in fRuns (FindColor()), otherwise return
1670// - Check for consistency between run numbers and number of files
1671// - Add fRuns to MReadMarsFile
1672// - Put into MParList:
1673// 1) MPedestalCam (pedcam)
1674// 2) MCalibrationQECam (fQECam)
1675// 3) MCalibrationChargeCam (fCalibrationCam)
1676// 4) MCalibrationRelTimeCam (fRelTimeCam) (only if flag IsRelTimes() is chosen)
1677// 5) MBadPixelsCam (fBadPixels)
1678// 6) MCalibrationChargePINDiode
1679// 7) MCalibrationBlindPix
1680// - Put into the MTaskList:
1681// 1) MReadMarsFile
1682// 2) MBadPixelsMerge
1683// 3) MGeomApply
1684// 4) MExtractor
1685// 5) MExtractPINDiode
1686// 6) MExtractBlindPixel
1687// 7) MExtractTime (only if flag IsRelTimes() is chosen)
1688// 8) MContinue(MFCosmics)
1689// 9) MFillH("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode", "FillPINDiode")
1690// 10) MFillH("MHCalibrationChargeBlindCam", "MExtractedSignalBlindPixel", "FillBlindCam")
1691// 11) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam", "FillChargeCam")
1692// 12) MFillH("MHCalibrationChargeCam", "MExtractedSignalCam", "FillRelTime")
1693// 13) MCalibrationChargeCalc
1694// 14) MFillH("MHCalibrationRelTimeCam", "MArrivalTimeCam") (only if flag IsRelTimes() is chosen)
1695// 15) MCalibrationRelTimeCalc
1696// - Execute MEvtLoop
1697// - DisplayResult()
1698// - WriteResult()
1699//
1700Bool_t MJCalibration::ProcessFile(MPedestalCam &pedcam)
1701{
1702 if (!fSequence.IsValid())
1703 {
1704 if (!fRuns)
1705 {
1706 *fLog << err << "No Runs choosen... abort." << endl;
1707 return kFALSE;
1708 }
1709
1710 if (fRuns->GetNumRuns() != fRuns->GetNumEntries())
1711 {
1712 *fLog << err << "Number of files found doesn't match number of runs... abort."
1713 << fRuns->GetNumRuns() << " vs. " << fRuns->GetNumEntries() << endl;
1714 return kFALSE;
1715 }
1716 }
1717
1718 // --------------------------------------------------------------------------------
1719
1720 *fLog << inf;
1721 fLog->Separator(GetDescriptor());
1722
1723 *fLog << "Calculate MCalibrationCam from ";
1724 if (fSequence.IsValid())
1725 *fLog << "Sequence #" << fSequence.GetSequence() << endl;
1726 else
1727 *fLog << "Runs " << fRuns->GetRunsAsString() << endl;
1728 *fLog << endl;
1729
1730 // --------------------------------------------------------------------------------
1731
1732 if (!CheckEnv())
1733 return kFALSE;
1734
1735 // --------------------------------------------------------------------------------
1736
1737 // Setup Tasklist
1738 MParList plist;
1739 MTaskList tlist;
1740 plist.AddToList(&tlist);
1741 plist.AddToList(this); // take care of fDisplay!
1742
1743 MDirIter iter;
1744 if (fSequence.IsValid())
1745 {
1746 const Int_t n0 = fSequence.SetupCalRuns(iter, fPathData, "C", IsUseRawData());
1747 const Int_t n1 = fSequence.GetNumCalRuns();
1748 if (n0==0)
1749 {
1750 *fLog << err << "ERROR - No input files of sequence found!" << endl;
1751 return kFALSE;
1752 }
1753 if (n0!=n1)
1754 {
1755 *fLog << err << "ERROR - Number of files found ("
1756 << n0 << ") doesn't match number of files in sequence ("
1757 << n1 << ")" << endl;
1758 return kFALSE;
1759 }
1760 }
1761
1762 //
1763 // Input containers
1764 //
1765 pedcam.SetName("MPedestalCam"); // MPedestalFundamental
1766 plist.AddToList(&pedcam);
1767 plist.AddToList(&fBadPixels);
1768
1769 //
1770 // Calibration Results containers
1771 //
1772 if (IsIntensity())
1773 {
1774 plist.AddToList(&fIntensQECam);
1775 plist.AddToList(&fIntensCalibCam);
1776 plist.AddToList(&fIntensBlindCam);
1777 // plist.AddToList(&fIntensCalibrationPINDiode);
1778 plist.AddToList(&fIntensRelTimeCam);
1779 plist.AddToList(&fIntensBadCam);
1780 }
1781 else
1782 {
1783 plist.AddToList(&fQECam);
1784 plist.AddToList(&fCalibrationCam);
1785 plist.AddToList(&fCalibrationBlindCam);
1786 plist.AddToList(&fCalibrationPINDiode);
1787 plist.AddToList(&fRelTimeCam);
1788 }
1789
1790 //
1791 // Initialize two histogram containers which could be modified in this class
1792 //
1793 MHCalibrationRelTimeCam reltimecam;
1794 MHCalibrationChargeCam chargecam;
1795 MHCalibrationChargeBlindCam blindcam;
1796 plist.AddToList(&chargecam);
1797
1798 if (IsUseBlindPixel())
1799 plist.AddToList(&blindcam);
1800 if (IsRelTimes())
1801 plist.AddToList(&reltimecam);
1802 //
1803 // Data Reading tasks
1804 //
1805 MReadMarsFile read("Events");
1806 MRawFileRead rawread(NULL);
1807
1808 if (IsUseRawData())
1809 {
1810 rawread.AddFiles(fSequence.IsValid() ? iter : *fRuns);
1811 tlist.AddToList(&rawread);
1812 }
1813 else
1814 {
1815 read.DisableAutoScheme();
1816 read.AddFiles(fSequence.IsValid() ? iter : *fRuns);
1817 tlist.AddToList(&read);
1818 }
1819
1820 //
1821 // Other Tasks
1822 //
1823 MCalibrationPatternDecode decode;
1824 MGeomApply apply;
1825 apply.SetGeometry(fGeometry);
1826
1827 MBadPixelsMerge merge(&fBadPixels);
1828 MExtractPINDiode pinext;
1829 MExtractBlindPixel blindext;
1830
1831 InitBlindPixel(blindext, blindcam);
1832
1833 MExtractSlidingWindow extract2;
1834 MExtractTimeHighestIntegral timehigh;
1835 MCalibrationChargeCalc calcalc;
1836 MCalibrationRelTimeCalc timecalc;
1837 calcalc.SetOutputFile("");
1838 timecalc.SetOutputFile("");
1839
1840 if (!fSequence.IsValid())
1841 {
1842 calcalc.SetOutputPath(fPathOut);
1843 calcalc.SetOutputFile(Form("%s-ChargeCalibStat.txt",(const char*)fRuns->GetRunsAsFileName()));
1844 timecalc.SetOutputPath(fPathOut);
1845 timecalc.SetOutputFile(Form("%s-ChargeCalibStat.txt",(const char*)fRuns->GetRunsAsFileName()));
1846 }
1847
1848 if (IsDebug())
1849 {
1850 chargecam.SetDebug();
1851 calcalc.SetDebug();
1852 }
1853
1854 //
1855 // Calibration histogramming
1856 //
1857 MFillH fillpin("MHCalibrationChargePINDiode", "MExtractedSignalPINDiode", "FillPINDiode");
1858 MFillH fillbnd("MHCalibrationChargeBlindCam", "MExtractedSignalBlindPixel", "FillBlindCam");
1859 MFillH fillcam("MHCalibrationChargeCam", "MExtractedSignalCam", "FillChargeCam");
1860 MFillH filltme("MHCalibrationRelTimeCam", "MArrivalTimeCam", "FillRelTime");
1861 fillpin.SetBit(MFillH::kDoNotDisplay);
1862 fillbnd.SetBit(MFillH::kDoNotDisplay);
1863 fillcam.SetBit(MFillH::kDoNotDisplay);
1864 filltme.SetBit(MFillH::kDoNotDisplay);
1865
1866 //
1867 // Set default extractors in case, none has been set...
1868 //
1869 if (!fExtractor)
1870 fExtractor = &extract2;
1871 if (!fTimeExtractor)
1872 fTimeExtractor = &timehigh;
1873
1874 const Bool_t istimecharge = fExtractor->InheritsFrom("MExtractTimeAndCharge");
1875 //
1876 // Look if the extractor is a pure charge or also a time extractor
1877 //
1878 if (istimecharge)
1879 {
1880 if (fExtractorCam.GetSize() == pedcam.GetSize())
1881 calcalc.SetPedestals(&fExtractorCam);
1882 else
1883 {
1884 *fLog << err << GetDescriptor() << "ERROR - ";
1885 *fLog << "Used Extractor derives from MExtractTimeAndCharge, " << endl;
1886 *fLog << "but MExtractorCam size " << fExtractorCam.GetSize() << " ";
1887 *fLog << "mismatch pedcam size " << pedcam.GetSize() << "! " << endl;
1888 return kFALSE;
1889 }
1890 }
1891
1892 //
1893 // Setup more tasks and tasklist
1894 //
1895 MTaskEnv taskenv("ExtractSignal");
1896 taskenv.SetDefault(fExtractor);
1897
1898 tlist.AddToList(&decode);
1899 tlist.AddToList(&merge);
1900 tlist.AddToList(&apply);
1901
1902 MPedCalcPedRun pedcalc;
1903 pedcalc.SetExtractWindow(fExtractor->GetHiGainFirst(),TMath::Nint(fExtractor->GetNumHiGainSamples()));
1904
1905 if (IsIntensity())
1906 tlist.AddToList(&pedcalc);
1907
1908 MCalibColorSet colorset;
1909 tlist.AddToList(&colorset);
1910
1911 tlist.AddToList(&taskenv);
1912
1913 if (IsUsePINDiode())
1914 tlist.AddToList(&pinext);
1915 if (IsUseBlindPixel())
1916 tlist.AddToList(&blindext);
1917
1918 MTaskEnv taskenv2("ExtractTime");
1919 if (!istimecharge)
1920 {
1921 taskenv2.SetDefault(fTimeExtractor);
1922
1923 if (IsRelTimes())
1924 tlist.AddToList(&taskenv2);
1925 }
1926
1927 //
1928 // Apply a filter against cosmics
1929 // (will have to be needed in the future
1930 // when the calibration hardware-trigger is working)
1931 //
1932 MFCosmics cosmics;
1933 MContinue cont(&cosmics);
1934
1935 if (fColor == MCalibrationCam::kCT1)
1936 tlist.AddToList(&cont);
1937
1938 MCalibColorSteer steer;
1939 if (IsIntensity())
1940 tlist.AddToList(&steer);
1941
1942 tlist.AddToList(&fillcam);
1943
1944 if (IsRelTimes())
1945 {
1946 tlist.AddToList(&filltme);
1947 tlist.AddToList(&timecalc);
1948 }
1949
1950 if (IsUseBlindPixel())
1951 tlist.AddToList(&fillbnd);
1952 if (IsUsePINDiode())
1953 tlist.AddToList(&fillpin);
1954
1955 tlist.AddToList(&calcalc);
1956
1957 // Create and setup the eventloop
1958 MEvtLoop evtloop(fName);
1959 evtloop.SetParList(&plist);
1960 evtloop.SetDisplay(fDisplay);
1961 evtloop.SetLogStream(fLog);
1962 if (!SetupEnv(evtloop))
1963 return kFALSE;
1964
1965 if (!taskenv.GetTask() && !taskenv2.GetTask())
1966 {
1967 *fLog << err << "ERROR - Neither ExtractSignal nor ExtractTime initialized or both '<dummy>'." << endl;
1968 return kFALSE;
1969 }
1970
1971 if (!WriteTasks(taskenv.GetTask(), istimecharge ? 0 : taskenv2.GetTask()))
1972 return kFALSE;
1973
1974 // Execute first analysis
1975 if (!evtloop.Eventloop())
1976 {
1977 DisplayResult(plist);
1978 WriteResult(plist);
1979 *fLog << err << GetDescriptor() << ": Failed." << endl;
1980 return kFALSE;
1981 }
1982
1983 tlist.PrintStatistics();
1984
1985 // --------------------------------------------------------------------------------
1986
1987 if (fIsPixelCheck)
1988 {
1989 MHCalibrationChargeCam *hcam = (MHCalibrationChargeCam*)plist.FindObject("MHCalibrationChargeCam");
1990 MHCalibrationChargePix &pix1 = (MHCalibrationChargePix&)(*hcam)[fCheckedPixId];
1991 pix1.DrawClone("");
1992
1993 MHCalibrationChargePix &pix2 = (MHCalibrationChargePix&)(*hcam)(fCheckedPixId);
1994 pix2.DrawClone("");
1995
1996 if (IsRelTimes())
1997 {
1998 MHCalibrationRelTimeCam *hccam = (MHCalibrationRelTimeCam*)plist.FindObject("MHCalibrationRelTimeCam");
1999 MHCalibrationPix &pix11 = (*hccam)[fCheckedPixId];
2000 pix11.DrawClone("");
2001 MHCalibrationPix &pix12 = (*hccam)(fCheckedPixId);
2002 pix12.DrawClone("");
2003 }
2004 }
2005
2006 if (!fCalibrationPINDiode.IsValid())
2007 SetUsePINDiode(kFALSE);
2008
2009 DisplayResult(plist);
2010
2011 if (!WriteResult(plist))
2012 return kFALSE;
2013
2014 *fLog << all << GetDescriptor() << ": Done." << endl;
2015
2016 return kTRUE;
2017}
2018
2019// --------------------------------------------------------------------------
2020//
2021// Read the following containers from GetOutputFile()
2022// - MCalibrationChargeCam
2023// - MCalibrationQECam
2024// - MBadPixelsCam
2025//
2026Bool_t MJCalibration::ReadCalibrationCam()
2027{
2028
2029 if (IsNoStorage())
2030 return kFALSE;
2031
2032 const TString fname = GetOutputFile();
2033
2034 if (gSystem->AccessPathName(fname, kFileExists))
2035 {
2036 *fLog << err << "Input file " << fname << " doesn't exist." << endl;
2037 return kFALSE;
2038 }
2039
2040 *fLog << inf << "Reading from file: " << fname << endl;
2041
2042 TFile file(fname, "READ");
2043 if (fCalibrationCam.Read()<=0)
2044 {
2045 *fLog << err << "Unable to read MCalibrationChargeCam from " << fname << endl;
2046 return kFALSE;
2047 }
2048
2049 if (fQECam.Read()<=0)
2050 {
2051 *fLog << err << "Unable to read MCalibrationQECam from " << fname << endl;
2052 return kFALSE;
2053 }
2054
2055
2056 if (file.FindKey("MCalibrationRelTimeCam"))
2057 if (fRelTimeCam.Read()<=0)
2058 {
2059 *fLog << err << "Unable to read MCalibrationRelTimeCam from " << fname << endl;
2060 return kFALSE;
2061 }
2062
2063 if (file.FindKey("MBadPixelsCam"))
2064 {
2065 MBadPixelsCam bad;
2066 if (bad.Read()<=0)
2067 {
2068 *fLog << err << "Unable to read MBadPixelsCam from " << fname << endl;
2069 return kFALSE;
2070 }
2071 fBadPixels.Merge(bad);
2072 }
2073
2074 if (fDisplay /*&& !fDisplay->GetCanvas("Pedestals")*/) // FIXME!
2075 fDisplay->Read();
2076
2077 return kTRUE;
2078}
2079
2080
2081// --------------------------------------------------------------------------
2082//
2083// Set the useage of the Blind Pixel device
2084//
2085void MJCalibration::SetUseBlindPixel(const Bool_t b)
2086{
2087 b ? SETBIT(fDevices,kUseBlindPixel) : CLRBIT(fDevices,kUseBlindPixel);
2088}
2089
2090// --------------------------------------------------------------------------
2091//
2092// Set the useage of the PIN Diode device
2093//
2094void MJCalibration::SetUsePINDiode(const Bool_t b)
2095{
2096 b ? SETBIT(fDevices,kUsePINDiode) : CLRBIT(fDevices,kUsePINDiode);
2097}
2098
2099/*
2100Bool_t MJCalibration::WriteEventloop(MEvtLoop &evtloop) const
2101{
2102 if (IsNoStorage())
2103 return kTRUE;
2104
2105 if (fPathOut.IsNull())
2106 return kTRUE;
2107
2108 const TString oname(GetOutputFile());
2109
2110 *fLog << inf << "Writing to file: " << oname << endl;
2111
2112 TFile file(oname, fOverwrite?"RECREATE":"NEW", "File created by MJCalibration", 9);
2113 if (!file.IsOpen())
2114 {
2115 *fLog << err << "ERROR - Couldn't open file " << oname << " for writing..." << endl;
2116 return kFALSE;
2117 }
2118
2119 if (evtloop.Write(fName)<=0)
2120 {
2121 *fLog << err << "Unable to write MEvtloop to " << oname << endl;
2122 return kFALSE;
2123 }
2124
2125 return kTRUE;
2126}
2127*/
2128
2129Bool_t MJCalibration::WriteTasks(MTask *t1, MTask *t2) const
2130{
2131 if (IsNoStorage())
2132 return kTRUE;
2133
2134 TObjArray cont;
2135 if (t1)
2136 cont.Add(t1);
2137 if (t2)
2138 cont.Add(t2);
2139
2140 return WriteContainer(cont, GetOutputFileName(), fOverwrite?"RECREATE":"NEW");
2141}
2142
2143// --------------------------------------------------------------------------
2144//
2145// Write the result into the output file GetOutputFile(), if fOutputPath exists.
2146//
2147// The following containers are written:
2148// - MStatusDisplay
2149// - MCalibrationChargeCam or MCalibrationIntensityChargeCam
2150// - MCalibrationBlindCam or MCalibrationIntensityBlindCam
2151// - MCalibrationQECam or MCalibrationIntensityQECam
2152// - MCalibrationChargePINDiode
2153// - MBadPixelsCam
2154// If the flag kRelTimes is set, then also:
2155// - MCalibrationRelTimeCam or MCalibrationIntensityRelTimeCam
2156//
2157Bool_t MJCalibration::WriteResult(MParList &plist)
2158{
2159 TObjArray cont;
2160
2161 if (!IsNoStorage())
2162 {
2163 if (IsIntensity())
2164 {
2165 cont.Add(&fIntensBadCam);
2166 cont.Add(&fIntensCalibCam);
2167 cont.Add(&fIntensQECam);
2168 cont.Add(&fIntensBlindCam);
2169 }
2170 else
2171 {
2172 cont.Add(&fBadPixels);
2173 cont.Add(&fCalibrationCam);
2174 cont.Add(&fQECam);
2175 cont.Add(&fCalibrationBlindCam);
2176 }
2177 cont.Add(&fCalibrationPINDiode);
2178
2179 //if (IsRelTimes())
2180 cont.Add(IsIntensity() ? (TObject*)&fIntensRelTimeCam : (TObject*)&fRelTimeCam);
2181
2182 if (fExtractorCam.GetSize() != 0)
2183 cont.Add(&fExtractorCam);
2184
2185 TObject *pedcam = plist.FindObject("MPedestalCam");
2186 if (!pedcam)
2187 *fLog << warn << " - WARNING - MPedestalCam (fundamental)... not found for writing!" << endl;
2188 else
2189 cont.Add(pedcam);
2190
2191 TObject *geom = plist.FindObject("MGeomCam");
2192 if (!geom)
2193 *fLog << warn << " - WARNING - MGeomCam... not found for writing!" << endl;
2194 else
2195 cont.Add(geom);
2196 }
2197
2198 if (IsHistsStorage())
2199 {
2200 cont.Add(plist.FindObject("MHCalibrationChargeCam"));
2201 cont.Add(plist.FindObject("MHCalibrationChargeBlindCam"));
2202 cont.Add(plist.FindObject("MHCalibrationChargePINDiode"));
2203 if (IsRelTimes())
2204 cont.Add(plist.FindObject("MHCalibrationRelTimeCam"));
2205 }
2206
2207 return WriteContainer(cont, GetOutputFileName(), "UPDATE");
2208}
2209
2210void MJCalibration::DisplayDoubleProject(MHCamera *cam, const char* whatsmall, const char* whatbig) const
2211{
2212
2213 TArrayI inner(1);
2214 inner[0] = 0;
2215
2216 TArrayI outer(1);
2217 outer[0] = 1;
2218
2219 TArrayI s1(3);
2220 s1[0] = 6;
2221 s1[1] = 1;
2222 s1[2] = 2;
2223
2224 TArrayI s2(3);
2225 s2[0] = 3;
2226 s2[1] = 4;
2227 s2[2] = 5;
2228
2229 TVirtualPad *pad = gPad;
2230 pad->Divide(2,1);
2231
2232 TH1D *inout[2];
2233
2234 for (int i=0; i<2; i++)
2235 {
2236 pad->cd(i+1);
2237 gPad->SetBorderMode(0);
2238 gPad->SetTicks();
2239
2240 inout[i] = cam->ProjectionS(TArrayI(), TArrayI(1,&i), i==0 ? "Inner" : "Outer");
2241 FixDataCheckHist(inout[i]);
2242 inout[i]->SetTitle(Form("%s %s",cam->GetTitle(),i==0 ? "Inner" : "Outer"));
2243 inout[i]->SetDirectory(NULL);
2244 inout[i]->SetLineColor(kRed+i);
2245 inout[i]->SetBit(kCanDelete);
2246 inout[i]->Draw();
2247 //
2248 // Display the outliers as dead and noisy pixels
2249 //
2250 if (!inout[i]->Fit("gaus","0Q"))
2251 DisplayOutliers(inout[i],whatsmall,whatbig);
2252
2253 gPad->Modified();
2254 gPad->Update();
2255 TPaveStats *st = (TPaveStats*)inout[i]->GetListOfFunctions()->FindObject("stats");
2256 st->SetY1NDC(0.6);
2257 st->SetY2NDC(0.9);
2258 st->SetX1NDC(0.55);
2259 st->SetX2NDC(0.99);
2260 gPad->Modified();
2261 gPad->Update();
2262
2263 TLegend *leg2 = new TLegend(0.55,0.4,0.99,0.6);
2264
2265 //
2266 // Display the two half of the camera separately
2267 //
2268 TH1D *half[2];
2269 half[0] = cam->ProjectionS(s1, TArrayI(1,&i), "Sector 6-1-2");
2270 half[1] = cam->ProjectionS(s2, TArrayI(1,&i), "Sector 3-4-5");
2271
2272 for (int j=0; j<2; j++)
2273 {
2274 half[j]->SetLineColor(kRed+i+2*j+1);
2275 half[j]->SetDirectory(NULL);
2276 half[j]->SetBit(kCanDelete);
2277 half[j]->Draw("same");
2278 leg2->AddEntry(half[j], half[j]->GetName(), "l");
2279 }
2280 leg2->Draw();
2281 }
2282}
Note: See TracBrowser for help on using the repository browser.