source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationRelTimeCam.cc@ 8300

Last change on this file since 8300 was 8192, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 25.7 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): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHCalibrationRelTimeCam
28//
29// Fills the extracted relative arrival times of MArrivalTimeCam into
30// the MHCalibrationPix-classes MHCalibrationPix for every:
31//
32// - Pixel, stored in the TObjArray's MHCalibrationCam::fHiGainArray
33// or MHCalibrationCam::fHiGainArray, respectively, depending if
34// MArrivalTimePix::IsLoGainUsed() is set.
35//
36// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
37// stored in the TObjArray's MHCalibrationCam::fAverageHiGainAreas and
38// MHCalibrationCam::fAverageHiGainAreas
39//
40// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
41// stored in the TObjArray's MHCalibrationCam::fAverageHiGainSectors
42// and MHCalibrationCam::fAverageHiGainSectors
43//
44// Every relative time is calculated as the difference between the individual
45// pixel arrival time and the one of pixel 1 (hardware number: 2).
46// The relative times are filled into a histogram and an array, in order to perform
47// a Fourier analysis (see MHGausEvents). The signals are moreover averaged on an
48// event-by-event basis and written into the corresponding average pixels.
49//
50// The histograms are fitted to a Gaussian, mean and sigma with its errors
51// and the fit probability are extracted. If none of these values are NaN's and
52// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
53// the fit is declared valid.
54// Otherwise, the fit is repeated within ranges of the previous mean
55// - MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
56// In case this does not make the fit valid, the histogram means and RMS's are
57// taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
58// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kRelTimeNotFitted ) and
59// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
60//
61// Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
62// from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
63//
64// The class also fills arrays with the signal vs. event number, creates a fourier
65// spectrum (see MHGausEvents::CreateFourierSpectrum()) and investigates if the
66// projected fourier components follow an exponential distribution.
67// In case that the probability of the exponential fit is less than
68// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
69// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kRelTimeOscillating ) and
70// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
71//
72// This same procedure is performed for the average pixels.
73//
74// The following results are written into MCalibrationRelTimeCam:
75//
76// - MCalibrationPix::SetMean()
77// - MCalibrationPix::SetMeanErr()
78// - MCalibrationPix::SetSigma()
79// - MCalibrationPix::SetSigmaErr()
80// - MCalibrationPix::SetProb()
81// - MCalibrationPix::SetNumPickup()
82//
83// For all averaged areas, the fitted sigma is multiplied with the square root of
84// the number involved pixels in order to be able to compare it to the average of
85// sigmas in the camera.
86//
87/////////////////////////////////////////////////////////////////////////////
88#include "MHCalibrationRelTimeCam.h"
89#include "MHCalibrationPix.h"
90
91#include "MLog.h"
92#include "MLogManip.h"
93
94#include "MParList.h"
95
96#include "MCalibrationIntensityRelTimeCam.h"
97
98#include "MCalibrationRelTimeCam.h"
99#include "MCalibrationRelTimePix.h"
100#include "MCalibrationPix.h"
101
102#include "MArrivalTimeCam.h"
103#include "MArrivalTimePix.h"
104
105#include "MGeomCam.h"
106#include "MGeomPix.h"
107
108#include "MBadPixelsIntensityCam.h"
109#include "MBadPixelsCam.h"
110#include "MBadPixelsPix.h"
111
112#include <TOrdCollection.h>
113#include <TPad.h>
114#include <TVirtualPad.h>
115#include <TCanvas.h>
116#include <TStyle.h>
117#include <TF1.h>
118#include <TLine.h>
119#include <TLatex.h>
120#include <TLegend.h>
121#include <TGraph.h>
122#include <TEnv.h>
123
124ClassImp(MHCalibrationRelTimeCam);
125
126using namespace std;
127
128const Float_t MHCalibrationRelTimeCam::fgNumHiGainSaturationLimit = 0.25;
129const UInt_t MHCalibrationRelTimeCam::fgReferencePixel = 1;
130const Int_t MHCalibrationRelTimeCam::fgNbins = 400;
131const Axis_t MHCalibrationRelTimeCam::fgFirst = -9.975;
132const Axis_t MHCalibrationRelTimeCam::fgLast = 10.025;
133const Float_t MHCalibrationRelTimeCam::fgProbLimit = 0.0;
134const TString MHCalibrationRelTimeCam::gsHistName = "RelTime";
135const TString MHCalibrationRelTimeCam::gsHistTitle = "Arr. Times";
136const TString MHCalibrationRelTimeCam::gsHistXTitle = "Arr. Time [FADC slices]";
137const TString MHCalibrationRelTimeCam::gsHistYTitle = "Nr. events";
138const TString MHCalibrationRelTimeCam::fgReferenceFile = "mjobs/calibrationref.rc";
139
140// --------------------------------------------------------------------------
141//
142// Default Constructor.
143//
144// Sets:
145// - fReferencePixel to fgReferencePixel
146// - fNbins to fgNbins
147// - fFirst to fgFirst
148// - fLast to fgLast
149//
150// - fHistName to gsHistName
151// - fHistTitle to gsHistTitle
152// - fHistXTitle to gsHistXTitle
153// - fHistYTitle to gsHistYTitle
154//
155MHCalibrationRelTimeCam::MHCalibrationRelTimeCam(const char *name, const char *title)
156{
157
158 fName = name ? name : "MHCalibrationRelTimeCam";
159 fTitle = title ? title : "Histogram class for the relative time calibration of the camera";
160
161 SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
162
163 SetReferencePixel();
164
165 SetNbins(fgNbins);
166 SetFirst(fgFirst);
167 SetLast (fgLast );
168
169 SetProbLimit(fgProbLimit);
170
171 SetHistName (gsHistName .Data());
172 SetHistTitle (gsHistTitle .Data());
173 SetHistXTitle(gsHistXTitle.Data());
174 SetHistYTitle(gsHistYTitle.Data());
175
176 SetReferenceFile();
177
178 fInnerRefTime = 2.95;
179 fOuterRefTime = 3.6;
180}
181
182// --------------------------------------------------------------------------
183//
184// Creates new MHCalibrationRelTimeCam only with the averaged areas:
185// the rest has to be retrieved directly, e.g. via:
186// MHCalibrationRelTimeCam *cam = MParList::FindObject("MHCalibrationRelTimeCam");
187// - cam->GetAverageSector(5).DrawClone();
188// - (*cam)[100].DrawClone()
189//
190TObject *MHCalibrationRelTimeCam::Clone(const char *) const
191{
192
193 MHCalibrationRelTimeCam *cam = new MHCalibrationRelTimeCam();
194
195 //
196 // Copy the data members
197 //
198 cam->fColor = fColor;
199 cam->fRunNumbers = fRunNumbers;
200 cam->fPulserFrequency = fPulserFrequency;
201 cam->fFlags = fFlags;
202 cam->fNbins = fNbins;
203 cam->fFirst = fFirst;
204 cam->fLast = fLast;
205
206 cam->fReferenceFile = fReferenceFile;
207 cam->fInnerRefTime = fInnerRefTime;
208 cam->fOuterRefTime = fOuterRefTime;
209
210 //
211 // Copy the MArrays
212 //
213 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
214 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
215 cam->fAverageAreaSat = fAverageAreaSat;
216 cam->fAverageAreaSigma = fAverageAreaSigma;
217 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
218 cam->fAverageAreaNum = fAverageAreaNum;
219 cam->fAverageSectorNum = fAverageSectorNum;
220
221 if (!IsAverageing())
222 return cam;
223
224 const Int_t navhi = fAverageHiGainAreas->GetSize();
225
226 for (int i=0; i<navhi; i++)
227 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
228
229 if (IsLoGain())
230 {
231
232 const Int_t navlo = fAverageLoGainAreas->GetSize();
233 for (int i=0; i<navlo; i++)
234 cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
235
236 }
237
238 return cam;
239}
240
241// --------------------------------------------------------------------------
242//
243// Gets or creates the pointers to:
244// - MCalibrationRelTimeCam
245//
246// Searches pointer to:
247// - MArrivalTimeCam
248//
249// Calls:
250// - MHCalibrationCam::InitHiGainArrays()
251// - MHCalibrationCam::InitLoGainArrays()
252//
253// Sets:
254// - fSumareahi to nareas
255// - fSumarealo to nareas
256// - fSumsectorhi to nareas
257// - fSumsectorlo to nareas
258// - fNumareahi to nareas
259// - fNumarealo to nareas
260// - fNumsectorhi to nareas
261// - fNumsectorlo to nareas
262//
263Bool_t MHCalibrationRelTimeCam::ReInitHists(MParList *pList)
264{
265
266 if (!InitCams(pList,"RelTime"))
267 return kFALSE;
268
269 MArrivalTimeCam *signal = (MArrivalTimeCam*)pList->FindObject("MArrivalTimeCam");
270 if (!signal)
271 {
272 *fLog << err << "MArrivalTimeCam not found... abort." << endl;
273 return kFALSE;
274 }
275
276 const Int_t npixels = fGeom->GetNumPixels();
277 const Int_t nsectors = fGeom->GetNumSectors();
278 const Int_t nareas = fGeom->GetNumAreas();
279
280 InitHiGainArrays(npixels,nareas,nsectors);
281 InitLoGainArrays(npixels,nareas,nsectors);
282
283 fSumareahi .Set(nareas);
284 fSumarealo .Set(nareas);
285 fSumsectorhi.Set(nsectors);
286 fSumsectorlo.Set(nsectors);
287 fNumareahi .Set(nareas);
288 fNumarealo .Set(nareas);
289 fNumsectorhi.Set(nsectors);
290 fNumsectorlo.Set(nsectors);
291
292 return kTRUE;
293}
294
295
296// -------------------------------------------------------------------------------
297//
298// Retrieves pointer to MArrivalTimeCam:
299//
300// Retrieves from MGeomCam:
301// - number of pixels
302// - number of pixel areas
303// - number of sectors
304//
305// Fills HiGain or LoGain histograms (MHGausEvents::FillHistAndArray()), respectively
306// depending on MArrivalTimePix::IsLoGainUsed(), with:
307// - MArrivalTimePix::GetArrivalTime(pixid) - MArrivalTimePix::GetArrivalTime(1);
308// (i.e. the time difference between pixel i and pixel 1 (hardware number: 2) )
309//
310Bool_t MHCalibrationRelTimeCam::FillHists(const MParContainer *par, const Stat_t w)
311{
312
313 MArrivalTimeCam *arrtime = (MArrivalTimeCam*)par;
314 if (!arrtime)
315 {
316 gLog << err << "No argument in MArrivalTime::Fill... abort." << endl;
317 return kFALSE;
318 }
319
320 const Int_t npixels = fGeom->GetNumPixels();
321 const Int_t nareas = fGeom->GetNumAreas();
322 const Int_t nsectors = fGeom->GetNumSectors();
323
324 fSumareahi .Reset();
325 fSumarealo .Reset();
326 fSumsectorhi.Reset();
327 fSumsectorlo.Reset();
328 fNumareahi .Reset();
329 fNumarealo .Reset();
330 fNumsectorhi.Reset();
331 fNumsectorlo.Reset();
332
333 const MArrivalTimePix &refpix = (*arrtime)[fReferencePixel];
334 // FIXME: What do we do if pixel is invalid?
335 if (!refpix.IsArrivalTimeValid())
336 {
337 gLog << warn << "WARNING - Arrival time in refrence pixel " << fReferencePixel << " invalid." << endl;
338 return kTRUE;
339 }
340
341 const Float_t reftime = refpix.GetArrivalTime();
342// refpix.IsHiGainSaturated()
343// ? refpix.GetArrivalTimeLoGain() : refpix.GetArrivalTimeHiGain();
344
345 for (Int_t i=0; i<npixels; i++)
346 {
347 MHCalibrationPix &histhi = (*this)[i];
348 if (histhi.IsExcluded())
349 continue;
350
351 const Int_t aidx = (*fGeom)[i].GetAidx();
352 const Int_t sector = (*fGeom)[i].GetSector();
353
354 const MArrivalTimePix &pix = (*arrtime)[i];
355
356 // If hi-gain arrival time has been extracted successfully
357 // fill hi-gain histograms and arrays
358 if (pix.IsHiGainValid() && !pix.IsHiGainSaturated())
359 {
360 const Float_t time = pix.GetArrivalTimeHiGain();
361
362 if (IsOscillations())
363 histhi.FillHistAndArray(time-reftime);
364 else
365 histhi.FillHist(time-reftime);
366
367 fSumareahi [aidx] += time;
368 fNumareahi [aidx] ++;
369 fSumsectorhi[sector] += time;
370 fNumsectorhi[sector] ++;
371 }
372
373 if (!pix.IsHiGainSaturated())
374 continue;
375
376 histhi.AddSaturated(1);
377
378 // If lo-gain arrival time has been extracted successfully,
379 // the hi-gain has saturateed and the lo-gain is switched on
380 // fill hi-gain histograms and arrays
381 if (pix.IsLoGainValid() && IsLoGain())
382 {
383 const Float_t time = pix.GetArrivalTimeLoGain();
384
385 MHCalibrationPix &histlo = (*this)(i);
386 if (IsOscillations())
387 histlo.FillHistAndArray(time-reftime);
388 else
389 histlo.FillHist(time-reftime);
390
391 fSumarealo [aidx] += time;
392 fNumarealo [aidx] ++;
393 fSumsectorlo[sector] += time;
394 fNumsectorlo[sector] ++;
395 }
396 }
397
398 for (Int_t j=0; j<nareas; j++)
399 {
400 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
401 if (IsOscillations())
402 histhi.FillHistAndArray(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
403 else
404 histhi.FillHist(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
405
406 if (IsLoGain())
407 {
408 MHCalibrationPix &histlo = GetAverageLoGainArea(j);
409 if (IsOscillations())
410 histlo.FillHistAndArray(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
411 else
412 histlo.FillHist(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
413 }
414 }
415
416 for (Int_t j=0; j<nsectors; j++)
417 {
418 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
419 if (IsOscillations())
420 histhi.FillHistAndArray(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
421 else
422 histhi.FillHist(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
423
424 if (IsLoGain())
425 {
426 MHCalibrationPix &histlo = GetAverageLoGainSector(j);
427 if (IsOscillations())
428 histlo.FillHistAndArray(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
429 else
430 histlo.FillHist(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
431 }
432 }
433
434 return kTRUE;
435}
436
437// --------------------------------------------------------------------------
438//
439// Calls:
440// - MHCalibrationCam::FitHiGainArrays() with flags:
441// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
442// - MHCalibrationCam::FitLoGainArrays() with flags:
443// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
444//
445Bool_t MHCalibrationRelTimeCam::FinalizeHists()
446{
447
448 *fLog << endl;
449
450 MCalibrationCam *relcam = fIntensCam ? fIntensCam->GetCam() : fCam;
451 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
452
453 const Int_t nareas = fAverageHiGainAreas->GetSize();
454 const Int_t nsectors = fAverageHiGainSectors->GetSize();
455
456 TArrayI satarea(nareas);
457 TArrayI satsect(nsectors);
458 fNumareahi .Reset();
459 fNumsectorhi.Reset();
460
461 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
462 {
463
464 MHCalibrationPix &histhi = (*this)[i];
465
466 if (histhi.IsExcluded())
467 continue;
468
469 const Int_t aidx = (*fGeom)[i].GetAidx();
470 const Int_t sector = (*fGeom)[i].GetSector();
471
472 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)(*relcam)[i] ;
473
474 fNumareahi[aidx]++;
475 fNumsectorhi[sector]++;
476 //
477 // Check saturation
478 //
479 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
480 {
481 pix.SetHiGainSaturation();
482 histhi.SetExcluded();
483 satarea[aidx]++;
484 satsect[sector]++;
485 }
486 else
487 if (IsLoGain())
488 (*this)(i).SetExcluded();
489
490 //
491 // Check histogram overflow
492 //
493 CheckOverflow(histhi);
494 if (IsLoGain())
495 CheckOverflow((*this)(i));
496
497 }
498
499 for (Int_t j=0; j<nareas; j++)
500 {
501
502 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
503 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageArea(j);
504
505 if (satarea[j] > 0.5*fNumareahi[j])
506 {
507 pix.SetHiGainSaturation();
508 histhi.SetExcluded();
509 }
510 else
511 if (IsLoGain())
512 GetAverageLoGainArea(j).SetExcluded();
513
514 //
515 // Check histogram overflow
516 //
517 CheckOverflow(histhi);
518 if (IsLoGain())
519 CheckOverflow(GetAverageLoGainArea(j));
520 }
521
522 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
523 {
524
525 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
526 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageSector(j) ;
527
528 if (satsect[j] > 0.5*fNumsectorhi[j])
529 {
530 pix.SetHiGainSaturation();
531 histhi.SetExcluded();
532 }
533 else
534 if (IsLoGain())
535 GetAverageLoGainSector(j).SetExcluded();
536
537 //
538 // Check histogram overflow
539 //
540 CheckOverflow(histhi);
541 if (IsLoGain())
542 CheckOverflow(GetAverageLoGainSector(j));
543 }
544
545 FitHiGainArrays(*relcam,*badcam,
546 MBadPixelsPix::kRelTimeNotFitted,
547 MBadPixelsPix::kRelTimeOscillating);
548
549 if (IsLoGain())
550 FitLoGainArrays(*relcam,*badcam,
551 MBadPixelsPix::kRelTimeNotFitted,
552 MBadPixelsPix::kRelTimeOscillating);
553
554 return kTRUE;
555}
556
557// --------------------------------------------------------------------------
558//
559// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
560// - MBadPixelsPix::kRelTimeNotFitted
561// - MBadPixelsPix::kRelTimeOscillating
562//
563void MHCalibrationRelTimeCam::FinalizeBadPixels()
564{
565
566 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
567
568 for (Int_t i=0; i<badcam->GetSize(); i++)
569 {
570 MBadPixelsPix &bad = (*badcam)[i];
571
572 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeNotFitted ))
573 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
574
575 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeOscillating))
576 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
577
578 }
579}
580
581// --------------------------------------------------------------------------
582//
583// The types are as follows:
584//
585// Fitted values:
586// ==============
587//
588// 0: Fitted Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMean()
589// 1: Error Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMeanErr()
590// 2: Sigma fitted Relative Arrival Time in FADC slices (MHGausEvents::GetSigma()
591// 3: Error Sigma Relative Arrival Time in FADC slices (MHGausEvents::GetSigmaErr()
592//
593// Useful variables derived from the fit results:
594// =============================================
595//
596// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
597//
598// Localized defects:
599// ==================
600//
601// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
602// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
603//
604Bool_t MHCalibrationRelTimeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
605{
606
607 if (fHiGainArray->GetSize() <= idx)
608 return kFALSE;
609
610 const MHCalibrationPix &pix = (*this)[idx];
611
612 switch (type)
613 {
614 case 0:
615 val = pix.GetMean();
616 break;
617 case 1:
618 val = pix.GetMeanErr();
619 break;
620 case 2:
621 val = pix.GetSigma();
622 break;
623 case 3:
624 val = pix.GetSigmaErr();
625 break;
626 case 4:
627 val = pix.GetProb();
628 break;
629 case 5:
630 if (!pix.IsGausFitOK())
631 val = 1.;
632 break;
633 case 6:
634 if (!pix.IsFourierSpectrumOK())
635 val = 1.;
636 break;
637 default:
638 return kFALSE;
639 }
640 return kTRUE;
641}
642
643// --------------------------------------------------------------------------
644//
645// Calls MHCalibrationPix::DrawClone() for pixel idx
646//
647void MHCalibrationRelTimeCam::DrawPixelContent(Int_t idx) const
648{
649 (*this)[idx].DrawClone();
650}
651
652// -----------------------------------------------------------------------------
653//
654// Default draw:
655//
656// Displays the averaged areas, both High Gain and Low Gain
657//
658// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
659//
660void MHCalibrationRelTimeCam::Draw(const Option_t *opt)
661{
662
663 const Int_t nareas = fAverageHiGainAreas->GetSize();
664 if (nareas == 0)
665 return;
666
667 TString option(opt);
668 option.ToLower();
669
670 if (!option.Contains("datacheck"))
671 {
672 MHCalibrationCam::Draw(opt);
673 return;
674 }
675
676 //
677 // From here on , the datacheck - Draw
678 //
679 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
680 pad->SetBorderMode(0);
681 pad->Divide(1,nareas);
682
683 //
684 // Loop over inner and outer pixels
685 //
686 for (Int_t i=0; i<nareas;i++)
687 {
688
689 pad->cd(i+1);
690
691 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
692 //
693 // Ask for Hi-Gain saturation
694 //
695 if (hipix.IsExcluded() && IsLoGain())
696 {
697 MHCalibrationPix &lopix = GetAverageLoGainArea(i);
698 DrawDataCheckPixel(lopix,i ? fOuterRefTime+1.5 : fInnerRefTime+1.5);
699 }
700 else
701 DrawDataCheckPixel(hipix,i ? fOuterRefTime : fInnerRefTime);
702 }
703}
704
705void MHCalibrationRelTimeCam::CheckOverflow( MHCalibrationPix &pix ) const
706{
707 if (pix.IsExcluded())
708 return;
709
710 const TH1F &hist = *pix.GetHGausHist();
711
712 const Int_t n = hist.GetNbinsX();
713 const Float_t max = fOverflowLimit*hist.GetEntries();
714
715 const Stat_t overflow = hist.GetBinContent(n+1);
716 if (overflow > max)
717 {
718 *fLog << warn << overflow << " overflows above " << hist.GetBinLowEdge(n);
719 *fLog << " in " << pix.GetName() << " (w/o saturation!) " << endl;
720 }
721
722 const Stat_t underflow = hist.GetBinContent(0);
723 if (underflow > max)
724 {
725 *fLog << warn << underflow << " underflows below " << hist.GetBinLowEdge(1);
726 *fLog << " in " << pix.GetName() << " (w/o saturation!) " << endl;
727 }
728}
729
730
731// -----------------------------------------------------------------------------
732//
733// Draw the average pixel for the datacheck:
734//
735// Displays the averaged areas, both High Gain and Low Gain
736//
737// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
738//
739void MHCalibrationRelTimeCam::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
740{
741
742 if (pix.IsEmpty())
743 return;
744
745 TVirtualPad *newpad = gPad;
746 newpad->Divide(1,2);
747 newpad->cd(1);
748
749 gPad->SetTicks();
750 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
751 gPad->SetLogy();
752
753 TH1F *hist = pix.GetHGausHist();
754
755 TH1F *null = new TH1F("Null",hist->GetTitle(),100,0.,pix.GetLast());
756
757 null->SetMaximum(1.1*hist->GetMaximum());
758 null->SetDirectory(NULL);
759 null->SetBit(kCanDelete);
760 null->SetStats(kFALSE);
761 //
762 // set the labels bigger
763 //
764 TAxis *xaxe = null->GetXaxis();
765 TAxis *yaxe = null->GetYaxis();
766 xaxe->CenterTitle();
767 yaxe->CenterTitle();
768 xaxe->SetTitleSize(0.07);
769 yaxe->SetTitleSize(0.07);
770 xaxe->SetTitleOffset(0.65);
771 yaxe->SetTitleOffset(0.55);
772 xaxe->SetLabelSize(0.06);
773 yaxe->SetLabelSize(0.06);
774 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
775 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
776
777 null->Draw();
778 hist->Draw("same");
779
780 gStyle->SetOptFit();
781
782 TF1 *fit = pix.GetFGausFit();
783
784 if (fit)
785 {
786 switch ( fColor )
787 {
788 case MCalibrationCam::kGREEN:
789 fit->SetLineColor(kGreen);
790 break;
791 case MCalibrationCam::kBLUE:
792 fit->SetLineColor(kBlue);
793 break;
794 case MCalibrationCam::kUV:
795 fit->SetLineColor(106);
796 break;
797 case MCalibrationCam::kCT1:
798 fit->SetLineColor(006);
799 break;
800 default:
801 fit->SetLineColor(kRed);
802 }
803 fit->Draw("same");
804 }
805
806 DisplayRefLines(null,refline);
807
808 newpad->cd(2);
809 gPad->SetTicks();
810
811 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
812
813 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
814 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
815 null2->SetDirectory(NULL);
816 null2->SetBit(kCanDelete);
817 null2->SetStats(kFALSE);
818 //
819 // set the labels bigger
820 //
821 TAxis *xaxe2 = null2->GetXaxis();
822 TAxis *yaxe2 = null2->GetYaxis();
823 xaxe2->CenterTitle();
824 yaxe2->CenterTitle();
825 xaxe2->SetTitleSize(0.07);
826 yaxe2->SetTitleSize(0.07);
827 xaxe2->SetTitleOffset(0.65);
828 yaxe2->SetTitleOffset(0.55);
829 xaxe2->SetLabelSize(0.06);
830 yaxe2->SetLabelSize(0.06);
831
832 pix.CreateGraphEvents();
833 TGraph *gr = pix.GetGraphEvents();
834 if (gr)
835 {
836 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
837 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
838 }
839
840 null2->Draw();
841
842 pix.DrawEvents("same");
843
844 return;
845}
846
847void MHCalibrationRelTimeCam::DisplayRefLines(const TH1F *hist, const Float_t refline) const
848{
849
850 TGraph *gr = new TGraph(2);
851 gr->SetPoint(0,refline,0.);
852 gr->SetPoint(1,refline,hist->GetMaximum());
853 gr->SetBit(kCanDelete);
854 gr->SetLineColor(kGreen);
855 gr->SetLineStyle(2);
856 gr->SetLineWidth(3);
857 gr->Draw("L");
858
859 TLegend *leg = new TLegend(0.75,0.7,0.99,0.99);
860 leg->SetBit(kCanDelete);
861 leg->AddEntry(gr,"Trigger Calibration","l");
862 leg->Draw();
863}
864
865Int_t MHCalibrationRelTimeCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
866{
867 Bool_t rc = kFALSE;
868
869 if (IsEnvDefined(env, prefix, "ReferenceFile", print))
870 {
871 SetReferenceFile(GetEnvValue(env, prefix, "ReferenceFile", fReferenceFile.Data()));
872 rc = kTRUE;
873 }
874
875 if (IsEnvDefined(env, prefix, "ReferencePixel", print))
876 {
877 SetReferencePixel(GetEnvValue(env, prefix, "ReferencePixel", (Int_t)fReferencePixel));
878 rc = kTRUE;
879 }
880
881 TEnv refenv(fReferenceFile);
882
883 fInnerRefTime = refenv.GetValue("InnerRefTime", fInnerRefTime);
884 fOuterRefTime = refenv.GetValue("OuterRefTime", fOuterRefTime);
885
886 return MHCalibrationCam::ReadEnv(env,prefix,print) ? kTRUE : rc;
887}
Note: See TracBrowser for help on using the repository browser.