source: trunk/Mars/mhcalib/MHCalibrationCam.cc@ 18750

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