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

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