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

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