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

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