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

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