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

Last change on this file since 5071 was 5014, checked in by gaug, 21 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//
177TObject *MHCalibrationRelTimeCam::Clone(const char *name) const
178{
179
180 const Int_t navhi = fAverageHiGainAreas->GetEntries();
181 const Int_t navlo = fAverageLoGainAreas->GetEntries();
182 const Int_t nsehi = fAverageHiGainSectors->GetEntries();
183 const Int_t nselo = fAverageLoGainSectors->GetEntries();
184
185 //
186 // FIXME, this might be done faster and more elegant, by direct copy.
187 //
188 MHCalibrationRelTimeCam *cam = new MHCalibrationRelTimeCam();
189
190 cam->fAverageHiGainAreas->Expand(navhi);
191 cam->fAverageHiGainSectors->Expand(nsehi);
192
193 for (int i=0; i<navhi; i++)
194 (*cam->fAverageHiGainAreas) [i] = (*fAverageHiGainAreas) [i]->Clone();
195 for (int i=0; i<nsehi; i++)
196 (*cam->fAverageHiGainSectors)[i] = (*fAverageHiGainSectors)[i]->Clone();
197
198 if (IsLoGain())
199 {
200 cam->fAverageLoGainAreas->Expand(navlo);
201 cam->fAverageLoGainSectors->Expand(nselo);
202
203 for (int i=0; i<navlo; i++)
204 (*cam->fAverageLoGainAreas) [i] = (*fAverageLoGainAreas) [i]->Clone();
205 for (int i=0; i<nselo; i++)
206 (*cam->fAverageLoGainSectors)[i] = (*fAverageLoGainSectors)[i]->Clone();
207 }
208
209 cam->fAverageAreaNum = fAverageAreaNum;
210 cam->fAverageAreaSat = fAverageAreaSat;
211 cam->fAverageAreaSigma = fAverageAreaSigma;
212 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
213 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
214 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
215 cam->fAverageSectorNum = fAverageSectorNum;
216 cam->fRunNumbers = fRunNumbers;
217
218 cam->fColor = fColor;
219 cam->fPulserFrequency = fPulserFrequency;
220 cam->fFlags = fFlags;
221
222 return cam;
223
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 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
410 {
411
412 MHCalibrationPix &histhi = (*this)[i];
413
414 if (histhi.IsExcluded())
415 continue;
416
417 MCalibrationRelTimePix &pix = fIntensCam
418 ? (MCalibrationRelTimePix&)(*fIntensCam)[i]
419 : (MCalibrationRelTimePix&)(*fCam)[i];
420
421 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
422 {
423 pix.SetHiGainSaturation();
424 histhi.SetExcluded();
425 }
426 else
427 if (IsLoGain())
428 (*this)(i).SetExcluded();
429
430 Stat_t overflow = histhi.GetHGausHist()->GetBinContent(histhi.GetHGausHist()->GetNbinsX()+1);
431 if (overflow > 0.1)
432 {
433 *fLog << warn << GetDescriptor()
434 << ": HiGain Histogram Overflow occurred " << overflow
435 << " times in pixel: " << i << " (without saturation!) " << endl;
436 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
437 }
438
439 overflow = histhi.GetHGausHist()->GetBinContent(0);
440 if (overflow > 0.1)
441 {
442 *fLog << warn << GetDescriptor()
443 << ": HiGain Histogram Underflow occurred " << overflow
444 << " times in pixel: " << i << " (without saturation!) " << endl;
445 // bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
446 }
447 }
448
449 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
450 {
451
452 MHCalibrationPix &histhi = GetAverageHiGainArea(j);
453
454 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
455 {
456 MCalibrationRelTimePix &pix = fIntensCam
457 ? (MCalibrationRelTimePix&)fIntensCam->GetAverageArea(j)
458 : (MCalibrationRelTimePix&)fCam->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 = fIntensCam
473 ? (MCalibrationRelTimePix&)fIntensCam->GetAverageSector(j)
474 : (MCalibrationRelTimePix&)fCam->GetAverageSector(j);
475
476 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
477 {
478 pix.SetHiGainSaturation();
479 histhi.SetExcluded();
480 }
481 else
482 if (IsLoGain())
483 GetAverageLoGainSector(j).SetExcluded();
484 }
485
486 FitHiGainArrays(fIntensCam ? (MCalibrationCam&)(*fIntensCam->GetCam()) : (MCalibrationCam&)(*fCam),
487 fIntensBad ? (*fIntensBad->GetCam()) : *fBadPixels,
488 MBadPixelsPix::kRelTimeNotFitted,
489 MBadPixelsPix::kRelTimeOscillating);
490
491 if (IsLoGain())
492 FitLoGainArrays(fIntensCam ? (MCalibrationCam&)(*fIntensCam->GetCam()) : (MCalibrationCam&)(*fCam),
493 fIntensBad ? (*fIntensBad->GetCam()) : *fBadPixels,
494 MBadPixelsPix::kRelTimeNotFitted,
495 MBadPixelsPix::kRelTimeOscillating);
496
497 return kTRUE;
498}
499
500// --------------------------------------------------------------------------
501//
502// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
503// - MBadPixelsPix::kRelTimeNotFitted
504// - MBadPixelsPix::kRelTimeOscillating
505//
506void MHCalibrationRelTimeCam::FinalizeBadPixels()
507{
508
509 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
510
511 for (Int_t i=0; i<badcam->GetSize(); i++)
512 {
513 MBadPixelsPix &bad = (*badcam)[i];
514
515 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeNotFitted ))
516 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
517
518 if (bad.IsUncalibrated( MBadPixelsPix::kRelTimeOscillating))
519 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
520
521 }
522}
523
524// --------------------------------------------------------------------------
525//
526// The types are as follows:
527//
528// Fitted values:
529// ==============
530//
531// 0: Fitted Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMean()
532// 1: Error Mean Relative Arrival Time in FADC slices (MHGausEvents::GetMeanErr()
533// 2: Sigma fitted Relative Arrival Time in FADC slices (MHGausEvents::GetSigma()
534// 3: Error Sigma Relative Arrival Time in FADC slices (MHGausEvents::GetSigmaErr()
535//
536// Useful variables derived from the fit results:
537// =============================================
538//
539// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
540//
541// Localized defects:
542// ==================
543//
544// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
545// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
546//
547Bool_t MHCalibrationRelTimeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
548{
549
550 if (fHiGainArray->GetSize() <= idx)
551 return kFALSE;
552
553 const MHCalibrationPix &pix = (*this)[idx];
554
555 switch (type)
556 {
557 case 0:
558 val = pix.GetMean();
559 break;
560 case 1:
561 val = pix.GetMeanErr();
562 break;
563 case 2:
564 val = pix.GetSigma();
565 break;
566 case 3:
567 val = pix.GetSigmaErr();
568 break;
569 case 4:
570 val = pix.GetProb();
571 break;
572 case 5:
573 if (!pix.IsGausFitOK())
574 val = 1.;
575 break;
576 case 6:
577 if (!pix.IsFourierSpectrumOK())
578 val = 1.;
579 break;
580 default:
581 return kFALSE;
582 }
583 return kTRUE;
584}
585
586// --------------------------------------------------------------------------
587//
588// Calls MHCalibrationPix::DrawClone() for pixel idx
589//
590void MHCalibrationRelTimeCam::DrawPixelContent(Int_t idx) const
591{
592 (*this)[idx].DrawClone();
593}
594
595// -----------------------------------------------------------------------------
596//
597// Default draw:
598//
599// Displays the averaged areas, both High Gain and Low Gain
600//
601// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
602//
603void MHCalibrationRelTimeCam::Draw(const Option_t *opt)
604{
605
606 const Int_t nareas = fAverageHiGainAreas->GetEntries();
607 if (nareas == 0)
608 return;
609
610 TString option(opt);
611 option.ToLower();
612
613 if (!option.Contains("datacheck"))
614 {
615 MHCalibrationCam::Draw(opt);
616 return;
617 }
618
619 //
620 // From here on , the datacheck - Draw
621 //
622 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
623 pad->SetBorderMode(0);
624 pad->Divide(1,nareas);
625
626 //
627 // Loop over inner and outer pixels
628 //
629 for (Int_t i=0; i<nareas;i++)
630 {
631
632 pad->cd(i+1);
633
634 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
635 //
636 // Ask for Hi-Gain saturation
637 //
638 if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
639 {
640 MHCalibrationPix &lopix = GetAverageLoGainArea(i);
641 DrawDataCheckPixel(lopix,0.);
642 }
643 else
644 DrawDataCheckPixel(hipix,0.);
645 }
646}
647
648// -----------------------------------------------------------------------------
649//
650// Draw the average pixel for the datacheck:
651//
652// Displays the averaged areas, both High Gain and Low Gain
653//
654// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
655//
656void MHCalibrationRelTimeCam::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
657{
658
659 TVirtualPad *newpad = gPad;
660 newpad->Divide(1,2);
661 newpad->cd(1);
662
663 gPad->SetTicks();
664 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
665 gPad->SetLogy();
666
667 TH1F *hist = pix.GetHGausHist();
668
669 TH1F *null = new TH1F("Null",hist->GetTitle(),100,pix.GetFirst() < -3.0 ? -3.0 : pix.GetFirst(),
670 pix.GetLast() > 3.0 ? 3.0 : pix.GetLast());
671
672 null->SetMaximum(1.1*hist->GetMaximum());
673 null->SetDirectory(NULL);
674 null->SetBit(kCanDelete);
675 null->SetStats(kFALSE);
676 //
677 // set the labels bigger
678 //
679 TAxis *xaxe = null->GetXaxis();
680 TAxis *yaxe = null->GetYaxis();
681 xaxe->CenterTitle();
682 yaxe->CenterTitle();
683 xaxe->SetTitleSize(0.07);
684 yaxe->SetTitleSize(0.07);
685 xaxe->SetTitleOffset(0.7);
686 yaxe->SetTitleOffset(0.55);
687 xaxe->SetLabelSize(0.06);
688 yaxe->SetLabelSize(0.06);
689
690 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
691 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
692
693 null->Draw();
694 hist->Draw("same");
695
696 gStyle->SetOptFit();
697
698 TF1 *fit = pix.GetFGausFit();
699
700 if (fit)
701 {
702 switch ( fColor )
703 {
704 case MCalibrationCam::kGREEN:
705 fit->SetLineColor(kGreen);
706 break;
707 case MCalibrationCam::kBLUE:
708 fit->SetLineColor(kBlue);
709 break;
710 case MCalibrationCam::kUV:
711 fit->SetLineColor(106);
712 break;
713 case MCalibrationCam::kCT1:
714 fit->SetLineColor(006);
715 break;
716 default:
717 fit->SetLineColor(kRed);
718 }
719 fit->Draw("same");
720 }
721
722 // DisplayRefLines(null,refline);
723
724 gPad->Modified();
725 gPad->Update();
726
727 newpad->cd(2);
728 gPad->SetTicks();
729
730 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
731
732 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
733 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
734 null2->SetDirectory(NULL);
735 null2->SetBit(kCanDelete);
736 null2->SetStats(kFALSE);
737 //
738 // set the labels bigger
739 //
740 TAxis *xaxe2 = null2->GetXaxis();
741 TAxis *yaxe2 = null2->GetYaxis();
742 xaxe2->CenterTitle();
743 yaxe2->CenterTitle();
744 xaxe2->SetTitleSize(0.07);
745 yaxe2->SetTitleSize(0.07);
746 xaxe2->SetTitleOffset(0.7);
747 yaxe2->SetTitleOffset(0.55);
748 xaxe2->SetLabelSize(0.06);
749 yaxe2->SetLabelSize(0.06);
750
751 pix.CreateGraphEvents();
752 TGraph *gr = pix.GetGraphEvents();
753
754 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
755 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
756
757 null2->Draw();
758
759 pix.DrawEvents("same");
760
761 return;
762}
763
Note: See TracBrowser for help on using the repository browser.