source: tags/Mars-V0.9.4.2/mhcalib/MHCalibrationCam.cc

Last change on this file was 7288, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 49.9 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 TOrdCollection'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 TOrdCollection'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// Class Version 5:
55// + Double_t fLowerFitLimitHiGain; // Lower limit for the fit range for the hi-gain hist
56// + Double_t fUpperFitLimitHiGain; // Upper limit for the fit range for the hi-gain hist
57// + Double_t fLowerFitLimitLoGain; // Lower limit for the fit range for the lo-gain hist
58// + Double_t fUpperFitLimitLoGain; // Upper limit for the fit range for the lo-gain hist
59// + Bool_t fIsHiGainFitRanges; // Are high-gain fit ranges defined?
60// + Bool_t fIsLoGainFitRanges; // Are low-gain fit ranges defined?
61//
62// Class Version 5:
63// + Int_t fMaxNumEvts; // Max Number of events
64//
65/////////////////////////////////////////////////////////////////////////////
66#include "MHCalibrationCam.h"
67#include "MHCalibrationPix.h"
68
69#include <TVirtualPad.h>
70#include <TCanvas.h>
71#include <TPad.h>
72#include <TText.h>
73#include <TPaveText.h>
74#include <TOrdCollection.h>
75#include <TROOT.h>
76
77#include "MLog.h"
78#include "MLogManip.h"
79
80#include "MCalibrationIntensityCam.h"
81#include "MCalibrationCam.h"
82#include "MCalibrationPix.h"
83
84#include "MBadPixelsIntensityCam.h"
85#include "MBadPixelsCam.h"
86#include "MBadPixelsPix.h"
87
88#include "MGeomCam.h"
89#include "MGeomPix.h"
90
91#include "MParList.h"
92
93#include "MRawRunHeader.h"
94
95ClassImp(MHCalibrationCam);
96
97using namespace std;
98
99const Double_t MHCalibrationCam::fgLowerFitLimitHiGain = 0;
100const Double_t MHCalibrationCam::fgUpperFitLimitHiGain = 0;
101const Double_t MHCalibrationCam::fgLowerFitLimitLoGain = 0;
102const Double_t MHCalibrationCam::fgUpperFitLimitLoGain = 0;
103
104const Int_t MHCalibrationCam::fgPulserFrequency = 500;
105const Float_t MHCalibrationCam::fgProbLimit = 0.0001;
106const Float_t MHCalibrationCam::fgOverflowLimit = 0.005;
107const Int_t MHCalibrationCam::fgMaxNumEvts = 4096;
108
109const TString MHCalibrationCam::gsHistName = "Hist";
110const TString MHCalibrationCam::gsHistTitle = "";
111const TString MHCalibrationCam::gsHistXTitle = "";
112const TString MHCalibrationCam::gsHistYTitle = "Nr. events";
113
114// --------------------------------------------------------------------------
115//
116// Default Constructor.
117//
118// Sets:
119// - all pointers to NULL
120//
121// Initializes and sets owner of:
122// - fHiGainArray, fLoGainArray
123// - fAverageHiGainAreas, fAverageLoGainAreas
124// - fAverageHiGainSectors, fAverageLoGainSectors
125//
126// Initializes:
127// - fPulserFrequency to fgPulserFrequency
128// - fProbLimit to fgProbLimit
129// - fOverflowLimit to fgOverflowLimit
130//
131// - SetAveregeing (kTRUE);
132// - SetDebug (kFALSE);
133// - SetLoGain (kTRUE);
134//- SetOscillations(kTRUE);
135//- SetSizeCheck (kTRUE);
136//- SetInterlaced (kFALSE);
137//- SetLowerFitLimitHiGain();
138//- SetUpperFitLimitHiGain();
139//- SetLowerFitLimitLoGain();
140//- SetUpperFitLimitLoGain();
141//
142MHCalibrationCam::MHCalibrationCam(const char *name, const char *title)
143 : fIsHiGainFitRanges(kFALSE), fIsLoGainFitRanges(kFALSE),
144 fHistName(gsHistName),fHistTitle(gsHistTitle),
145 fHistXTitle(gsHistXTitle),fHistYTitle(gsHistYTitle),
146 fCurrentNumEvts(0),
147 fColor(MCalibrationCam::kNONE), fIntensBad(NULL),
148 fBadPixels(NULL), fIntensCam(NULL), fCam(NULL), fGeom(NULL),
149 fRunHeader(NULL)
150{
151
152 fHiGainArray = new TOrdCollection;
153 fHiGainArray->SetOwner();
154
155 fLoGainArray = new TOrdCollection;
156 fLoGainArray->SetOwner();
157
158 fAverageHiGainAreas = new TOrdCollection;
159 fAverageHiGainAreas->SetOwner();
160
161 fAverageLoGainAreas = new TOrdCollection;
162 fAverageLoGainAreas->SetOwner();
163
164 fAverageHiGainSectors = new TOrdCollection;
165 fAverageHiGainSectors->SetOwner();
166
167 fAverageLoGainSectors = new TOrdCollection;
168 fAverageLoGainSectors->SetOwner();
169
170 SetPulserFrequency();
171 SetProbLimit();
172 SetOverflowLimit();
173 SetMaxNumEvts();
174
175 SetAverageing (kTRUE);
176 SetDebug (kFALSE);
177 SetLoGain (kTRUE);
178 SetOscillations(kTRUE);
179 SetSizeCheck (kTRUE);
180 SetIsReset (kTRUE);
181
182 SetLowerFitLimitHiGain();
183 SetUpperFitLimitHiGain();
184 SetLowerFitLimitLoGain();
185 SetUpperFitLimitLoGain();
186}
187
188// --------------------------------------------------------------------------
189//
190// Deletes the TOrdCollection of:
191// - fHiGainArray, fLoGainArray
192// - fAverageHiGainAreas, fAverageLoGainAreas
193// - fAverageHiGainSectors, fAverageLoGainSectors
194//
195MHCalibrationCam::~MHCalibrationCam()
196{
197
198 delete fHiGainArray;
199 delete fLoGainArray;
200
201 delete fAverageHiGainAreas;
202 delete fAverageLoGainAreas;
203
204 delete fAverageHiGainSectors;
205 delete fAverageLoGainSectors;
206 /*
207
208 Remove ( fHiGainArray );
209 Remove ( fLoGainArray );
210
211 Remove ( fAverageHiGainAreas );
212 Remove ( fAverageLoGainAreas );
213
214 Remove ( fAverageHiGainSectors );
215 Remove ( fAverageLoGainSectors );
216 */
217
218}
219
220void MHCalibrationCam::Remove(TOrdCollection *col)
221{
222
223 if (!col)
224 return;
225
226 TOrdCollectionIter Next(col);
227
228 Int_t count = 0;
229
230 while (MHCalibrationPix *obj = (MHCalibrationPix*)Next())
231 {
232 *fLog << ++count << " " << obj << flush;
233 if (obj && obj->IsOnHeap())
234 {
235 obj->Draw();
236 delete obj;
237 }
238 }
239
240 delete col;
241}
242
243// --------------------------------------------------------------------------
244//
245// Returns size of fHiGainArray
246//
247const Int_t MHCalibrationCam::GetSize() const
248{
249 return fHiGainArray->GetSize();
250}
251
252// --------------------------------------------------------------------------
253//
254// Get i-th High Gain pixel (pixel number)
255//
256MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i)
257{
258 return *static_cast<MHCalibrationPix*>(fHiGainArray->At(i));
259}
260
261// --------------------------------------------------------------------------
262//
263// Get i-th High Gain pixel (pixel number)
264//
265const MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i) const
266{
267 return *static_cast<MHCalibrationPix*>(fHiGainArray->At(i));
268}
269
270// --------------------------------------------------------------------------
271//
272// Get i-th Low Gain pixel (pixel number)
273//
274MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i)
275{
276 return *static_cast<MHCalibrationPix*>(fLoGainArray->At(i));
277}
278
279// --------------------------------------------------------------------------
280//
281// Get i-th Low Gain pixel (pixel number)
282//
283const MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i) const
284{
285 return *static_cast<MHCalibrationPix*>(fLoGainArray->At(i));
286}
287
288// --------------------------------------------------------------------------
289//
290// Returns the current size of the TOrdCollection fAverageHiGainAreas
291// independently if the MHCalibrationPix is filled with values or not.
292//
293const Int_t MHCalibrationCam::GetAverageAreas() const
294{
295 return fAverageHiGainAreas->GetSize();
296}
297
298// --------------------------------------------------------------------------
299//
300// Get i-th High Gain pixel Area (area number)
301//
302MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i)
303{
304 return *static_cast<MHCalibrationPix*>(fAverageHiGainAreas->At(i));
305}
306
307// --------------------------------------------------------------------------
308//
309// Get i-th High Gain pixel Area (area number)
310//
311const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i) const
312{
313 return *static_cast<MHCalibrationPix *>(fAverageHiGainAreas->At(i));
314}
315
316// --------------------------------------------------------------------------
317//
318// Get i-th Low Gain pixel Area (area number)
319//
320MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i)
321{
322 return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->At(i));
323}
324
325// --------------------------------------------------------------------------
326//
327// Get i-th Low Gain pixel Area (area number)
328//
329const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i) const
330{
331 return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->At(i));
332}
333
334// --------------------------------------------------------------------------
335//
336// Returns the current size of the TOrdCollection fAverageHiGainSectors
337// independently if the MHCalibrationPix is filled with values or not.
338//
339const Int_t MHCalibrationCam::GetAverageSectors() const
340{
341 return fAverageHiGainSectors->GetSize();
342}
343
344// --------------------------------------------------------------------------
345//
346// Get i-th High Gain Sector (sector number)
347//
348MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i)
349{
350 return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->At(i));
351}
352
353// --------------------------------------------------------------------------
354//
355// Get i-th High Gain Sector (sector number)
356//
357const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i) const
358{
359 return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->At(i));
360}
361
362// --------------------------------------------------------------------------
363//
364// Get i-th Low Gain Sector (sector number)
365//
366MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i)
367{
368 return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->At(i));
369}
370
371// --------------------------------------------------------------------------
372//
373// Get i-th Low Gain Sector (sector number)
374//
375const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i) const
376{
377 return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->At(i));
378}
379
380// --------------------------------------------------------------------------
381//
382// Calls ResetHistTitles()
383//
384// Calls Reset() for each entry in:
385// - fHiGainArray, fLoGainArray
386// - fAverageHiGainAreas, fAverageLoGainAreas
387// - fAverageHiGainSectors, fAverageLoGainSectors
388//
389void MHCalibrationCam::ResetHists()
390{
391 SetIsReset();
392
393 fCurrentNumEvts = 0;
394
395 ResetHistTitles();
396
397 if (fHiGainArray)
398 { fHiGainArray->ForEach(MHCalibrationPix,Reset)(); }
399
400 if (IsAverageing())
401 {
402 if (fAverageHiGainAreas)
403 { fAverageHiGainAreas->ForEach(MHCalibrationPix,Reset)(); }
404 if (fAverageHiGainSectors)
405 { fAverageHiGainSectors->ForEach(MHCalibrationPix,Reset)(); }
406 }
407
408 if (!IsLoGain())
409 return;
410
411 if (fLoGainArray)
412 { fLoGainArray->ForEach(MHCalibrationPix,Reset)(); }
413 if (IsAverageing())
414 {
415 if (fAverageLoGainAreas)
416 { fAverageLoGainAreas->ForEach(MHCalibrationPix,Reset)(); }
417 if (fAverageLoGainSectors)
418 { fAverageLoGainSectors->ForEach(MHCalibrationPix,Reset)(); }
419 }
420}
421
422// --------------------------------------------------------------------------
423//
424// Resets the histogram titles for each entry in:
425// - fHiGainArray, fLoGainArray
426// - fAverageHiGainAreas, fAverageLoGainAreas
427// - fAverageHiGainSectors, fAverageLoGainSectors
428//
429void MHCalibrationCam::ResetHistTitles()
430{
431
432 TH1F *h;
433
434 if (fHiGainArray)
435 for (Int_t i=0;i<fHiGainArray->GetSize(); i++)
436 {
437 MHCalibrationPix &pix = (*this)[i];
438 h = pix.GetHGausHist();
439 h->SetName (Form("H%sHiGainPix%04d",fHistName.Data(),i));
440 h->SetTitle(Form("%s High Gain Pixel %04d Runs: ",fHistTitle.Data(),i));
441 h->SetXTitle(fHistXTitle.Data());
442 h->SetYTitle(fHistYTitle.Data());
443 }
444
445 if (IsAverageing())
446 {
447 if (fAverageHiGainAreas)
448 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
449 {
450 MHCalibrationPix &pix = GetAverageHiGainArea(j);
451 h = pix.GetHGausHist();
452 h->SetName (Form("H%sHiGainArea%d",fHistName.Data(),j));
453 h->SetXTitle(fHistXTitle.Data());
454 h->SetYTitle(fHistYTitle.Data());
455 if (fGeom->InheritsFrom("MGeomCamMagic"))
456 h->SetTitle(Form("%s averaged on event-by-event basis %s High Gain Runs: ",
457 fHistTitle.Data(), j==0 ? "Inner Pixels" : "Outer Pixels"));
458 else
459 h->SetTitle(Form("%s averaged on event-by-event basis High Gain Area Idx %d Runs: ",
460 fHistTitle.Data(), j));
461 }
462
463 if (fAverageHiGainSectors)
464 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
465 {
466 MHCalibrationPix &pix = GetAverageHiGainSector(j);
467 h = pix.GetHGausHist();
468 h->SetName (Form("H%sHiGainSector%02d",fHistName.Data(),j));
469 h->SetTitle(Form("%s averaged on event-by-event basis High Gain Sector %02d Runs: ",
470 fHistTitle.Data(), j));
471 h->SetXTitle(fHistXTitle.Data());
472 h->SetYTitle(fHistYTitle.Data());
473 }
474 }
475
476 if (!IsLoGain())
477 return;
478
479 if (fLoGainArray)
480 for (Int_t i=0;i<fLoGainArray->GetSize(); i++)
481 {
482 MHCalibrationPix &pix = (*this)(i);
483 h = pix.GetHGausHist();
484 h->SetName (Form("H%sLoGainPix%04d",fHistName.Data(),i));
485 h->SetTitle(Form("%s Low Gain Pixel %04d Runs: ", fHistTitle.Data(),i));
486 h->SetXTitle(fHistXTitle.Data());
487 h->SetYTitle(fHistYTitle.Data());
488 }
489
490 if (IsAverageing())
491 {
492 if (fAverageLoGainAreas)
493 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
494 {
495 MHCalibrationPix &pix = GetAverageLoGainArea(j);
496 h = pix.GetHGausHist();
497 h->SetName (Form("H%sLoGainArea%d", fHistName.Data(), j));
498 h->SetXTitle(fHistXTitle.Data());
499 h->SetYTitle(fHistYTitle.Data());
500 if (fGeom->InheritsFrom("MGeomCamMagic"))
501 h->SetTitle(Form("%s averaged on event-by-event basis %s Low Gain Runs: ",
502 fHistTitle.Data(), j==0 ? "Inner Pixels" : "Outer Pixels"));
503 else
504 h->SetTitle(Form("%s averaged on event-by-event basis Low Gain Area Idx %d Runs: ",
505 fHistTitle.Data(), j));
506 }
507
508 if (fAverageLoGainSectors)
509 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
510 {
511 MHCalibrationPix &pix = GetAverageLoGainSector(j);
512 h = pix.GetHGausHist();
513 h->SetName (Form("H%sLoGainSector%02d",fHistName.Data(),j));
514 h->SetTitle(Form("%s averaged on event-by-event basis Low Gain Sector %02d Runs: ",
515 fHistTitle.Data(), j));
516 h->SetXTitle(fHistXTitle.Data());
517 h->SetYTitle(fHistYTitle.Data());
518 }
519 }
520}
521
522// --------------------------------------------------------------------------
523//
524// Gets the pointers to:
525// - MGeomCam
526//
527// Calls SetupHists(const MParList *pList)
528//
529// Calls Delete-Function of:
530// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
531// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
532// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
533//
534Bool_t MHCalibrationCam::SetupFill(const MParList *pList)
535{
536
537 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
538 if (!fGeom)
539 {
540 *fLog << err << GetDescriptor()
541 << ": MGeomCam not found... aborting." << endl;
542 return kFALSE;
543 }
544
545 fRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
546 if (!fRunHeader)
547 {
548 *fLog << warn << GetDescriptor()
549 << ": MRawRunHeader not found... will not store run numbers." << endl;
550 }
551
552 fCurrentNumEvts = 0;
553
554 return SetupHists(pList);
555}
556
557
558// --------------------------------------------------------------------------
559//
560// Searches MRawEvtHeader to find the correct pulser colour
561//
562// Gets or creates the pointers to:
563// - MBadPixelsIntensityCam
564// - MBadPixelsCam
565//
566// Searches pointer to:
567// - MArrivalTimeCam
568//
569// Initializes, if empty to MArrivalTimeCam::GetSize() for:
570// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
571//
572// Initializes, if empty to MGeomCam::GetNumAreas() for:
573// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
574//
575// Initializes, if empty to MGeomCam::GetNumSectors() for:
576// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
577//
578// Initializes TArray's to MGeomCam::GetNumAreas and MGeomCam::GetNumSectors, respectively
579// Fills with number of valid pixels (if !MBadPixelsPix::IsBad()):
580// - MHCalibrationCam::fAverageAreaNum[area index]
581// - MHCalibrationCam::fAverageSectorNum[area index]
582//
583// Calls InitializeHists() for every entry in:
584// - MHCalibrationCam::fHiGainArray
585// - MHCalibrationCam::fAverageHiGainAreas
586// - MHCalibrationCam::fAverageHiGainSectors
587//
588// Sets Titles and Names for the Histograms
589// - MHCalibrationCam::fAverageHiGainAreas
590// - MHCalibrationCam::fAverageHiGainSectors
591//
592// Retrieves the run numbers from MRawRunHeader and stores them in fRunNumbers
593//
594Bool_t MHCalibrationCam::ReInit(MParList *pList)
595{
596
597 const Int_t npixels = fGeom->GetNumPixels();
598 const Int_t nsectors = fGeom->GetNumSectors();
599 const Int_t nareas = fGeom->GetNumAreas();
600
601 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
602 if (fIntensBad)
603 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
604 else
605 {
606 fBadPixels = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
607 if (!fBadPixels)
608 {
609
610 fBadPixels = (MBadPixelsCam*)pList->FindCreateObj(AddSerialNumber("MBadPixelsCam"));
611 if (!fBadPixels)
612 return kFALSE;
613
614 fBadPixels->InitSize(npixels);
615 }
616 }
617
618 if (IsAverageing())
619 {
620 //
621 // The function TArrayF::Set() already sets all entries to 0.
622 //
623 fAverageAreaNum. Set(nareas);
624 fAverageAreaSat. Set(nareas);
625 fAverageAreaSigma. Set(nareas);
626 fAverageAreaSigmaVar. Set(nareas);
627 fAverageAreaRelSigma. Set(nareas);
628 fAverageAreaRelSigmaVar.Set(nareas);
629 fAverageSectorNum. Set(nsectors);
630
631 for (Int_t aidx=0; aidx<nareas; aidx++)
632 fAverageAreaNum[aidx] = 0;
633
634 for (Int_t sector=0; sector<nsectors; sector++)
635 fAverageSectorNum[sector] = 0;
636
637 for (Int_t i=0; i<npixels; i++)
638 {
639
640 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
641 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
642 continue;
643
644 fAverageAreaNum [(*fGeom)[i].GetAidx() ]++;
645 fAverageSectorNum[(*fGeom)[i].GetSector()]++;
646 }
647 }
648
649 //
650 // Because ReInit has been called, a new run number is added
651 //
652 fRunNumbers.Set(fRunNumbers.GetSize()+1);
653
654 if (fRunHeader)
655 {
656 fRunNumbers[fRunNumbers.GetSize()-1] = fRunHeader->GetRunNumber();
657 if (IsLoGain())
658 SetLoGain(fRunHeader->GetNumSamplesLoGain());
659 }
660
661 fCurrentNumEvts = 0;
662
663 if (!ReInitHists(pList))
664 return kFALSE;
665
666 ResetHistTitles();
667
668 if (!fRunHeader)
669 return kTRUE;
670
671 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
672 {
673 TH1F *h = (*this)[i].GetHGausHist();
674 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
675 }
676
677 if (IsLoGain())
678 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
679 {
680 TH1F *h = (*this)(i).GetHGausHist();
681 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
682 }
683
684 if (!IsAverageing())
685 return kTRUE;
686
687 for (Int_t j=0; j<nareas; j++)
688 {
689 TH1F *h = GetAverageHiGainArea(j).GetHGausHist();
690 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
691 }
692
693 if (IsLoGain())
694 for (Int_t j=0; j<nareas; j++)
695 {
696 TH1F *h = GetAverageLoGainArea(j).GetHGausHist();
697 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
698 }
699
700 for (Int_t j=0; j<nsectors; j++)
701 {
702 TH1F *h = GetAverageHiGainSector(j).GetHGausHist();
703 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
704 }
705
706 if (IsLoGain())
707 for (Int_t j=0; j<nsectors; j++)
708 {
709 TH1F *h = GetAverageLoGainSector(j).GetHGausHist();
710 h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
711 }
712
713 return kTRUE;
714}
715
716//--------------------------------------------------------------------------------------
717//
718// Initializes the High Gain Arrays:
719//
720// - For every entry in the expanded arrays:
721// * Initialize an MHCalibrationPix
722// * Set Binning from fNbins, fFirst and fLast
723// * Set Histgram names and titles from fHistName and fHistTitle
724// * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
725// * Call InitHists
726//
727void MHCalibrationCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
728{
729
730 if (fHiGainArray->GetSize()==0)
731 {
732 for (Int_t i=0; i<npixels; i++)
733 {
734 fHiGainArray->AddAt(new MHCalibrationPix(Form("%sHiGainPix%04d",fHistName.Data(),i),
735 Form("%s High Gain Pixel %4d",fHistTitle.Data(),i)),i);
736
737 MHCalibrationPix &pix = (*this)[i];
738 pix.SetNbins(fNbins);
739 pix.SetFirst(fFirst);
740 pix.SetLast (fLast);
741 pix.SetProbLimit(fProbLimit);
742
743 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
744 InitHists(pix,bad,i);
745
746 if (fCam)
747 (*fCam)[i].SetPixId(i);
748 }
749 }
750
751 if (!IsAverageing())
752 return;
753
754 if (fAverageHiGainAreas->GetSize()==0)
755 {
756 for (Int_t j=0; j<nareas; j++)
757 {
758 fAverageHiGainAreas->AddAt(new MHCalibrationPix(Form("%sHiGainArea%d",fHistName.Data(),j),
759 Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
760
761 MHCalibrationPix &pix = GetAverageHiGainArea(j);
762
763 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
764 pix.SetFirst(fFirst);
765 pix.SetLast (fLast);
766
767 if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
768 {
769 pix.InitBins();
770 pix.SetEventFrequency(fPulserFrequency);
771 }
772 else
773 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
774 }
775 }
776
777 if (fAverageHiGainSectors->GetSize()==0)
778 {
779 for (Int_t j=0; j<nsectors; j++)
780 {
781 fAverageHiGainSectors->AddAt(new MHCalibrationPix(Form("%sHiGainSector%02d",fHistName.Data(),j),
782 Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
783 MHCalibrationPix &pix = GetAverageHiGainSector(j);
784
785 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
786 pix.SetFirst(fFirst);
787 pix.SetLast (fLast);
788
789 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
790 }
791 }
792}
793
794//--------------------------------------------------------------------------------------
795//
796// Return, if IsLoGain() is kFALSE
797//
798// Initializes the Low Gain Arrays:
799//
800// - For every entry in the expanded arrays:
801// * Initialize an MHCalibrationPix
802// * Set Binning from fNbins, fFirst and fLast
803// * Set Histgram names and titles from fHistName and fHistTitle
804// * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
805// * Call InitHists
806//
807void MHCalibrationCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
808{
809
810 if (!IsLoGain())
811 return;
812
813 if (fLoGainArray->GetSize()==0)
814 {
815 for (Int_t i=0; i<npixels; i++)
816 {
817 fLoGainArray->AddAt(new MHCalibrationPix(Form("%sLoGainPix%04d",fHistName.Data(),i),
818 Form("%s Low Gain Pixel%04d",fHistTitle.Data(),i)),i);
819
820 MHCalibrationPix &pix = (*this)(i);
821 pix.SetNbins(fNbins);
822 pix.SetFirst(fFirst);
823 pix.SetLast (fLast);
824 pix.SetProbLimit(fProbLimit);
825
826 MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
827 InitHists(pix,bad,i);
828 }
829 }
830
831 if (!IsAverageing())
832 return;
833
834 if (fAverageLoGainAreas->GetSize()==0)
835 {
836 for (Int_t j=0; j<nareas; j++)
837 {
838 fAverageLoGainAreas->AddAt(new MHCalibrationPix(Form("%sLoGainArea%d",fHistName.Data(),j),
839 Form("%s Low Gain Area Idx %d",fHistTitle.Data(),j)),j);
840
841 MHCalibrationPix &pix = GetAverageLoGainArea(j);
842
843 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
844 pix.SetFirst(fFirst);
845 pix.SetLast (fLast);
846
847 if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
848 {
849 pix.InitBins();
850 pix.SetEventFrequency(fPulserFrequency);
851 }
852 else
853 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
854 }
855 }
856
857 if (fAverageLoGainSectors->GetSize()==0)
858 {
859 for (Int_t j=0; j<nsectors; j++)
860 {
861 fAverageLoGainSectors->AddAt(new MHCalibrationPix(Form("%sLoGainSector%02d",fHistName.Data(),j),
862 Form("%s Low Gain Sector %02d",fHistTitle.Data(),j)),j);
863 MHCalibrationPix &pix = GetAverageLoGainSector(j);
864
865 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
866 pix.SetFirst(fFirst);
867 pix.SetLast (fLast);
868
869 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
870 }
871 }
872}
873
874//--------------------------------------------------------------------------------
875//
876// Retrieves from MGeomCam:
877// - number of pixels
878// - number of pixel areas
879// - number of sectors
880//
881// Return kFALSE, if sizes of the TOrdCollections do not match npixels, nareas or nsectors
882//
883// Call FillHists()
884//
885Bool_t MHCalibrationCam::Fill(const MParContainer *par, const Stat_t w)
886{
887 if (fCurrentNumEvts >= fMaxNumEvts)
888 return kTRUE;
889
890 fCurrentNumEvts++;
891
892 SetIsReset(kFALSE);
893
894 if (!IsSizeCheck())
895 return FillHists(par,w);
896
897 const Int_t npixels = fGeom->GetNumPixels();
898 const Int_t nareas = fGeom->GetNumAreas();
899 const Int_t nsectors = fGeom->GetNumSectors();
900
901 //
902 // Hi-Gain OrdCollections
903 //
904 if (fHiGainArray->GetSize() != npixels)
905 {
906 *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
907 return kFALSE;
908 }
909
910 if (IsLoGain())
911 {
912 if (fLoGainArray->GetSize() != npixels)
913 {
914 *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
915 return kFALSE;
916 }
917 }
918
919 if (!IsAverageing())
920 return FillHists(par,w);
921
922 if (fAverageHiGainAreas->GetSize() != nareas)
923 {
924 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
925 return kFALSE;
926 }
927
928 if (fAverageHiGainSectors->GetSize() != nsectors)
929 {
930 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
931 return kFALSE;
932 }
933
934 if (IsLoGain())
935 {
936
937 if (fAverageLoGainAreas->GetSize() != nareas)
938 {
939 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
940 return kFALSE;
941 }
942
943 if (fAverageLoGainSectors->GetSize() != nsectors)
944 {
945 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
946 return kFALSE;
947 }
948 }
949
950 return FillHists(par, w);
951}
952
953// --------------------------------------------------------------------------
954//
955// 0) Ask if fHiGainArray and fLoGainArray have been initialized,
956// otherwise return kFALSE.
957// 1) FinalizeHists()
958// 2) FinalizeBadPixels()
959// 3) CalcAverageSigma()
960//
961Bool_t MHCalibrationCam::Finalize()
962{
963 if (IsReset())
964 return kTRUE;
965
966 if (GetNumExecutions() < 2)
967 return kTRUE;
968
969 *fLog << inf << GetDescriptor() << ": Number of event used to fill histograms == " << fCurrentNumEvts << endl;
970
971 if (fHiGainArray->GetSize() == 0 && fLoGainArray->GetSize() == 0)
972 {
973 *fLog << err << GetDescriptor()
974 << ": ERROR - Both (HiGain and LoGain) histogram arrays have not been initialized... abort." << endl;
975 return kFALSE;
976 }
977
978 for (Int_t i=0; i<fAverageHiGainAreas->GetSize(); i++)
979 {
980 TH1F *h = GetAverageHiGainArea(i).GetHGausHist();
981 switch (fColor)
982 {
983 case MCalibrationCam::kNONE:
984 break;
985 case MCalibrationCam::kBLUE:
986 h->SetTitle( Form("%s%s", h->GetTitle(),"BLUE "));
987 break;
988 case MCalibrationCam::kGREEN:
989 h->SetTitle( Form("%s%s", h->GetTitle(),"GREEN "));
990 break;
991 case MCalibrationCam::kUV:
992 h->SetTitle( Form("%s%s", h->GetTitle(),"UV "));
993 break;
994 case MCalibrationCam::kCT1:
995 h->SetTitle( Form("%s%s", h->GetTitle(),"CT1-Pulser "));
996 break;
997 }
998 }
999
1000 for (Int_t i=0; i<fAverageLoGainAreas->GetSize(); i++)
1001 {
1002 TH1F *h = GetAverageLoGainArea(i).GetHGausHist();
1003 switch (fColor)
1004 {
1005 case MCalibrationCam::kNONE:
1006 break;
1007 case MCalibrationCam::kBLUE:
1008 h->SetTitle( Form("%s%s", h->GetTitle(),"BLUE "));
1009 break;
1010 case MCalibrationCam::kGREEN:
1011 h->SetTitle( Form("%s%s", h->GetTitle(),"GREEN "));
1012 break;
1013 case MCalibrationCam::kUV:
1014 h->SetTitle( Form("%s%s", h->GetTitle(),"UV "));
1015 break;
1016 case MCalibrationCam::kCT1:
1017 h->SetTitle( Form("%s%s", h->GetTitle(),"CT1-Pulser "));
1018 break;
1019 }
1020 }
1021
1022 if (!FinalizeHists())
1023 return kFALSE;
1024
1025
1026 FinalizeBadPixels();
1027 CalcAverageSigma();
1028
1029 return kTRUE;
1030}
1031
1032// -------------------------------------------------------------
1033//
1034// If MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnsuitableRun):
1035// - calls MHCalibrationPix::SetExcluded()
1036//
1037// Calls:
1038// - MHGausEvents::InitBins()
1039// - MHCalibrationPix::ChangeHistId(i)
1040// - MHCalibrationPix::SetEventFrequency(fPulserFrequency)
1041//
1042void MHCalibrationCam::InitHists(MHCalibrationPix &hist, MBadPixelsPix &bad, const Int_t i)
1043{
1044
1045 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1046 hist.SetExcluded();
1047
1048 hist.InitBins();
1049 hist.SetEventFrequency(fPulserFrequency);
1050}
1051
1052// -------------------------------------------------------------
1053//
1054// - Searches for the CalibrationIntensity*Cam corresponding to 'name'.
1055// - In case, it does not exist in the parameter list, it searches
1056// for the corresponding MCalibration*Cam.
1057// - Initializes the MCalibration*Cam, if not yet done.
1058//
1059Bool_t MHCalibrationCam::InitCams( MParList *plist, const TString name )
1060{
1061
1062 TString intensname = "MCalibrationIntensity";
1063 intensname += name;
1064 intensname += "Cam";
1065
1066 TString ordname = "MCalibration";
1067 ordname += name;
1068 ordname += "Cam";
1069
1070 fIntensCam = (MCalibrationIntensityCam*)plist->FindObject(AddSerialNumber(intensname));
1071 if (fIntensCam)
1072 *fLog << inf << "Found " << intensname << "... " << endl;
1073 else
1074 {
1075 fCam = (MCalibrationCam*)plist->FindObject(AddSerialNumber(ordname));
1076 if (!fCam)
1077 {
1078 fCam = (MCalibrationCam*)plist->FindCreateObj(AddSerialNumber(ordname));
1079 if (!fCam)
1080 return kFALSE;
1081
1082 fCam->Init(*fGeom);
1083 }
1084 }
1085 return kTRUE;
1086}
1087
1088// --------------------------------------------------------------------------
1089//
1090// Calls FitHiGainHists for every entry in:
1091// - fHiGainArray
1092// - fAverageHiGainAreas
1093// - fAverageHiGainSectors
1094//
1095void MHCalibrationCam::FitHiGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
1096 MBadPixelsPix::UncalibratedType_t fittyp,
1097 MBadPixelsPix::UncalibratedType_t osctyp)
1098{
1099 fIsHiGainFitRanges = TMath::Abs(fUpperFitLimitHiGain - fLowerFitLimitHiGain) > 1E-5;
1100
1101 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
1102 {
1103
1104 MHCalibrationPix &hist = (*this)[i];
1105
1106 if (hist.IsExcluded())
1107 continue;
1108
1109 MCalibrationPix &pix = calcam[i];
1110 MBadPixelsPix &bad = badcam[i];
1111
1112 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
1113 }
1114
1115 if (!IsAverageing())
1116 return;
1117
1118 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
1119 {
1120
1121 MHCalibrationPix &hist = GetAverageHiGainArea(j);
1122 MCalibrationPix &pix = calcam.GetAverageArea(j);
1123 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
1124
1125 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
1126 }
1127
1128 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
1129 {
1130 MHCalibrationPix &hist = GetAverageHiGainSector(j);
1131 MCalibrationPix &pix = calcam.GetAverageSector(j);
1132 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
1133
1134 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
1135 }
1136}
1137
1138// --------------------------------------------------------------------------
1139//
1140// Calls FitLoGainHists for every entry in:
1141// - fLoGainArray
1142// - fAverageLoGainAreas
1143// - fAverageLoGainSectors
1144//
1145void MHCalibrationCam::FitLoGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
1146 MBadPixelsPix::UncalibratedType_t fittyp,
1147 MBadPixelsPix::UncalibratedType_t osctyp)
1148{
1149 fIsLoGainFitRanges = TMath::Abs(fUpperFitLimitLoGain - fLowerFitLimitLoGain) > 1E-5;
1150
1151 if (!IsLoGain())
1152 return;
1153
1154 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
1155 {
1156
1157 MHCalibrationPix &hist = (*this)(i);
1158
1159 if (hist.IsExcluded())
1160 continue;
1161
1162 MCalibrationPix &pix = calcam[i];
1163 MBadPixelsPix &bad = badcam[i];
1164
1165 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1166
1167 }
1168
1169 if (!IsAverageing())
1170 return;
1171
1172 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
1173 {
1174
1175 MHCalibrationPix &hist = GetAverageLoGainArea(j);
1176 MCalibrationPix &pix = calcam.GetAverageArea(j);
1177 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
1178
1179 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1180 }
1181
1182 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
1183 {
1184
1185 MHCalibrationPix &hist = GetAverageLoGainSector(j);
1186 MCalibrationPix &pix = calcam.GetAverageSector(j);
1187 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
1188
1189 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
1190 }
1191}
1192
1193//------------------------------------------------------------
1194//
1195// For all averaged areas, the fitted sigma is multiplied with the square root of
1196// the number involved pixels
1197//
1198void MHCalibrationCam::CalcAverageSigma()
1199{
1200 if (!IsAverageing())
1201 return;
1202
1203 MCalibrationCam *cam = fIntensCam ? fIntensCam->GetCam() : fCam;
1204 if (!cam)
1205 return;
1206
1207 for (UInt_t j=0; j<fGeom->GetNumAreas(); j++)
1208 {
1209
1210 MCalibrationPix &pix = cam->GetAverageArea(j);
1211
1212 const Float_t numsqr = TMath::Sqrt((Float_t)fAverageAreaNum[j]);
1213 fAverageAreaSigma[j] = pix.GetSigma () * numsqr;
1214 fAverageAreaSigmaVar[j] = pix.GetSigmaErr () * pix.GetSigmaErr() * numsqr;
1215
1216 pix.SetSigma (fAverageAreaSigma[j]);
1217 pix.SetSigmaVar(fAverageAreaSigmaVar[j]);
1218
1219 fAverageAreaRelSigma [j] = fAverageAreaSigma[j] / pix.GetMean();
1220 fAverageAreaRelSigmaVar[j] = fAverageAreaSigmaVar[j] / (fAverageAreaSigma[j]*fAverageAreaSigma[j]);
1221 fAverageAreaRelSigmaVar[j] += pix.GetMeanRelVar();
1222 fAverageAreaRelSigmaVar[j] *= fAverageAreaRelSigma[j];
1223 }
1224
1225 for (UInt_t j=0; j<fGeom->GetNumSectors(); j++)
1226 {
1227 MCalibrationPix &pix = cam->GetAverageSector(j);
1228
1229 const Float_t numsqr = TMath::Sqrt((Float_t)fAverageSectorNum[j]);
1230 pix.SetSigma (pix.GetSigma() * numsqr);
1231 pix.SetSigmaVar(pix.GetSigmaErr() * pix.GetSigmaErr() * numsqr);
1232 }
1233}
1234
1235// ---------------------------------------------------------------------------
1236//
1237// Returns if the histogram is empty and sets the following flag:
1238// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
1239//
1240// Fits the histograms with a Gaussian, in case of failure
1241// calls MHCalibrationPix::RepeatFit(), in case of repeated failure
1242// calls MHCalibrationPix::BypassFit() and sets the following flags:
1243// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
1244// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1245//
1246// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
1247// In case no, sets the following flags:
1248// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
1249// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1250//
1251// Retrieves the results and store them in MCalibrationPix
1252//
1253void MHCalibrationCam::FitHiGainHists(MHCalibrationPix &hist,
1254 MCalibrationPix &pix,
1255 MBadPixelsPix &bad,
1256 MBadPixelsPix::UncalibratedType_t fittyp,
1257 MBadPixelsPix::UncalibratedType_t osctyp)
1258{
1259
1260 if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
1261 {
1262 *fLog << warn << GetDescriptor()
1263 << ": Only over- or underflow in hi-gain pixel: " << pix.GetPixId() << endl;
1264 return;
1265 }
1266 //
1267 // 2) Fit the Hi Gain histograms with a Gaussian
1268 //
1269 if (fIsHiGainFitRanges)
1270 {
1271 if (!hist.FitGaus("R",fLowerFitLimitHiGain,fUpperFitLimitHiGain))
1272 bad.SetUncalibrated( fittyp );
1273 }
1274 else
1275 if (!hist.FitGaus())
1276 //
1277 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
1278 //
1279 if (!hist.RepeatFit())
1280 {
1281 hist.BypassFit();
1282 bad.SetUncalibrated( fittyp );
1283 }
1284
1285
1286 //
1287 // 4) Check for oscillations
1288 //
1289 if (IsOscillations())
1290 {
1291 hist.CreateFourierSpectrum();
1292
1293 if (!hist.IsFourierSpectrumOK())
1294 bad.SetUncalibrated( osctyp );
1295 }
1296
1297 //
1298 // 5) Retrieve the results and store them in this class
1299 //
1300 pix.SetHiGainMean ( hist.GetMean() );
1301 pix.SetHiGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
1302 pix.SetHiGainRms ( hist.GetHistRms() );
1303 pix.SetHiGainSigma ( hist.GetSigma() );
1304 pix.SetHiGainSigmaVar ( hist.GetSigmaErr()* hist.GetSigmaErr() );
1305 pix.SetHiGainProb ( hist.GetProb() );
1306 pix.SetHiGainNumBlackout( hist.GetBlackout() );
1307 pix.SetHiGainNumPickup ( hist.GetPickup() );
1308
1309 if (IsDebug())
1310 {
1311 *fLog << dbginf << GetDescriptor() << ": ID " << GetName()
1312 << " "<<pix.GetPixId()
1313 << " HiGainSaturation: " << pix.IsHiGainSaturation()
1314 << " HiGainMean: " << hist.GetMean ()
1315 << " HiGainMeanErr: " << hist.GetMeanErr ()
1316 << " HiGainMeanSigma: " << hist.GetSigma ()
1317 << " HiGainMeanSigmaErr: " << hist.GetSigmaErr()
1318 << " HiGainMeanProb: " << hist.GetProb ()
1319 << " HiGainNumBlackout: " << hist.GetBlackout()
1320 << " HiGainNumPickup : " << hist.GetPickup ()
1321 << endl;
1322 }
1323
1324}
1325
1326
1327// ---------------------------------------------------------------------------
1328//
1329// Returns if the histogram is empty and sets the following flag:
1330// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
1331//
1332// Fits the histograms with a Gaussian, in case of failure
1333// calls MHCalibrationPix::RepeatFit(), in case of repeated failure
1334// calls MHCalibrationPix::BypassFit() and sets the following flags:
1335// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
1336// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1337//
1338// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
1339// In case no, sets the following flags:
1340// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
1341// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
1342//
1343// Retrieves the results and store them in MCalibrationPix
1344//
1345void MHCalibrationCam::FitLoGainHists(MHCalibrationPix &hist,
1346 MCalibrationPix &pix,
1347 MBadPixelsPix &bad,
1348 MBadPixelsPix::UncalibratedType_t fittyp,
1349 MBadPixelsPix::UncalibratedType_t osctyp)
1350{
1351
1352 if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
1353 {
1354 *fLog << warn << GetDescriptor()
1355 << ": Only over- or underflow in lo-gain pixel: " << pix.GetPixId() << endl;
1356 return;
1357 }
1358 //
1359 // 2) Fit the Hi Gain histograms with a Gaussian
1360 //
1361 if (fIsLoGainFitRanges)
1362 {
1363 if (!hist.FitGaus("R",fLowerFitLimitLoGain,fUpperFitLimitLoGain))
1364 bad.SetUncalibrated( fittyp );
1365 }
1366 else
1367 if (!hist.FitGaus())
1368 //
1369 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
1370 //
1371 if (!hist.RepeatFit())
1372 {
1373 hist.BypassFit();
1374 if (pix.IsHiGainSaturation())
1375 bad.SetUncalibrated( fittyp );
1376 }
1377
1378 //
1379 // 4) Check for oscillations
1380 //
1381 if (IsOscillations())
1382 {
1383 hist.CreateFourierSpectrum();
1384
1385 if (!hist.IsFourierSpectrumOK())
1386 bad.SetUncalibrated( osctyp );
1387 }
1388
1389 //
1390 // 5) Retrieve the results and store them in this class
1391 //
1392 pix.SetLoGainMean ( hist.GetMean() );
1393 pix.SetLoGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
1394 pix.SetLoGainRms ( hist.GetHistRms() );
1395 pix.SetLoGainSigma ( hist.GetSigma() );
1396 pix.SetLoGainSigmaVar ( hist.GetSigmaErr() * hist.GetSigmaErr() );
1397 pix.SetLoGainProb ( hist.GetProb() );
1398 pix.SetLoGainNumBlackout( hist.GetBlackout() );
1399 pix.SetLoGainNumPickup ( hist.GetPickup() );
1400
1401 if (IsDebug())
1402 {
1403 *fLog << dbginf << GetDescriptor() << "ID: " << hist.GetName()
1404 << " "<<pix.GetPixId()
1405 << " HiGainSaturation: " << pix.IsHiGainSaturation()
1406 << " LoGainMean: " << hist.GetMean ()
1407 << " LoGainMeanErr: " << hist.GetMeanErr ()
1408 << " LoGainMeanSigma: " << hist.GetSigma ()
1409 << " LoGainMeanSigmaErr: " << hist.GetSigmaErr()
1410 << " LoGainMeanProb: " << hist.GetProb ()
1411 << " LoGainNumBlackout: " << hist.GetBlackout()
1412 << " LoGainNumPickup : " << hist.GetPickup ()
1413 << endl;
1414 }
1415
1416}
1417
1418
1419
1420// -----------------------------------------------------------------------------
1421//
1422// Default draw:
1423//
1424// Displays the averaged areas, both High Gain and Low Gain
1425//
1426// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1427//
1428void MHCalibrationCam::Draw(const Option_t *opt)
1429{
1430
1431 if (!IsAverageing())
1432 return;
1433
1434 const Int_t nareas = fAverageHiGainAreas->GetSize();
1435 if (nareas == 0)
1436 return;
1437
1438 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
1439 pad->SetBorderMode(0);
1440
1441 pad->Divide(IsLoGain() ? 2 : 1,nareas);
1442
1443 for (Int_t i=0; i<nareas;i++)
1444 {
1445
1446 pad->cd(IsLoGain() ? 2*(i+1)-1 : i+1);
1447 GetAverageHiGainArea(i).Draw(opt);
1448
1449 if (!fAverageAreaSat[i])
1450 DrawAverageSigma(fAverageAreaSat[i], i,
1451 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
1452 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
1453
1454 if (IsLoGain())
1455 {
1456 pad->cd(2*(i+1));
1457 GetAverageLoGainArea(i).Draw(opt);
1458 }
1459
1460 if (fAverageAreaSat[i])
1461 DrawAverageSigma(fAverageAreaSat[i], i,
1462 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
1463 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
1464 }
1465
1466}
1467
1468// -----------------------------------------------------------------------------
1469//
1470// Default draw:
1471//
1472// Displays a TPaveText with the re-normalized sigmas of the average area
1473//
1474void MHCalibrationCam::DrawAverageSigma(Bool_t sat, Bool_t inner,
1475 Float_t sigma, Float_t sigmavar,
1476 Float_t relsigma, Float_t relsigmavar) const
1477{
1478
1479 if (sigma != 0 && sigmavar >= 0 && relsigmavar >= 0.)
1480 {
1481
1482 TPad *newpad = new TPad("newpad","transparent",0,0,1,1);
1483 newpad->SetFillStyle(4000);
1484 newpad->Draw();
1485 newpad->cd();
1486
1487 TPaveText *text = new TPaveText(sat? 0.1 : 0.35,0.7,sat ? 0.4 : 0.7,1.0);
1488 text->SetTextSize(0.05);
1489 const TString line1 = Form("%s%s%s",inner ? "Outer" : "Inner",
1490 " Pixels ", sat ? "Low Gain" : "High Gain");
1491 TText *txt1 = text->AddText(line1.Data());
1492 const TString line2 = Form("#sigma per pix: %2.2f #pm %2.2f",sigma,TMath::Sqrt(sigmavar));
1493 TText *txt2 = text->AddText(line2.Data());
1494 const TString line3 = Form("Rel. #sigma per pix: %2.2f #pm %2.2f",relsigma,TMath::Sqrt(relsigmavar));
1495 TText *txt3 = text->AddText(line3.Data());
1496 text->Draw("");
1497
1498 text->SetBit(kCanDelete);
1499 txt1->SetBit(kCanDelete);
1500 txt2->SetBit(kCanDelete);
1501 txt3->SetBit(kCanDelete);
1502 newpad->SetBit(kCanDelete);
1503 }
1504}
1505
1506// -----------------------------------------------------------------------------
1507//
1508// Available options
1509// Debug
1510// LoGain
1511// Oscillations
1512// SizeCheck
1513// Averageing
1514// Nbins
1515// First
1516// Last
1517// ProbLimit
1518// OverflowLimit
1519// PulserFrequency
1520// LowerFitLimitHiGain
1521// UpperFitLimitHiGain
1522// LowerFitLimitLoGain
1523// UpperFitLimitLoGain
1524//
1525Int_t MHCalibrationCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1526{
1527
1528 Bool_t rc = kFALSE;
1529 if (IsEnvDefined(env, prefix, "Debug", print))
1530 {
1531 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
1532 rc = kTRUE;
1533 }
1534 if (IsEnvDefined(env, prefix, "LoGain", print))
1535 {
1536 SetDebug(GetEnvValue(env, prefix, "LoGain", IsLoGain()));
1537 rc = kTRUE;
1538 }
1539 if (IsEnvDefined(env, prefix, "Oscillations", print))
1540 {
1541 SetOscillations(GetEnvValue(env, prefix, "Oscillations", IsOscillations()));
1542 rc = kTRUE;
1543 }
1544 if (IsEnvDefined(env, prefix, "SizeCheck", print))
1545 {
1546 SetSizeCheck(GetEnvValue(env, prefix, "SizeCheck", IsSizeCheck()));
1547 rc = kTRUE;
1548 }
1549 if (IsEnvDefined(env, prefix, "Averageing", print))
1550 {
1551 SetAverageing(GetEnvValue(env, prefix, "Averageing", IsAverageing()));
1552 rc = kTRUE;
1553 }
1554
1555 if (IsEnvDefined(env, prefix, "Nbins", print))
1556 {
1557 SetNbins(GetEnvValue(env, prefix, "Nbins", fNbins));
1558 rc = kTRUE;
1559 }
1560 if (IsEnvDefined(env, prefix, "First", print))
1561 {
1562 SetFirst(GetEnvValue(env, prefix, "First", fFirst));
1563 rc = kTRUE;
1564 }
1565
1566 if (IsEnvDefined(env, prefix, "Last", print))
1567 {
1568 SetLast(GetEnvValue(env, prefix, "Last", fLast));
1569 rc = kTRUE;
1570 }
1571
1572 if (IsEnvDefined(env, prefix, "ProbLimit", print))
1573 {
1574 SetProbLimit(GetEnvValue(env, prefix, "ProbLimit", fProbLimit));
1575 rc = kTRUE;
1576 }
1577
1578 if (IsEnvDefined(env, prefix, "MaxNumEvts", print))
1579 {
1580 SetMaxNumEvts(GetEnvValue(env, prefix, "MaxNumEvts", fMaxNumEvts));
1581 rc = kTRUE;
1582 }
1583
1584 if (IsEnvDefined(env, prefix, "OverflowLimit", print))
1585 {
1586 SetOverflowLimit(GetEnvValue(env, prefix, "OverflowLimit", fOverflowLimit));
1587 rc = kTRUE;
1588 }
1589
1590 if (IsEnvDefined(env, prefix, "PulserFrequency", print))
1591 {
1592 SetPulserFrequency(GetEnvValue(env, prefix, "PulserFrequency", fPulserFrequency));
1593 rc = kTRUE;
1594 }
1595
1596 if (IsEnvDefined(env, prefix, "LowerFitLimitHiGain", print))
1597 {
1598 SetLowerFitLimitHiGain(GetEnvValue(env, prefix, "LowerFitLimitHiGain", fLowerFitLimitHiGain));
1599 rc = kTRUE;
1600 }
1601
1602 if (IsEnvDefined(env, prefix, "UpperFitLimitHiGain", print))
1603 {
1604 SetUpperFitLimitHiGain(GetEnvValue(env, prefix, "UpperFitLimitHiGain", fUpperFitLimitHiGain));
1605 rc = kTRUE;
1606 }
1607
1608 if (IsEnvDefined(env, prefix, "LowerFitLimitLoGain", print))
1609 {
1610 SetLowerFitLimitLoGain(GetEnvValue(env, prefix, "LowerFitLimitLoGain", fLowerFitLimitLoGain));
1611 rc = kTRUE;
1612 }
1613
1614 if (IsEnvDefined(env, prefix, "UpperFitLimitLoGain", print))
1615 {
1616 SetUpperFitLimitLoGain(GetEnvValue(env, prefix, "UpperFitLimitLoGain", fUpperFitLimitLoGain));
1617 rc = kTRUE;
1618 }
1619
1620
1621 return rc;
1622}
1623
Note: See TracBrowser for help on using the repository browser.