source: trunk/Mars/mhcalib/MHCalibrationRelTimeCam.cc@ 10083

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