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

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