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

Last change on this file since 5138 was 5137, checked in by gaug, 20 years ago
*** empty log message ***
File size: 23.0 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 = 900;
129const Axis_t MHCalibrationRelTimeCam::fgFirst = -5.;
130const Axis_t MHCalibrationRelTimeCam::fgLast = 5.;
131const TString MHCalibrationRelTimeCam::gsHistName = "RelTime";
132const TString MHCalibrationRelTimeCam::gsHistTitle = "Rel. Arr. Times";
133const TString MHCalibrationRelTimeCam::gsHistXTitle = "Rel. Arr. Time [FADC slices]";
134const TString MHCalibrationRelTimeCam::gsHistYTitle = "Nr. events";
135// --------------------------------------------------------------------------
136//
137// Default Constructor.
138//
139// Sets:
140// - fReferencePixel to fgReferencePixel
141// - fNbins to fgNbins
142// - fFirst to fgFirst
143// - fLast to fgLast
144//
145// - fHistName to gsHistName
146// - fHistTitle to gsHistTitle
147// - fHistXTitle to gsHistXTitle
148// - fHistYTitle to gsHistYTitle
149//
150MHCalibrationRelTimeCam::MHCalibrationRelTimeCam(const char *name, const char *title)
151{
152
153 fName = name ? name : "MHCalibrationRelTimeCam";
154 fTitle = title ? title : "Histogram class for the relative time calibration of the camera";
155
156 SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
157
158 SetReferencePixel();
159
160 SetNbins(fgNbins);
161 SetFirst(fgFirst);
162 SetLast (fgLast );
163
164 SetHistName (gsHistName .Data());
165 SetHistTitle (gsHistTitle .Data());
166 SetHistXTitle(gsHistXTitle.Data());
167 SetHistYTitle(gsHistYTitle.Data());
168
169}
170
171// --------------------------------------------------------------------------
172//
173// Creates new MHCalibrationRelTimeCam only with the averaged areas:
174// the rest has to be retrieved directly, e.g. via:
175// MHCalibrationRelTimeCam *cam = MParList::FindObject("MHCalibrationRelTimeCam");
176// - cam->GetAverageSector(5).DrawClone();
177// - (*cam)[100].DrawClone()
178//
179#if 0
180TObject *MHCalibrationRelTimeCam::Clone(const char *) const
181{
182
183 MHCalibrationRelTimeCam *cam = new MHCalibrationRelTimeCam();
184
185 //
186 // Copy the data members
187 //
188 cam->fRunNumbers = fRunNumbers;
189 cam->fPulserFrequency = fPulserFrequency;
190 cam->fFlags = fFlags;
191 cam->fNbins = fNbins;
192 cam->fFirst = fFirst;
193 cam->fLast = fLast;
194
195 //
196 // Copy the MArrays
197 //
198 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
199 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
200 cam->fAverageAreaSat = fAverageAreaSat;
201 cam->fAverageAreaSigma = fAverageAreaSigma;
202 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
203 cam->fAverageAreaNum = fAverageAreaNum;
204 cam->fAverageSectorNum = fAverageSectorNum;
205
206 if (!IsAverageing())
207 return cam;
208
209 const Int_t navhi = fAverageHiGainAreas->GetSize();
210
211 for (int i=0; i<navhi; i++)
212 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
213
214 if (IsLoGain())
215 {
216
217 const Int_t navlo = fAverageLoGainAreas->GetSize();
218 for (int i=0; i<navlo; i++)
219 cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
220
221 }
222
223 return cam;
224}
225#endif
226
227// --------------------------------------------------------------------------
228//
229// Gets or creates the pointers to:
230// - MCalibrationRelTimeCam
231//
232// Searches pointer to:
233// - MArrivalTimeCam
234//
235// Calls:
236// - MHCalibrationCam::InitHiGainArrays()
237// - MHCalibrationCam::InitLoGainArrays()
238//
239// Sets:
240// - fSumareahi to nareas
241// - fSumarealo to nareas
242// - fSumsectorhi to nareas
243// - fSumsectorlo to nareas
244// - fNumareahi to nareas
245// - fNumarealo to nareas
246// - fNumsectorhi to nareas
247// - fNumsectorlo to nareas
248//
249Bool_t MHCalibrationRelTimeCam::ReInitHists(MParList *pList)
250{
251
252 fIntensCam = (MCalibrationIntensityCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityRelTimeCam"));
253 if (fIntensCam)
254 *fLog << inf << "Found MCalibrationIntensityRelTimeCam ... " << endl;
255 else
256 {
257 fCam = (MCalibrationCam*)pList->FindObject(AddSerialNumber("MCalibrationRelTimeCam"));
258 if (!fCam)
259 {
260 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationRelTimeCam"));
261 if (!fCam)
262 {
263 *fLog << err << "Cannot find nor create MCalibrationRelTimeCam ... abort." << endl;
264 return kFALSE;
265 }
266 fCam->Init(*fGeom);
267 }
268 }
269
270 MArrivalTimeCam *signal = (MArrivalTimeCam*)pList->FindObject("MArrivalTimeCam");
271 if (!signal)
272 {
273 *fLog << err << "MArrivalTimeCam not found... abort." << endl;
274 return kFALSE;
275 }
276
277 const Int_t npixels = fGeom->GetNumPixels();
278 const Int_t nsectors = fGeom->GetNumSectors();
279 const Int_t nareas = fGeom->GetNumAreas();
280
281 InitHiGainArrays(npixels,nareas,nsectors);
282 InitLoGainArrays(npixels,nareas,nsectors);
283
284 fSumareahi .Set(nareas);
285 fSumarealo .Set(nareas);
286 fSumsectorhi.Set(nsectors);
287 fSumsectorlo.Set(nsectors);
288 fNumareahi .Set(nareas);
289 fNumarealo .Set(nareas);
290 fNumsectorhi.Set(nsectors);
291 fNumsectorlo.Set(nsectors);
292
293 return kTRUE;
294}
295
296
297// -------------------------------------------------------------------------------
298//
299// Retrieves pointer to MArrivalTimeCam:
300//
301// Retrieves from MGeomCam:
302// - number of pixels
303// - number of pixel areas
304// - number of sectors
305//
306// Fills HiGain or LoGain histograms (MHGausEvents::FillHistAndArray()), respectively
307// depending on MArrivalTimePix::IsLoGainUsed(), with:
308// - MArrivalTimePix::GetArrivalTime(pixid) - MArrivalTimePix::GetArrivalTime(1);
309// (i.e. the time difference between pixel i and pixel 1 (hardware number: 2) )
310//
311Bool_t MHCalibrationRelTimeCam::FillHists(const MParContainer *par, const Stat_t w)
312{
313
314 MArrivalTimeCam *arrtime = (MArrivalTimeCam*)par;
315 if (!arrtime)
316 {
317 gLog << err << "No argument in MArrivalTime::Fill... abort." << endl;
318 return kFALSE;
319 }
320
321 const Int_t npixels = fGeom->GetNumPixels();
322 const Int_t nareas = fGeom->GetNumAreas();
323 const Int_t nsectors = fGeom->GetNumSectors();
324
325 fSumareahi .Reset();
326 fSumarealo .Reset();
327 fSumsectorhi.Reset();
328 fSumsectorlo.Reset();
329 fNumareahi .Reset();
330 fNumarealo .Reset();
331 fNumsectorhi.Reset();
332 fNumsectorlo.Reset();
333
334 const MArrivalTimePix &refpix = (*arrtime)[fReferencePixel];
335 const Float_t reftime = refpix.IsLoGainUsed()
336 ? refpix.GetArrivalTimeLoGain() : refpix.GetArrivalTimeHiGain();
337
338 for (Int_t i=0; i<npixels; i++)
339 {
340
341 MHCalibrationPix &histhi = (*this)[i];
342 MHCalibrationPix &histlo = (*this)(i);
343
344 if (histhi.IsExcluded())
345 continue;
346
347 const MArrivalTimePix &pix = (*arrtime)[i];
348 const Int_t aidx = (*fGeom)[i].GetAidx();
349 const Int_t sector = (*fGeom)[i].GetSector();
350
351 if (pix.IsLoGainUsed() && IsLoGain())
352 {
353 const Float_t reltime = pix.GetArrivalTimeLoGain() - reftime;
354 histhi.AddSaturated(1);
355 histlo.FillHistAndArray(reltime);
356 fSumarealo [aidx] += reltime;
357 fNumarealo [aidx] ++;
358 fSumsectorlo[sector] += reltime;
359 fNumsectorlo[sector] ++;
360 }
361 else
362 {
363 const Float_t reltime = pix.GetArrivalTimeHiGain() - reftime;
364
365 histhi.FillHistAndArray(reltime) ;
366 fSumareahi [aidx] += reltime;
367 fNumareahi [aidx] ++;
368 fSumsectorhi[sector] += reltime;
369 fNumsectorhi[sector] ++;
370 }
371 }
372
373 for (Int_t j=0; j<nareas; j++)
374 {
375 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
376 histhi.FillHistAndArray(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
377
378 if (IsLoGain())
379 {
380 MHCalibrationPix &histlo = GetAverageLoGainArea(j);
381 histlo.FillHistAndArray(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
382 }
383 }
384
385 for (Int_t j=0; j<nsectors; j++)
386 {
387 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
388 histhi.FillHistAndArray(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
389
390 if (IsLoGain())
391 {
392 MHCalibrationPix &histlo = GetAverageLoGainSector(j);
393 histlo.FillHistAndArray(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
394 }
395 }
396
397 return kTRUE;
398}
399
400// --------------------------------------------------------------------------
401//
402// Calls:
403// - MHCalibrationCam::FitHiGainArrays() with flags:
404// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
405// - MHCalibrationCam::FitLoGainArrays() with flags:
406// MBadPixelsPix::kRelTimeNotFitted and MBadPixelsPix::kRelTimeOscillating
407//
408Bool_t MHCalibrationRelTimeCam::FinalizeHists()
409{
410
411 *fLog << endl;
412
413 MCalibrationCam *relcam = fIntensCam ? fIntensCam->GetCam() : fCam;
414 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
415
416 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
417 {
418
419 MHCalibrationPix &histhi = (*this)[i];
420
421 if (histhi.IsExcluded())
422 continue;
423
424 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)(*relcam)[i] ;
425
426 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
427 {
428 pix.SetHiGainSaturation();
429 histhi.SetExcluded();
430 }
431 else
432 if (IsLoGain())
433 (*this)(i).SetExcluded();
434
435 Stat_t overflow = histhi.GetHGausHist()->GetBinContent(histhi.GetHGausHist()->GetNbinsX()+1);
436 if (overflow > 0.1)
437 {
438 *fLog << warn << GetDescriptor()
439 << ": HiGain Histogram Overflow occurred " << overflow
440 << " times in pixel: " << i << " (without saturation!) " << endl;
441 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
442 }
443
444 overflow = histhi.GetHGausHist()->GetBinContent(0);
445 if (overflow > 0.1)
446 {
447 *fLog << warn << GetDescriptor()
448 << ": HiGain Histogram Underflow occurred " << overflow
449 << " times in pixel: " << i << " (without saturation!) " << endl;
450 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
451 }
452 }
453
454 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
455 {
456
457 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
458
459 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
460 {
461 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageArea(j);
462 pix.SetHiGainSaturation();
463 histhi.SetExcluded();
464 }
465 else
466 if (IsLoGain())
467 GetAverageLoGainArea(j).SetExcluded();
468
469 }
470
471 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
472 {
473
474 MHCalibrationPix &histhi = GetAverageHiGainSector(j);
475 MCalibrationRelTimePix &pix = (MCalibrationRelTimePix&)relcam->GetAverageSector(j) ;
476
477 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
478 {
479 pix.SetHiGainSaturation();
480 histhi.SetExcluded();
481 }
482 else
483 if (IsLoGain())
484 GetAverageLoGainSector(j).SetExcluded();
485 }
486
487 FitHiGainArrays(*relcam,*badcam,
488 MBadPixelsPix::kRelTimeNotFitted,
489 MBadPixelsPix::kRelTimeOscillating);
490
491 if (IsLoGain())
492 FitLoGainArrays(*relcam,*badcam,
493 MBadPixelsPix::kRelTimeNotFitted,
494 MBadPixelsPix::kRelTimeOscillating);
495
496 return kTRUE;
497}
498
499// --------------------------------------------------------------------------
500//
501// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
502// - MBadPixelsPix::kRelTimeNotFitted
503// - MBadPixelsPix::kRelTimeOscillating
504//
505void MHCalibrationRelTimeCam::FinalizeBadPixels()
506{
507
508 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
509
510 for (Int_t i=0; i<badcam->GetSize(); i++)
511 {
512 MBadPixelsPix &bad = (*badcam)[i];
513
514 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeNotFitted ))
515 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
516
517 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeOscillating))
518 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
519
520 }
521}
522
523// --------------------------------------------------------------------------
524//
525// The types are as follows:
526//
527// Fitted values:
528// ==============
529//
530// 0: Fitted Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMean()
531// 1: Error Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMeanErr()
532// 2: Sigma fitted Relative Arrival Time in FADC slices (MHGausEvents::GetSigma()
533// 3: Error Sigma Relative Arrival Time in FADC slices (MHGausEvents::GetSigmaErr()
534//
535// Useful variables derived from the fit results:
536// =============================================
537//
538// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
539//
540// Localized defects:
541// ==================
542//
543// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
544// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
545//
546Bool_t MHCalibrationRelTimeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
547{
548
549 if (fHiGainArray->GetSize() <= idx)
550 return kFALSE;
551
552 const MHCalibrationPix &pix = (*this)[idx];
553
554 switch (type)
555 {
556 case 0:
557 val = pix.GetMean();
558 break;
559 case 1:
560 val = pix.GetMeanErr();
561 break;
562 case 2:
563 val = pix.GetSigma();
564 break;
565 case 3:
566 val = pix.GetSigmaErr();
567 break;
568 case 4:
569 val = pix.GetProb();
570 break;
571 case 5:
572 if (!pix.IsGausFitOK())
573 val = 1.;
574 break;
575 case 6:
576 if (!pix.IsFourierSpectrumOK())
577 val = 1.;
578 break;
579 default:
580 return kFALSE;
581 }
582 return kTRUE;
583}
584
585// --------------------------------------------------------------------------
586//
587// Calls MHCalibrationPix::DrawClone() for pixel idx
588//
589void MHCalibrationRelTimeCam::DrawPixelContent(Int_t idx) const
590{
591 (*this)[idx].DrawClone();
592}
593
594// -----------------------------------------------------------------------------
595//
596// Default draw:
597//
598// Displays the averaged areas, both High Gain and Low Gain
599//
600// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
601//
602void MHCalibrationRelTimeCam::Draw(const Option_t *opt)
603{
604
605 const Int_t nareas = fAverageHiGainAreas->GetSize();
606 if (nareas == 0)
607 return;
608
609 TString option(opt);
610 option.ToLower();
611
612 if (!option.Contains("datacheck"))
613 {
614 MHCalibrationCam::Draw(opt);
615 return;
616 }
617
618 //
619 // From here on , the datacheck - Draw
620 //
621 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
622 pad->SetBorderMode(0);
623 pad->Divide(1,nareas);
624
625 //
626 // Loop over inner and outer pixels
627 //
628 for (Int_t i=0; i<nareas;i++)
629 {
630
631 pad->cd(i+1);
632
633 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
634 //
635 // Ask for Hi-Gain saturation
636 //
637 if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
638 {
639 MHCalibrationPix &lopix = GetAverageLoGainArea(i);
640 DrawDataCheckPixel(lopix,0.);
641 }
642 else
643 DrawDataCheckPixel(hipix,0.);
644 }
645}
646
647// -----------------------------------------------------------------------------
648//
649// Draw the average pixel for the datacheck:
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::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
656{
657
658 TVirtualPad *newpad = gPad;
659 newpad->Divide(1,2);
660 newpad->cd(1);
661
662 gPad->SetTicks();
663 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
664 gPad->SetLogy();
665
666 TH1F *hist = pix.GetHGausHist();
667
668 TH1F *null = new TH1F("Null",hist->GetTitle(),100,pix.GetFirst() < -3.0 ? -3.0 : pix.GetFirst(),
669 pix.GetLast() > 3.0 ? 3.0 : pix.GetLast());
670
671 null->SetMaximum(1.1*hist->GetMaximum());
672 null->SetDirectory(NULL);
673 null->SetBit(kCanDelete);
674 null->SetStats(kFALSE);
675 //
676 // set the labels bigger
677 //
678 TAxis *xaxe = null->GetXaxis();
679 TAxis *yaxe = null->GetYaxis();
680 xaxe->CenterTitle();
681 yaxe->CenterTitle();
682 xaxe->SetTitleSize(0.07);
683 yaxe->SetTitleSize(0.07);
684 xaxe->SetTitleOffset(0.7);
685 yaxe->SetTitleOffset(0.55);
686 xaxe->SetLabelSize(0.06);
687 yaxe->SetLabelSize(0.06);
688
689 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
690 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
691
692 null->Draw();
693 hist->Draw("same");
694
695 gStyle->SetOptFit();
696
697 TF1 *fit = pix.GetFGausFit();
698
699 if (fit)
700 {
701 switch ( fColor )
702 {
703 case MCalibrationCam::kGREEN:
704 fit->SetLineColor(kGreen);
705 break;
706 case MCalibrationCam::kBLUE:
707 fit->SetLineColor(kBlue);
708 break;
709 case MCalibrationCam::kUV:
710 fit->SetLineColor(106);
711 break;
712 case MCalibrationCam::kCT1:
713 fit->SetLineColor(006);
714 break;
715 default:
716 fit->SetLineColor(kRed);
717 }
718 fit->Draw("same");
719 }
720
721 // DisplayRefLines(null,refline);
722
723 gPad->Modified();
724 gPad->Update();
725
726 newpad->cd(2);
727 gPad->SetTicks();
728
729 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
730
731 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
732 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
733 null2->SetDirectory(NULL);
734 null2->SetBit(kCanDelete);
735 null2->SetStats(kFALSE);
736 //
737 // set the labels bigger
738 //
739 TAxis *xaxe2 = null2->GetXaxis();
740 TAxis *yaxe2 = null2->GetYaxis();
741 xaxe2->CenterTitle();
742 yaxe2->CenterTitle();
743 xaxe2->SetTitleSize(0.07);
744 yaxe2->SetTitleSize(0.07);
745 xaxe2->SetTitleOffset(0.7);
746 yaxe2->SetTitleOffset(0.55);
747 xaxe2->SetLabelSize(0.06);
748 yaxe2->SetLabelSize(0.06);
749
750 pix.CreateGraphEvents();
751 TGraph *gr = pix.GetGraphEvents();
752
753 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
754 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
755
756 null2->Draw();
757
758 pix.DrawEvents("same");
759
760 return;
761}
762
Note: See TracBrowser for help on using the repository browser.