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

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