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

Last change on this file since 5910 was 5864, checked in by gaug, 20 years ago
*** empty log message ***
File size: 23.1 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 if (!InitCams(pList,"RelTime"))
255 return kFALSE;
256
257 MArrivalTimeCam *signal = (MArrivalTimeCam*)pList->FindObject("MArrivalTimeCam");
258 if (!signal)
259 {
260 *fLog << err << "MArrivalTimeCam not found... abort." << endl;
261 return kFALSE;
262 }
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 const MArrivalTimePix &refpix = (*arrtime)[fReferencePixel];
322 const Float_t reftime = refpix.IsLoGainUsed()
323 ? refpix.GetArrivalTimeLoGain() : refpix.GetArrivalTimeHiGain();
324
325 for (Int_t i=0; i<npixels; i++)
326 {
327
328 MHCalibrationPix &histhi = (*this)[i];
329
330 if (histhi.IsExcluded())
331 continue;
332
333 const MArrivalTimePix &pix = (*arrtime)[i];
334 const Int_t aidx = (*fGeom)[i].GetAidx();
335 const Int_t sector = (*fGeom)[i].GetSector();
336
337 if (pix.IsLoGainUsed() && IsLoGain())
338 {
339 const Float_t reltime = pix.GetArrivalTimeLoGain() - reftime;
340 histhi.AddSaturated(1);
341
342 MHCalibrationPix &histlo = (*this)(i);
343 histlo.FillHistAndArray(reltime);
344
345 fSumarealo [aidx] += reltime;
346 fNumarealo [aidx] ++;
347 fSumsectorlo[sector] += reltime;
348 fNumsectorlo[sector] ++;
349 }
350 else
351 {
352 const Float_t reltime = pix.GetArrivalTimeHiGain() - reftime;
353
354 histhi.FillHistAndArray(reltime) ;
355 fSumareahi [aidx] += reltime;
356 fNumareahi [aidx] ++;
357 fSumsectorhi[sector] += reltime;
358 fNumsectorhi[sector] ++;
359 }
360 }
361
362 for (Int_t j=0; j<nareas; j++)
363 {
364 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
365 histhi.FillHistAndArray(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
366
367 if (IsLoGain())
368 {
369 MHCalibrationPix &histlo = GetAverageLoGainArea(j);
370 histlo.FillHistAndArray(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
371 }
372 }
373
374 for (Int_t j=0; j<nsectors; j++)
375 {
376 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
377 histhi.FillHistAndArray(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
378
379 if (IsLoGain())
380 {
381 MHCalibrationPix &histlo = GetAverageLoGainSector(j);
382 histlo.FillHistAndArray(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
383 }
384 }
385
386 return kTRUE;
387}
388
389// --------------------------------------------------------------------------
390//
391// Calls:
392// - MHCalibrationCam::FitHiGainArrays() with flags:
393// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
394// - MHCalibrationCam::FitLoGainArrays() with flags:
395// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
396//
397Bool_t MHCalibrationRelTimeCam::FinalizeHists()
398{
399
400 *fLog << endl;
401
402 MCalibrationCam *relcam = fIntensCam ? fIntensCam->GetCam() : fCam;
403 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
404
405 const Int_t nareas = fAverageHiGainAreas->GetSize();
406 const Int_t nsectors = fAverageHiGainSectors->GetSize();
407
408 TArrayI satarea(nareas);
409 TArrayI satsect(nsectors);
410 fNumareahi .Reset();
411 fNumsectorhi.Reset();
412
413 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
414 {
415
416 MHCalibrationPix &histhi = (*this)[i];
417
418 if (histhi.IsExcluded())
419 continue;
420
421 const Int_t aidx = (*fGeom)[i].GetAidx();
422 const Int_t sector = (*fGeom)[i].GetSector();
423
424 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)(*relcam)[i] ;
425
426 fNumareahi[aidx]++;
427 fNumsectorhi[sector]++;
428 //
429 // Check saturation
430 //
431 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
432 {
433 pix.SetHiGainSaturation();
434 histhi.SetExcluded();
435 satarea[aidx]++;
436 satsect[sector]++;
437 }
438 else
439 if (IsLoGain())
440 (*this)(i).SetExcluded();
441
442 //
443 // Check histogram overflow
444 //
445 CheckOverflow(histhi);
446 if (IsLoGain())
447 CheckOverflow((*this)(i));
448
449 }
450
451 for (Int_t j=0; j<nareas; j++)
452 {
453
454 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
455 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageArea(j);
456
457 if (satarea[j] > 0.5*fNumareahi[j])
458 {
459 pix.SetHiGainSaturation();
460 histhi.SetExcluded();
461 }
462 else
463 if (IsLoGain())
464 GetAverageLoGainArea(j).SetExcluded();
465
466 //
467 // Check histogram overflow
468 //
469 CheckOverflow(histhi);
470 if (IsLoGain())
471 CheckOverflow(GetAverageLoGainArea(j));
472 }
473
474 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
475 {
476
477 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
478 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageSector(j) ;
479
480 if (satsect[j] > 0.5*fNumsectorhi[j])
481 {
482 pix.SetHiGainSaturation();
483 histhi.SetExcluded();
484 }
485 else
486 if (IsLoGain())
487 GetAverageLoGainSector(j).SetExcluded();
488
489 //
490 // Check histogram overflow
491 //
492 CheckOverflow(histhi);
493 if (IsLoGain())
494 CheckOverflow(GetAverageLoGainSector(j));
495 }
496
497 FitHiGainArrays(*relcam,*badcam,
498 MBadPixelsPix::kRelTimeNotFitted,
499 MBadPixelsPix::kRelTimeOscillating);
500
501 if (IsLoGain())
502 FitLoGainArrays(*relcam,*badcam,
503 MBadPixelsPix::kRelTimeNotFitted,
504 MBadPixelsPix::kRelTimeOscillating);
505
506 return kTRUE;
507}
508
509// --------------------------------------------------------------------------
510//
511// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
512// - MBadPixelsPix::kRelTimeNotFitted
513// - MBadPixelsPix::kRelTimeOscillating
514//
515void MHCalibrationRelTimeCam::FinalizeBadPixels()
516{
517
518 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
519
520 for (Int_t i=0; i<badcam->GetSize(); i++)
521 {
522 MBadPixelsPix &bad = (*badcam)[i];
523
524 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeNotFitted ))
525 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
526
527 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeOscillating))
528 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
529
530 }
531}
532
533// --------------------------------------------------------------------------
534//
535// The types are as follows:
536//
537// Fitted values:
538// ==============
539//
540// 0: Fitted Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMean()
541// 1: Error Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMeanErr()
542// 2: Sigma fitted Relative Arrival Time in FADC slices (MHGausEvents::GetSigma()
543// 3: Error Sigma Relative Arrival Time in FADC slices (MHGausEvents::GetSigmaErr()
544//
545// Useful variables derived from the fit results:
546// =============================================
547//
548// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
549//
550// Localized defects:
551// ==================
552//
553// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
554// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
555//
556Bool_t MHCalibrationRelTimeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
557{
558
559 if (fHiGainArray->GetSize() <= idx)
560 return kFALSE;
561
562 const MHCalibrationPix &pix = (*this)[idx];
563
564 switch (type)
565 {
566 case 0:
567 val = pix.GetMean();
568 break;
569 case 1:
570 val = pix.GetMeanErr();
571 break;
572 case 2:
573 val = pix.GetSigma();
574 break;
575 case 3:
576 val = pix.GetSigmaErr();
577 break;
578 case 4:
579 val = pix.GetProb();
580 break;
581 case 5:
582 if (!pix.IsGausFitOK())
583 val = 1.;
584 break;
585 case 6:
586 if (!pix.IsFourierSpectrumOK())
587 val = 1.;
588 break;
589 default:
590 return kFALSE;
591 }
592 return kTRUE;
593}
594
595// --------------------------------------------------------------------------
596//
597// Calls MHCalibrationPix::DrawClone() for pixel idx
598//
599void MHCalibrationRelTimeCam::DrawPixelContent(Int_t idx) const
600{
601 (*this)[idx].DrawClone();
602}
603
604// -----------------------------------------------------------------------------
605//
606// Default draw:
607//
608// Displays the averaged areas, both High Gain and Low Gain
609//
610// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
611//
612void MHCalibrationRelTimeCam::Draw(const Option_t *opt)
613{
614
615 const Int_t nareas = fAverageHiGainAreas->GetSize();
616 if (nareas == 0)
617 return;
618
619 TString option(opt);
620 option.ToLower();
621
622 if (!option.Contains("datacheck"))
623 {
624 MHCalibrationCam::Draw(opt);
625 return;
626 }
627
628 //
629 // From here on , the datacheck - Draw
630 //
631 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
632 pad->SetBorderMode(0);
633 pad->Divide(1,nareas);
634
635 //
636 // Loop over inner and outer pixels
637 //
638 for (Int_t i=0; i<nareas;i++)
639 {
640
641 pad->cd(i+1);
642
643 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
644 //
645 // Ask for Hi-Gain saturation
646 //
647 if (hipix.IsExcluded() && IsLoGain())
648 {
649 MHCalibrationPix &lopix = GetAverageLoGainArea(i);
650 DrawDataCheckPixel(lopix,0.);
651 }
652 else
653 DrawDataCheckPixel(hipix,0.);
654 }
655}
656
657void MHCalibrationRelTimeCam::CheckOverflow( MHCalibrationPix &pix )
658{
659
660 if (pix.IsExcluded())
661 return;
662
663 TH1F *hist = pix.GetHGausHist();
664
665 Stat_t overflow = hist->GetBinContent(hist->GetNbinsX()+1);
666 if (overflow > fOverflowLimit*hist->GetEntries())
667 {
668 *fLog << warn << "HiGain Hist-overflow " << overflow
669 << " times in " << pix.GetName() << " (w/o saturation!) " << endl;
670 }
671
672 overflow = hist->GetBinContent(0);
673 if (overflow > fOverflowLimit*hist->GetEntries())
674 {
675 *fLog << warn << "HiGain Hist-underflow " << overflow
676 << " times in " << pix.GetName() << " (w/o saturation!) " << endl;
677 }
678}
679
680
681// -----------------------------------------------------------------------------
682//
683// Draw the average pixel for the datacheck:
684//
685// Displays the averaged areas, both High Gain and Low Gain
686//
687// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
688//
689void MHCalibrationRelTimeCam::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
690{
691
692 TVirtualPad *newpad = gPad;
693 newpad->Divide(1,2);
694 newpad->cd(1);
695
696 gPad->SetTicks();
697 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
698 gPad->SetLogy();
699
700 TH1F *hist = pix.GetHGausHist();
701
702 TH1F *null = new TH1F("Null",hist->GetTitle(),100,pix.GetFirst() < -3.0 ? -3.0 : pix.GetFirst(),
703 pix.GetLast() > 3.0 ? 3.0 : pix.GetLast());
704
705 null->SetMaximum(1.1*hist->GetMaximum());
706 null->SetDirectory(NULL);
707 null->SetBit(kCanDelete);
708 null->SetStats(kFALSE);
709 //
710 // set the labels bigger
711 //
712 TAxis *xaxe = null->GetXaxis();
713 TAxis *yaxe = null->GetYaxis();
714 xaxe->CenterTitle();
715 yaxe->CenterTitle();
716 xaxe->SetTitleSize(0.07);
717 yaxe->SetTitleSize(0.07);
718 xaxe->SetTitleOffset(0.7);
719 yaxe->SetTitleOffset(0.55);
720 xaxe->SetLabelSize(0.06);
721 yaxe->SetLabelSize(0.06);
722
723 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
724 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
725
726 null->Draw();
727 hist->Draw("same");
728
729 gStyle->SetOptFit();
730
731 TF1 *fit = pix.GetFGausFit();
732
733 if (fit)
734 {
735 switch ( fColor )
736 {
737 case MCalibrationCam::kGREEN:
738 fit->SetLineColor(kGreen);
739 break;
740 case MCalibrationCam::kBLUE:
741 fit->SetLineColor(kBlue);
742 break;
743 case MCalibrationCam::kUV:
744 fit->SetLineColor(106);
745 break;
746 case MCalibrationCam::kCT1:
747 fit->SetLineColor(006);
748 break;
749 default:
750 fit->SetLineColor(kRed);
751 }
752 fit->Draw("same");
753 }
754
755 // DisplayRefLines(null,refline);
756
757 gPad->Modified();
758 gPad->Update();
759
760 newpad->cd(2);
761 gPad->SetTicks();
762
763 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
764
765 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
766 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
767 null2->SetDirectory(NULL);
768 null2->SetBit(kCanDelete);
769 null2->SetStats(kFALSE);
770 //
771 // set the labels bigger
772 //
773 TAxis *xaxe2 = null2->GetXaxis();
774 TAxis *yaxe2 = null2->GetYaxis();
775 xaxe2->CenterTitle();
776 yaxe2->CenterTitle();
777 xaxe2->SetTitleSize(0.07);
778 yaxe2->SetTitleSize(0.07);
779 xaxe2->SetTitleOffset(0.7);
780 yaxe2->SetTitleOffset(0.55);
781 xaxe2->SetLabelSize(0.06);
782 yaxe2->SetLabelSize(0.06);
783
784 pix.CreateGraphEvents();
785 TGraph *gr = pix.GetGraphEvents();
786
787 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
788 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
789
790 null2->Draw();
791
792 pix.DrawEvents("same");
793
794 return;
795}
796
Note: See TracBrowser for help on using the repository browser.