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

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