source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationCam.cc@ 5026

Last change on this file since 5026 was 5026, checked in by gaug, 20 years ago
*** empty log message ***
File size: 43.5 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// MHCalibrationCam
27//
28// Base class for camera calibration classes. Incorporates the TObjArray's:
29// - fHiGainArray (for calibrated High Gains per pixel)
30// - fLoGainArray (for calibrated Low Gains per pixel)
31// - fAverageHiGainAreas (for averaged High Gains events per camera area index)
32// - fAverageLoGainAreas (for averaged High Gains events per camera area index)
33// - fAverageHiGainSectors (for averaged High Gains events per camera sector )
34// - fAverageLoGainSectors (for averaged High Gains events per camera sector )
35// These TObjArray's are called by their default constructors, thus no objects
36// are created, until the derived class does so.
37//
38// The corresponding operators: [],() and the operators GetAverageHiGainArea(),
39// GetAverageLoGainArea(), GetAverageHiGainSector() and GetAverageLoGainSector()
40// have to be cast to the corresponding class. It is assumed that all classes
41// dealing with calibration pixels derive from MHCalibrationPix.
42//
43// The following flag can be set:
44// - SetAverageing() for calculating the event-by-event averages.
45// - SetDebug() for debug output.
46// - SetLoGain() for the case that low-gain slices are available, but
47// MRawRunHeader::GetNumLoGainSlices() gives still 0.
48// - SetCheckSize() for testing the sizes of the f*Arrays in ReInit() if they
49// are coinciding with the equiv. sizes in MGeomCam
50//
51// The flag kLoGain steers if the low-gain signal is treated at all or not.
52// The flag kAverageing steers if the event-by-event averages are treated at all.
53//
54/////////////////////////////////////////////////////////////////////////////
55#include "MHCalibrationCam.h"
56#include "MHCalibrationPix.h"
57
58#include <TVirtualPad.h>
59#include <TCanvas.h>
60#include <TPad.h>
61#include <TText.h>
62#include <TPaveText.h>
63
64#include "MLog.h"
65#include "MLogManip.h"
66
67#include "MCalibrationIntensityCam.h"
68#include "MCalibrationCam.h"
69#include "MCalibrationPix.h"
70
71#include "MBadPixelsIntensityCam.h"
72#include "MBadPixelsCam.h"
73#include "MBadPixelsPix.h"
74
75#include "MGeomCam.h"
76#include "MGeomPix.h"
77
78#include "MParList.h"
79
80#include "MRawRunHeader.h"
81
82ClassImp(MHCalibrationCam);
83
84using namespace std;
85
86const Int_t MHCalibrationCam::fgPulserFrequency = 500;
87const TString MHCalibrationCam::gsHistName = "Hist";
88const TString MHCalibrationCam::gsHistTitle = "";
89const TString MHCalibrationCam::gsHistXTitle = "";
90const TString MHCalibrationCam::gsHistYTitle = "Nr. events";
91// --------------------------------------------------------------------------
92//
93// Default Constructor.
94//
95// Sets:
96// - all pointers to NULL
97//
98// Initializes and sets owner of:
99// - fHiGainArray, fLoGainArray
100// - fAverageHiGainAreas, fAverageLoGainAreas
101// - fAverageHiGainSectors, fAverageLoGainSectors
102//
103// Initializes:
104// - fPulserFrequency to fgPulserFrequency
105//
106// - SetAveregeing (kTRUE);
107// - SetDebug (kFALSE);
108// - SetLoGain (kTRUE);
109//- SetOscillations(kTRUE);
110//- SetSizeCheck (kTRUE);
111//
112MHCalibrationCam::MHCalibrationCam(const char *name, const char *title)
113 : fHistName(gsHistName),fHistTitle(gsHistTitle),
114 fHistXTitle(gsHistXTitle),fHistYTitle(gsHistYTitle),
115 fColor(MCalibrationCam::kNONE), fIntensBad(NULL),
116 fBadPixels(NULL), fIntensCam(NULL), fCam(NULL), fGeom(NULL), fRunHeader(NULL)
117{
118
119 fHiGainArray = new TObjArray;
120 fHiGainArray->SetOwner();
121
122 fLoGainArray = new TObjArray;
123 fLoGainArray->SetOwner();
124
125 fAverageHiGainAreas = new TObjArray;
126 fAverageHiGainAreas->SetOwner();
127
128 fAverageLoGainAreas = new TObjArray;
129 fAverageLoGainAreas->SetOwner();
130
131 fAverageHiGainSectors = new TObjArray;
132 fAverageHiGainSectors->SetOwner();
133
134 fAverageLoGainSectors = new TObjArray;
135 fAverageLoGainSectors->SetOwner();
136
137 SetPulserFrequency();
138
139 SetAverageing (kTRUE);
140 SetDebug (kFALSE);
141 SetLoGain (kTRUE);
142 SetOscillations(kTRUE);
143 SetSizeCheck (kTRUE);
144}
145
146// --------------------------------------------------------------------------
147//
148// Returns size of fHiGainArray
149//
150const Int_t MHCalibrationCam::GetSize() const
151{
152 return fHiGainArray->GetSize();
153}
154
155// --------------------------------------------------------------------------
156//
157// Deletes the TClonesArray of:
158// - fHiGainArray, fLoGainArray
159// - fAverageHiGainAreas, fAverageLoGainAreas
160// - fAverageHiGainSectors, fAverageLoGainSectors
161//
162MHCalibrationCam::~MHCalibrationCam()
163{
164
165 delete fHiGainArray;
166 delete fLoGainArray;
167
168 delete fAverageHiGainAreas;
169 delete fAverageLoGainAreas;
170
171 delete fAverageHiGainSectors;
172 delete fAverageLoGainSectors;
173}
174
175// --------------------------------------------------------------------------
176//
177// Calls Reset() for each entry in:
178// - fHiGainArray, fLoGainArray
179// - fAverageHiGainAreas, fAverageLoGainAreas
180// - fAverageHiGainSectors, fAverageLoGainSectors
181//
182void MHCalibrationCam::ResetHists()
183{
184
185 if (fHiGainArray)
186 { fHiGainArray->ForEach(MHCalibrationPix,Reset)(); }
187 if (IsAverageing())
188 {
189 if (fAverageHiGainAreas)
190 { fAverageHiGainAreas->ForEach(MHCalibrationPix,Reset)(); }
191 if (fAverageHiGainSectors)
192 { fAverageHiGainSectors->ForEach(MHCalibrationPix,Reset)(); }
193 }
194
195 if (!IsLoGain())
196 return;
197
198 if (fLoGainArray)
199 { fLoGainArray->ForEach(MHCalibrationPix,Reset)(); }
200 if (IsAverageing())
201 {
202 if (fAverageLoGainAreas)
203 { fAverageLoGainAreas->ForEach(MHCalibrationPix,Reset)(); }
204 if (fAverageLoGainSectors)
205 { fAverageLoGainSectors->ForEach(MHCalibrationPix,Reset)(); }
206 }
207}
208
209// --------------------------------------------------------------------------
210//
211// Get i-th High Gain pixel (pixel number)
212//
213MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i)
214{
215 return *static_cast<MHCalibrationPix*>(fHiGainArray->UncheckedAt(i));
216}
217
218// --------------------------------------------------------------------------
219//
220// Get i-th High Gain pixel (pixel number)
221//
222const MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i) const
223{
224 return *static_cast<MHCalibrationPix*>(fHiGainArray->UncheckedAt(i));
225}
226
227// --------------------------------------------------------------------------
228//
229// Get i-th Low Gain pixel (pixel number)
230//
231MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i)
232{
233 return *static_cast<MHCalibrationPix*>(fLoGainArray->UncheckedAt(i));
234}
235
236// --------------------------------------------------------------------------
237//
238// Get i-th Low Gain pixel (pixel number)
239//
240const MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i) const
241{
242 return *static_cast<MHCalibrationPix*>(fLoGainArray->UncheckedAt(i));
243}
244
245// --------------------------------------------------------------------------
246//
247// Returns the current size of the TObjArray fAverageHiGainAreas
248// independently if the MHCalibrationPix is filled with values or not.
249//
250const Int_t MHCalibrationCam::GetAverageAreas() const
251{
252 return fAverageHiGainAreas->GetEntries();
253}
254
255// --------------------------------------------------------------------------
256//
257// Get i-th High Gain pixel Area (area number)
258//
259MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i)
260{
261 return *static_cast<MHCalibrationPix*>(fAverageHiGainAreas->UncheckedAt(i));
262}
263
264// --------------------------------------------------------------------------
265//
266// Get i-th High Gain pixel Area (area number)
267//
268const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i) const
269{
270 return *static_cast<MHCalibrationPix *>(fAverageHiGainAreas->UncheckedAt(i));
271}
272
273// --------------------------------------------------------------------------
274//
275// Get i-th Low Gain pixel Area (area number)
276//
277MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i)
278{
279 return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->UncheckedAt(i));
280}
281
282// --------------------------------------------------------------------------
283//
284// Get i-th Low Gain pixel Area (area number)
285//
286const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i) const
287{
288 return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->UncheckedAt(i));
289}
290
291// --------------------------------------------------------------------------
292//
293// Returns the current size of the TObjArray fAverageHiGainSectors
294// independently if the MHCalibrationPix is filled with values or not.
295//
296const Int_t MHCalibrationCam::GetAverageSectors() const
297{
298 return fAverageHiGainSectors->GetEntries();
299}
300
301// --------------------------------------------------------------------------
302//
303// Get i-th High Gain Sector (sector number)
304//
305MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i)
306{
307 return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->UncheckedAt(i));
308}
309
310// --------------------------------------------------------------------------
311//
312// Get i-th High Gain Sector (sector number)
313//
314const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i) const
315{
316 return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->UncheckedAt(i));
317}
318
319// --------------------------------------------------------------------------
320//
321// Get i-th Low Gain Sector (sector number)
322//
323MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i)
324{
325 return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->UncheckedAt(i));
326}
327
328// --------------------------------------------------------------------------
329//
330// Get i-th Low Gain Sector (sector number)
331//
332const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i) const
333{
334 return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->UncheckedAt(i));
335}
336
337
338// --------------------------------------------------------------------------
339//
340// Creates new MHCalibrationCam only for the averaged areas:
341// the rest has to be retrieved directly, e.g. via:
342// MHCalibrationCam *cam = MParList::FindObject("MHCalibrationCam");
343// - cam->GetAverageSector(5).DrawClone();
344// - (*cam)[100].DrawClone()
345//
346TObject *MHCalibrationCam::Clone(const char *) const
347{
348
349 //
350 // FIXME, this might be done faster and more elegant, by direct copy.
351 //
352 MHCalibrationCam *cam = new MHCalibrationCam();
353
354 cam->fRunNumbers = fRunNumbers;
355 cam->fPulserFrequency = fPulserFrequency;
356 cam->fFlags = fFlags;
357 cam->fNbins = fNbins;
358 cam->fFirst = fFirst;
359 cam->fLast = fLast;
360
361 if (!IsAverageing())
362 return cam;
363
364 // const Int_t nhi = fHiGainArray->GetEntries();
365 // const Int_t nlo = fLoGainArray->GetEntries();
366 const Int_t navhi = fAverageHiGainAreas->GetEntries();
367 // const Int_t nsehi = fAverageHiGainSectors->GetEntries();
368
369 // cam->fHiGainArray->Expand(nhi);
370 // cam->fLoGainArray->Expand(nlo);
371 cam->fAverageHiGainAreas->Expand(navhi);
372 // cam->fAverageHiGainSectors->Expand(nsehi);
373 // cam->fAverageLoGainSectors->Expand(nselo);
374
375 /*
376 for (int i=0; i<nhi; i++)
377 {
378 delete (*cam->fHiGainArray)[i];
379 (*cam->fHiGainArray)[i] = (*fHiGainArray)[i]->Clone();
380 }
381 for (int i=0; i<nlo; i++)
382 {
383 delete (*cam->fLoGainArray)[i];
384 (*cam->fLoGainArray)[i] = (*fLoGainArray)[i]->Clone();
385 }
386 */
387
388 for (int i=0; i<navhi; i++)
389 (*cam->fAverageHiGainAreas)[i] = (*fAverageHiGainAreas)[i]->Clone();
390
391 // for (int i=0; i<nsehi; i++)
392 // (*cam->fAverageHiGainSectors)[i] = (*fAverageHiGainSectors)[i]->Clone();
393
394 if (IsLoGain())
395 {
396
397 const Int_t navlo = fAverageLoGainAreas->GetEntries();
398
399 cam->fAverageLoGainAreas->Expand(navlo);
400 // cam->fAverageLoGainSectors->Expand(nselo);
401 // const Int_t nselo = fAverageLoGainSectors->GetEntries();
402
403 for (int i=0; i<navlo; i++)
404 (*cam->fAverageLoGainAreas)[i] = (*fAverageLoGainAreas)[i]->Clone();
405
406 // for (int i=0; i<nselo; i++)
407 // (*cam->fAverageLoGainSectors)[i] = (*fAverageLoGainSectors)[i]->Clone();
408 }
409
410 cam->fAverageAreaNum = fAverageAreaNum;
411 cam->fAverageAreaSat = fAverageAreaSat;
412 cam->fAverageAreaSigma = fAverageAreaSigma;
413 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
414 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
415 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
416 cam->fAverageSectorNum = fAverageSectorNum;
417
418 return cam;
419}
420
421// --------------------------------------------------------------------------
422//
423// Gets the pointers to:
424// - MGeomCam
425//
426// Calls SetupHists(const MParList *pList)
427//
428// Calls Delete-Function of:
429// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
430// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
431// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
432//
433Bool_t MHCalibrationCam::SetupFill(const MParList *pList)
434{
435
436 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
437 if (!fGeom)
438 {
439 *fLog << err << GetDescriptor()
440 << ": MGeomCam not found... aborting." << endl;
441 return kFALSE;
442 }
443
444 fRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
445 if (!fRunHeader)
446 {
447 *fLog << warn << GetDescriptor()
448 << ": MRawRunHeader not found... will not store run numbers." << endl;
449 }
450
451 /* not needed -> Class is never used twice
452 fHiGainArray->Delete();
453 fLoGainArray->Delete();
454
455 fAverageHiGainAreas->Delete();
456 fAverageLoGainAreas->Delete();
457
458 fAverageHiGainSectors->Delete();
459 fAverageLoGainSectors->Delete();
460 */
461
462 return SetupHists(pList);
463}
464
465
466// --------------------------------------------------------------------------
467//
468// Gets or creates the pointers to:
469// - MBadPixelsIntensityCam
470// - MBadPixelsCam
471//
472// Searches pointer to:
473// - MArrivalTimeCam
474//
475// Initializes, if empty to MArrivalTimeCam::GetSize() for:
476// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
477//
478// Initializes, if empty to MGeomCam::GetNumAreas() for:
479// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
480//
481// Initializes, if empty to MGeomCam::GetNumSectors() for:
482// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
483//
484// Initializes TArray's to MGeomCam::GetNumAreas and MGeomCam::GetNumSectors, respectively
485// Fills with number of valid pixels (if !MBadPixelsPix::IsBad()):
486// - MHCalibrationCam::fAverageAreaNum[area index]
487// - MHCalibrationCam::fAverageSectorNum[area index]
488//
489// Calls InitializeHists() for every entry in:
490// - MHCalibrationCam::fHiGainArray
491// - MHCalibrationCam::fAverageHiGainAreas
492// - MHCalibrationCam::fAverageHiGainSectors
493//
494// Sets Titles and Names for the Histograms
495// - MHCalibrationCam::fAverageHiGainAreas
496// - MHCalibrationCam::fAverageHiGainSectors
497//
498// Retrieves the run numbers from MRawRunHeader and stores them in fRunNumbers
499//
500Bool_t MHCalibrationCam::ReInit(MParList *pList)
501{
502
503 const Int_t npixels = fGeom->GetNumPixels();
504 const Int_t nsectors = fGeom->GetNumSectors();
505 const Int_t nareas = fGeom->GetNumAreas();
506
507 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
508 if (fIntensBad)
509 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
510 else
511 {
512 fBadPixels = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
513 if (!fBadPixels)
514 {
515
516 fBadPixels = (MBadPixelsCam*)pList->FindCreateObj(AddSerialNumber("MBadPixelsCam"));
517 if (!fBadPixels)
518 {
519 *fLog << err << "Cannot find nor create MBadPixelsCam ... abort." << endl;
520 return kFALSE;
521 }
522 else
523 fBadPixels->InitSize(npixels);
524 }
525 }
526
527 if (IsAverageing())
528 {
529 //
530 // The function TArrayF::Set() already sets all entries to 0.
531 //
532 fAverageAreaNum. Set(nareas);
533 fAverageAreaSat. Set(nareas);
534 fAverageAreaSigma. Set(nareas);
535 fAverageAreaSigmaVar. Set(nareas);
536 fAverageAreaRelSigma. Set(nareas);
537 fAverageAreaRelSigmaVar.Set(nareas);
538 fAverageSectorNum. Set(nsectors);
539
540 for (Int_t aidx=0; aidx<nareas; aidx++)
541 fAverageAreaNum[aidx] = 0;
542
543 for (Int_t sector=0; sector<nsectors; sector++)
544 fAverageSectorNum[sector] = 0;
545
546 for (Int_t i=0; i<npixels; i++)
547 {
548
549 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
550 if (bad.IsBad())
551 continue;
552
553 fAverageAreaNum [(*fGeom)[i].GetAidx() ]++;
554 fAverageSectorNum[(*fGeom)[i].GetSector()]++;
555 }
556 }
557
558 //
559 // Because ReInit has been called, a new run number is added
560 //
561 fRunNumbers.Set(fRunNumbers.GetSize()+1);
562
563 if (fRunHeader)
564 {
565 fRunNumbers[fRunNumbers.GetSize()-1] = fRunHeader->GetRunNumber();
566 if (IsLoGain())
567 SetLoGain(fRunHeader->GetNumSamplesLoGain());
568 }
569
570 if (!ReInitHists(pList))
571 return kFALSE;
572
573 if (!fRunHeader)
574 return kTRUE;
575
576 for (Int_t i=0; i<fHiGainArray->GetEntries(); i++)
577 {
578 TH1F *h = (*this)[i].GetHGausHist();
579 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
580 }
581
582 if (IsLoGain())
583 for (Int_t i=0; i<fLoGainArray->GetEntries(); i++)
584 {
585 TH1F *h = (*this)(i).GetHGausHist();
586 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
587 }
588
589 if (!IsAverageing())
590 return kTRUE;
591
592 for (Int_t j=0; j<nareas; j++)
593 {
594 TH1F *h = GetAverageHiGainArea(j).GetHGausHist();
595 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
596 }
597
598 if (IsLoGain())
599 for (Int_t j=0; j<nareas; j++)
600 {
601 TH1F *h = GetAverageLoGainArea(j).GetHGausHist();
602 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
603 }
604
605 for (Int_t j=0; j<nsectors; j++)
606 {
607 TH1F *h = GetAverageHiGainSector(j).GetHGausHist();
608 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
609 }
610
611 if (IsLoGain())
612 for (Int_t j=0; j<nsectors; j++)
613 {
614 TH1F *h = GetAverageLoGainSector(j).GetHGausHist();
615 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
616 }
617
618 return kTRUE;
619}
620
621//--------------------------------------------------------------------------------------
622//
623// Initializes the High Gain Arrays:
624//
625// - Expand fHiGainArrays to npixels
626// - Expand fAverageHiGainAreas to nareas
627// - Expand fAverageHiGainSectors to nsectors
628//
629// - For every entry in the expanded arrays:
630// * Initialize an MHCalibrationPix
631// * Set Binning from fNbins, fFirst and fLast
632// * Set Histgram names and titles from fHistName and fHistTitle
633// * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
634// * Call InitHists
635//
636void MHCalibrationCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
637{
638
639 if (fHiGainArray->GetEntries()==0)
640 {
641 fHiGainArray->Expand(npixels);
642 for (Int_t i=0; i<npixels; i++)
643 {
644 (*fHiGainArray)[i] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"HiGainPix"),
645 Form("%s%s",fHistTitle.Data()," High Gain Pixel"));
646
647 MHCalibrationPix &pix = (*this)[i];
648 pix.SetNbins(fNbins);
649 pix.SetFirst(fFirst);
650 pix.SetLast (fLast);
651
652 TH1F *h = pix.GetHGausHist();
653
654 h->SetName (Form("%s%s%s","H",fHistName.Data(),"HiGainPix"));
655 h->SetTitle(Form("%s%s",fHistTitle.Data()," High Gain Pixel "));
656 h->SetXTitle(fHistXTitle.Data());
657 h->SetYTitle(fHistYTitle.Data());
658
659 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
660 InitHists(pix,bad,i);
661 }
662 }
663
664 if (!IsAverageing())
665 return;
666
667 if (fAverageHiGainAreas->GetEntries()==0)
668 {
669 fAverageHiGainAreas->Expand(nareas);
670
671 for (Int_t j=0; j<nareas; j++)
672 {
673 (*fAverageHiGainAreas)[j] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"HiGainArea"),
674 Form("%s%s",fHistTitle.Data()," High Gain Area Idx "));
675
676 MHCalibrationPix &pix = GetAverageHiGainArea(j);
677
678 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
679 pix.SetFirst(fFirst);
680 pix.SetLast (fLast);
681
682 TH1F *h = pix.GetHGausHist();
683
684 h->SetName (Form("%s%s%s","H",fHistName.Data(),"HiGainArea"));
685 h->SetXTitle(fHistXTitle.Data());
686 h->SetYTitle(fHistYTitle.Data());
687
688 if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
689 {
690 h->SetTitle(Form("%s%s%s%s",fHistTitle.Data()," averaged on event-by-event basis ",
691 j==0 ? "Inner Pixels " : "Outer Pixels ","High Gain Runs: "));
692 pix.InitBins();
693 pix.SetEventFrequency(fPulserFrequency);
694 }
695 else
696 {
697 h->SetTitle(Form("%s%s",fHistTitle.Data()," averaged on event-by-event basis High Gain Area Idx "));
698 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
699 }
700 }
701 }
702
703 if (fAverageHiGainSectors->GetEntries()==0)
704 {
705 fAverageHiGainSectors->Expand(nsectors);
706
707 for (Int_t j=0; j<nsectors; j++)
708 {
709 (*fAverageHiGainSectors)[j] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"HiGainSector"),
710 Form("%s%s",fHistTitle.Data()," High Gain Sector "));
711 MHCalibrationPix &pix = GetAverageHiGainSector(j);
712
713 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
714 pix.SetFirst(fFirst);
715 pix.SetLast (fLast);
716
717 TH1F *h = pix.GetHGausHist();
718
719 h->SetName (Form("%s%s%s","H",fHistName.Data(),"HiGainSector"));
720 h->SetTitle(Form("%s%s",fHistTitle.Data()," High Gain Sector "));
721 h->SetXTitle(fHistXTitle.Data());
722 h->SetYTitle(fHistYTitle.Data());
723
724 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
725 }
726 }
727}
728
729//--------------------------------------------------------------------------------------
730//
731// Return, if IsLoGain() is kFALSE
732//
733// Initializes the Low Gain Arrays:
734//
735// - Expand fLoGainArrays to npixels
736// - Expand fAverageLoGainAreas to nareas
737// - Expand fAverageLoGainSectors to nsectors
738//
739// - For every entry in the expanded arrays:
740// * Initialize an MHCalibrationPix
741// * Set Binning from fNbins, fFirst and fLast
742// * Set Histgram names and titles from fHistName and fHistTitle
743// * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
744// * Call InitHists
745//
746void MHCalibrationCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
747{
748
749 if (!IsLoGain())
750 return;
751
752 if (fLoGainArray->GetEntries()==0)
753 {
754 fLoGainArray->Expand(npixels);
755 for (Int_t i=0; i<npixels; i++)
756 {
757 (*fLoGainArray)[i] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"LoGainPix"),
758 Form("%s%s",fHistTitle.Data()," Low Gain Pixel"));
759
760 MHCalibrationPix &pix = (*this)[i];
761 pix.SetNbins(fNbins);
762 pix.SetFirst(fFirst);
763 pix.SetLast (fLast);
764
765 TH1F *h = pix.GetHGausHist();
766
767 h->SetName (Form("%s%s%s","H",fHistName.Data(),"LoGainPix"));
768 h->SetTitle(Form("%s%s",fHistTitle.Data()," Low Gain Pixel "));
769 h->SetXTitle(fHistXTitle.Data());
770 h->SetYTitle(fHistYTitle.Data());
771
772 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
773 InitHists(pix,bad,i);
774 }
775 }
776
777 if (!IsAverageing())
778 return;
779
780 if (fAverageLoGainAreas->GetEntries()==0)
781 {
782 fAverageLoGainAreas->Expand(nareas);
783
784 for (Int_t j=0; j<nareas; j++)
785 {
786 (*fAverageLoGainAreas)[j] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"LoGainArea"),
787 Form("%s%s",fHistTitle.Data()," Low Gain Area Idx "));
788
789 MHCalibrationPix &pix = GetAverageLoGainArea(j);
790
791 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
792 pix.SetFirst(fFirst);
793 pix.SetLast (fLast);
794
795 TH1F *h = pix.GetHGausHist();
796
797 h->SetName (Form("%s%s%s","H",fHistName.Data(),"LoGainArea"));
798 h->SetXTitle(fHistXTitle.Data());
799 h->SetYTitle(fHistYTitle.Data());
800
801 if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
802 {
803 h->SetTitle(Form("%s%s%s%s",fHistTitle.Data()," averaged on event-by-event basis ",
804 j==0 ? "Inner Pixels " : "Outer Pixels ","Low Gain Runs: "));
805 pix.InitBins();
806 pix.SetEventFrequency(fPulserFrequency);
807 }
808 else
809 {
810 h->SetTitle(Form("%s%s",fHistTitle.Data()," averaged on event-by-event basis Low Gain Area Idx "));
811 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
812 }
813 }
814 }
815
816 if (fAverageLoGainSectors->GetEntries()==0)
817 {
818 fAverageLoGainSectors->Expand(nsectors);
819
820 for (Int_t j=0; j<nsectors; j++)
821 {
822 (*fAverageLoGainSectors)[j] = new MHCalibrationPix(Form("%s%s",fHistName.Data(),"LoGainSector"),
823 Form("%s%s",fHistTitle.Data()," Low Gain Sector "));
824 MHCalibrationPix &pix = GetAverageLoGainSector(j);
825
826 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
827 pix.SetFirst(fFirst);
828 pix.SetLast (fLast);
829
830 TH1F *h = pix.GetHGausHist();
831
832 h->SetName (Form("%s%s%s","H",fHistName.Data(),"LoGainSector"));
833 h->SetTitle(Form("%s%s",fHistTitle.Data()," Low Gain Sector "));
834 h->SetXTitle(fHistXTitle.Data());
835 h->SetYTitle(fHistYTitle.Data());
836
837 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
838 }
839 }
840}
841
842//--------------------------------------------------------------------------------
843//
844// Retrieves from MGeomCam:
845// - number of pixels
846// - number of pixel areas
847// - number of sectors
848//
849// Return kFALSE, if sizes of the TObjArrays do not match npixels, nareas or nsectors
850//
851// Call FillHists()
852//
853Bool_t MHCalibrationCam::Fill(const MParContainer *par, const Stat_t w)
854{
855
856 if (!IsSizeCheck())
857 return FillHists(par,w);
858
859 const Int_t npixels = fGeom->GetNumPixels();
860 const Int_t nareas = fGeom->GetNumAreas();
861 const Int_t nsectors = fGeom->GetNumSectors();
862
863 //
864 // Hi-Gain ObjArrays
865 //
866 if (fHiGainArray->GetEntries() != npixels)
867 {
868 *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
869 return kFALSE;
870 }
871
872 if (IsLoGain())
873 {
874 if (fLoGainArray->GetEntries() != npixels)
875 {
876 *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
877 return kFALSE;
878 }
879 }
880
881 if (!IsAverageing())
882 return FillHists(par,w);
883
884 if (fAverageHiGainAreas->GetEntries() != nareas)
885 {
886 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
887 return kFALSE;
888 }
889
890 if (fAverageHiGainSectors->GetEntries() != nsectors)
891 {
892 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
893 return kFALSE;
894 }
895
896 if (IsLoGain())
897 {
898
899 if (fAverageLoGainAreas->GetEntries() != nareas)
900 {
901 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
902 return kFALSE;
903 }
904
905 if (fAverageLoGainSectors->GetEntries() != nsectors)
906 {
907 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
908 return kFALSE;
909 }
910 }
911
912 return FillHists(par, w);
913}
914
915// --------------------------------------------------------------------------
916//
917// 0) Ask if fHiGainArray and fLoGainArray have been initialized,
918// otherwise return kFALSE.
919// 1) FinalizeHists()
920// 2) FinalizeBadPixels()
921// 3) CalcAverageSigma()
922//
923Bool_t MHCalibrationCam::Finalize()
924{
925
926 if (fHiGainArray->GetEntries() == 0 && fLoGainArray->GetEntries() == 0)
927 {
928 *fLog << err << GetDescriptor()
929 << ": ERROR: Both (HiGain and LoGain) histogram arrays have not been initialized... abort." << endl;
930 return kFALSE;
931 }
932
933 if (!FinalizeHists())
934 return kFALSE;
935
936 FinalizeBadPixels();
937 CalcAverageSigma();
938
939 return kTRUE;
940}
941
942// -------------------------------------------------------------
943//
944// If MBadPixelsPix::IsBad():
945// - calls MHCalibrationPix::SetExcluded()
946//
947// Calls:
948// - MHGausEvents::InitBins()
949// - MHCalibrationPix::ChangeHistId(i)
950// - MHCalibrationPix::SetEventFrequency(fPulserFrequency)
951//
952void MHCalibrationCam::InitHists(MHCalibrationPix &hist, MBadPixelsPix &bad, const Int_t i)
953{
954
955 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
956 hist.SetExcluded();
957
958 hist.InitBins();
959 hist.ChangeHistId(i);
960 hist.SetEventFrequency(fPulserFrequency);
961
962 TH1F *h = hist.GetHGausHist();
963 h->SetTitle( Form("%s%s", h->GetTitle()," Runs: "));
964}
965
966// --------------------------------------------------------------------------
967//
968// Calls FitHiGainHists for every entry in:
969// - fHiGainArray
970// - fAverageHiGainAreas
971// - fAverageHiGainSectors
972//
973void MHCalibrationCam::FitHiGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
974 MBadPixelsPix::UncalibratedType_t fittyp,
975 MBadPixelsPix::UncalibratedType_t osctyp)
976{
977
978 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
979 {
980
981 MHCalibrationPix &hist = (*this)[i];
982
983 if (hist.IsExcluded())
984 continue;
985
986 MCalibrationPix &pix = calcam[i];
987 MBadPixelsPix &bad = badcam[i];
988
989 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
990
991 }
992
993 if (!IsAverageing())
994 return;
995
996 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
997 {
998
999 MHCalibrationPix &hist = GetAverageHiGainArea(j);
1000 MCalibrationPix &pix = calcam.GetAverageArea(j);
1001 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
1002
1003 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
1004 }
1005
1006
1007 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
1008 {
1009
1010 MHCalibrationPix &hist = GetAverageHiGainSector(j);
1011 MCalibrationPix &pix = calcam.GetAverageSector(j);
1012 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
1013
1014 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
1015 }
1016
1017}
1018
1019// --------------------------------------------------------------------------
1020//
1021// Calls FitLoGainHists for every entry in:
1022// - fLoGainArray
1023// - fAverageLoGainAreas
1024// - fAverageLoGainSectors
1025//
1026void MHCalibrationCam::FitLoGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
1027 MBadPixelsPix::UncalibratedType_t fittyp,
1028 MBadPixelsPix::UncalibratedType_t osctyp)
1029{
1030
1031 if (!IsLoGain())
1032 return;
1033
1034 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
1035 {
1036
1037 MHCalibrationPix &hist = (*this)(i);
1038
1039 if (hist.IsExcluded())
1040 continue;
1041
1042 MCalibrationPix &pix = calcam[i];
1043 MBadPixelsPix &bad = badcam[i];
1044
1045 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1046
1047 }
1048
1049 if (!IsAverageing())
1050 return;
1051
1052 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
1053 {
1054
1055 MHCalibrationPix &hist = GetAverageLoGainArea(j);
1056 MCalibrationPix &pix = calcam.GetAverageArea(j);
1057 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
1058
1059 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1060 }
1061
1062
1063 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
1064 {
1065
1066 MHCalibrationPix &hist = GetAverageLoGainSector(j);
1067 MCalibrationPix &pix = calcam.GetAverageSector(j);
1068 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
1069
1070 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1071 }
1072}
1073
1074//------------------------------------------------------------
1075//
1076// For all averaged areas, the fitted sigma is multiplied with the square root of
1077// the number involved pixels
1078//
1079void MHCalibrationCam::CalcAverageSigma()
1080{
1081
1082 if (!fCam)
1083 return;
1084
1085 if (!IsAverageing())
1086 return;
1087
1088 for (UInt_t j=0; j<fGeom->GetNumAreas(); j++)
1089 {
1090
1091 MCalibrationPix &pix = fCam->GetAverageArea(j);
1092
1093 const Float_t numsqr = TMath::Sqrt((Float_t)fAverageAreaNum[j]);
1094 fAverageAreaSigma[j] = pix.GetSigma () * numsqr;
1095 fAverageAreaSigmaVar[j] = pix.GetSigmaErr () * pix.GetSigmaErr() * numsqr;
1096
1097 pix.SetSigma (fAverageAreaSigma[j]);
1098 pix.SetSigmaVar(fAverageAreaSigmaVar[j]);
1099
1100 fAverageAreaRelSigma [j] = fAverageAreaSigma[j] / pix.GetMean();
1101 fAverageAreaRelSigmaVar[j] = fAverageAreaSigmaVar[j] / (fAverageAreaSigma[j]*fAverageAreaSigma[j]);
1102 fAverageAreaRelSigmaVar[j] += pix.GetMeanRelVar();
1103 fAverageAreaRelSigmaVar[j] *= fAverageAreaRelSigma[j];
1104 }
1105}
1106
1107// ---------------------------------------------------------------------------
1108//
1109// Returns if the histogram is empty and sets the following flag:
1110// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
1111//
1112// Fits the histograms with a Gaussian, in case of failure
1113// calls MHCalibrationPix::RepeatFit(), in case of repeated failure
1114// calls MHCalibrationPix::BypassFit() and sets the following flags:
1115// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
1116// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1117//
1118// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
1119// In case no, sets the following flags:
1120// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
1121// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1122//
1123// Retrieves the results and store them in MCalibrationPix
1124//
1125void MHCalibrationCam::FitHiGainHists(MHCalibrationPix &hist,
1126 MCalibrationPix &pix,
1127 MBadPixelsPix &bad,
1128 MBadPixelsPix::UncalibratedType_t fittyp,
1129 MBadPixelsPix::UncalibratedType_t osctyp)
1130{
1131
1132
1133 if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
1134 return;
1135
1136 //
1137 // 2) Fit the Hi Gain histograms with a Gaussian
1138 //
1139 if (!hist.FitGaus())
1140 //
1141 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
1142 //
1143 if (!hist.RepeatFit())
1144 {
1145 hist.BypassFit();
1146 bad.SetUncalibrated( fittyp );
1147 }
1148
1149 hist.Renorm();
1150 //
1151 // 4) Check for oscillations
1152 //
1153 if (IsOscillations())
1154 {
1155 hist.CreateFourierSpectrum();
1156
1157 if (!hist.IsFourierSpectrumOK())
1158 bad.SetUncalibrated( osctyp );
1159 }
1160
1161 //
1162 // 5) Retrieve the results and store them in this class
1163 //
1164 pix.SetHiGainMean ( hist.GetMean() );
1165 pix.SetHiGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
1166 pix.SetHiGainSigma ( hist.GetSigma() );
1167 pix.SetHiGainSigmaVar ( hist.GetSigmaErr()* hist.GetSigmaErr() );
1168 pix.SetHiGainProb ( hist.GetProb() );
1169 pix.SetHiGainNumBlackout( hist.GetBlackout() );
1170 pix.SetHiGainNumPickup ( hist.GetPickup() );
1171
1172 if (IsDebug())
1173 {
1174 *fLog << dbginf << GetDescriptor() << ": ID " << hist.GetPixId()
1175 << " HiGainSaturation: " << pix.IsHiGainSaturation()
1176 << " HiGainMean: " << hist.GetMean ()
1177 << " HiGainMeanErr: " << hist.GetMeanErr ()
1178 << " HiGainMeanSigma: " << hist.GetSigma ()
1179 << " HiGainMeanSigmaErr: " << hist.GetSigmaErr()
1180 << " HiGainMeanProb: " << hist.GetProb ()
1181 << " HiGainNumBlackout: " << hist.GetBlackout()
1182 << " HiGainNumPickup : " << hist.GetPickup ()
1183 << endl;
1184 }
1185
1186}
1187
1188
1189// ---------------------------------------------------------------------------
1190//
1191// Returns if the histogram is empty and sets the following flag:
1192// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
1193//
1194// Fits the histograms with a Gaussian, in case of failure
1195// calls MHCalibrationPix::RepeatFit(), in case of repeated failure
1196// calls MHCalibrationPix::BypassFit() and sets the following flags:
1197// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
1198// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1199//
1200// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
1201// In case no, sets the following flags:
1202// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
1203// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1204//
1205// Retrieves the results and store them in MCalibrationPix
1206//
1207void MHCalibrationCam::FitLoGainHists(MHCalibrationPix &hist,
1208 MCalibrationPix &pix,
1209 MBadPixelsPix &bad,
1210 MBadPixelsPix::UncalibratedType_t fittyp,
1211 MBadPixelsPix::UncalibratedType_t osctyp)
1212{
1213
1214 if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
1215 return;
1216
1217
1218 //
1219 // 2) Fit the Hi Gain histograms with a Gaussian
1220 //
1221 if (!hist.FitGaus())
1222 //
1223 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
1224 //
1225 if (!hist.RepeatFit())
1226 {
1227 hist.BypassFit();
1228 bad.SetUncalibrated( fittyp );
1229 }
1230
1231 //
1232 // 4) Check for oscillations
1233 //
1234 if (IsOscillations())
1235 {
1236 hist.CreateFourierSpectrum();
1237
1238 if (!hist.IsFourierSpectrumOK())
1239 bad.SetUncalibrated( osctyp );
1240 }
1241
1242 //
1243 // 5) Retrieve the results and store them in this class
1244 //
1245 pix.SetLoGainMean ( hist.GetMean() );
1246 pix.SetLoGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
1247 pix.SetLoGainSigma ( hist.GetSigma() );
1248 pix.SetLoGainSigmaVar ( hist.GetSigmaErr() * hist.GetSigmaErr() );
1249 pix.SetLoGainProb ( hist.GetProb() );
1250 pix.SetLoGainNumBlackout( hist.GetBlackout() );
1251 pix.SetLoGainNumPickup ( hist.GetPickup() );
1252
1253 if (IsDebug())
1254 {
1255 *fLog << dbginf << GetDescriptor() << "ID: " << hist.GetPixId()
1256 << " HiGainSaturation: " << pix.IsHiGainSaturation()
1257 << " LoGainMean: " << hist.GetMean ()
1258 << " LoGainMeanErr: " << hist.GetMeanErr ()
1259 << " LoGainMeanSigma: " << hist.GetSigma ()
1260 << " LoGainMeanSigmaErr: " << hist.GetSigmaErr()
1261 << " LoGainMeanProb: " << hist.GetProb ()
1262 << " LoGainNumBlackout: " << hist.GetBlackout()
1263 << " LoGainNumPickup : " << hist.GetPickup ()
1264 << endl;
1265 }
1266
1267}
1268
1269
1270
1271// -----------------------------------------------------------------------------
1272//
1273// Default draw:
1274//
1275// Displays the averaged areas, both High Gain and Low Gain
1276//
1277// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1278//
1279void MHCalibrationCam::Draw(const Option_t *opt)
1280{
1281
1282 if (!IsAverageing())
1283 return;
1284
1285 const Int_t nareas = fAverageHiGainAreas->GetEntries();
1286 if (nareas == 0)
1287 return;
1288
1289 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
1290 pad->SetBorderMode(0);
1291
1292 pad->Divide(IsLoGain() ? 2 : 1,nareas);
1293
1294 for (Int_t i=0; i<nareas;i++)
1295 {
1296
1297 pad->cd(IsLoGain() ? 2*(i+1)-1 : i+1);
1298 GetAverageHiGainArea(i).Draw(opt);
1299
1300 if (!fAverageAreaSat[i])
1301 DrawAverageSigma(fAverageAreaSat[i], i,
1302 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
1303 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
1304
1305 if (IsLoGain())
1306 {
1307 pad->cd(2*(i+1));
1308 GetAverageLoGainArea(i).Draw(opt);
1309 }
1310
1311 if (fAverageAreaSat[i])
1312 DrawAverageSigma(fAverageAreaSat[i], i,
1313 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
1314 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
1315 }
1316
1317}
1318
1319// -----------------------------------------------------------------------------
1320//
1321// Default draw:
1322//
1323// Displays a TPaveText with the re-normalized sigmas of the average area
1324//
1325void MHCalibrationCam::DrawAverageSigma(Bool_t sat, Bool_t inner,
1326 Float_t sigma, Float_t sigmavar,
1327 Float_t relsigma, Float_t relsigmavar) const
1328{
1329
1330 if (sigma != 0 && sigmavar >= 0 && relsigmavar >= 0.)
1331 {
1332
1333 TPad *newpad = new TPad("newpad","transparent",0,0,1,1);
1334 newpad->SetFillStyle(4000);
1335 newpad->Draw();
1336 newpad->cd();
1337
1338 TPaveText *text = new TPaveText(sat? 0.1 : 0.35,0.7,sat ? 0.4 : 0.7,1.0);
1339 text->SetTextSize(0.07);
1340 const TString line1 = Form("%s%s%s",inner ? "Outer" : "Inner",
1341 " Pixels ", sat ? "Low Gain" : "High Gain");
1342 TText *txt1 = text->AddText(line1.Data());
1343 const TString line2 = Form("#sigma per pix: %2.2f #pm %2.2f",sigma,TMath::Sqrt(sigmavar));
1344 TText *txt2 = text->AddText(line2.Data());
1345 const TString line3 = Form("Rel. #sigma per pix: %2.2f #pm %2.2f",relsigma,TMath::Sqrt(relsigmavar));
1346 TText *txt3 = text->AddText(line3.Data());
1347 text->Draw("");
1348
1349 text->SetBit(kCanDelete);
1350 txt1->SetBit(kCanDelete);
1351 txt2->SetBit(kCanDelete);
1352 txt3->SetBit(kCanDelete);
1353 newpad->SetBit(kCanDelete);
1354 }
1355}
1356
1357Int_t MHCalibrationCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1358{
1359
1360 Bool_t rc = kFALSE;
1361 if (IsEnvDefined(env, prefix, "Debug", print))
1362 {
1363 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
1364 rc = kTRUE;
1365 }
1366 if (IsEnvDefined(env, prefix, "LoGain", print))
1367 {
1368 SetDebug(GetEnvValue(env, prefix, "LoGain", IsLoGain()));
1369 rc = kTRUE;
1370 }
1371 if (IsEnvDefined(env, prefix, "Oscillations", print))
1372 {
1373 SetOscillations(GetEnvValue(env, prefix, "Oscillations", IsOscillations()));
1374 rc = kTRUE;
1375 }
1376 if (IsEnvDefined(env, prefix, "SizeCheck", print))
1377 {
1378 SetSizeCheck(GetEnvValue(env, prefix, "SizeCheck", IsSizeCheck()));
1379 rc = kTRUE;
1380 }
1381 if (IsEnvDefined(env, prefix, "Averageing", print))
1382 {
1383 SetAverageing(GetEnvValue(env, prefix, "Averageing", IsAverageing()));
1384 rc = kTRUE;
1385 }
1386
1387 if (IsEnvDefined(env, prefix, "Nbins", print))
1388 {
1389 SetNbins(GetEnvValue(env, prefix, "Nbins", fNbins));
1390 rc = kTRUE;
1391 }
1392 if (IsEnvDefined(env, prefix, "First", print))
1393 {
1394 SetFirst(GetEnvValue(env, prefix, "First", fFirst));
1395 rc = kTRUE;
1396 }
1397 if (IsEnvDefined(env, prefix, "Last", print))
1398 {
1399 SetLast(GetEnvValue(env, prefix, "Last", fLast));
1400 rc = kTRUE;
1401 }
1402
1403 return rc;
1404}
Note: See TracBrowser for help on using the repository browser.