source: trunk/MagicSoft/Mars/mcalib/MHCalibrationCam.cc@ 4693

Last change on this file since 4693 was 4630, checked in by gaug, 20 years ago
*** empty log message ***
File size: 31.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24/////////////////////////////////////////////////////////////////////////////
25//
26// MHCalibrationCam
27//
28// Base class for camera calibration classes. Incorporates the TObjArray's:
29// - fHiGainArray (for calibrated High Gains per pixel)
30// - fLoGainArray (for calibrated Low Gains per pixel)
31// - fAverageHiGainAreas (for averaged High Gains events per camera area index)
32// - fAverageLoGainAreas (for averaged High Gains events per camera area index)
33// - fAverageHiGainSectors (for averaged High Gains events per camera sector )
34// - fAverageLoGainSectors (for averaged High Gains events per camera sector )
35// These TObjArray's are called by their default constructors, thus no objects
36// are created, until the derived class does so.
37//
38// The corresponding operators: [],() and the operators GetAverageHiGainArea(),
39// GetAverageLoGainArea(), GetAverageHiGainSector() and GetAverageLoGainSector()
40// have to be cast to the corresponding class. It is assumed that all classes
41// dealing with calibration pixels derive from MHGausEvents.
42//
43/////////////////////////////////////////////////////////////////////////////
44#include "MHCalibrationCam.h"
45
46#include <TVirtualPad.h>
47#include <TCanvas.h>
48#include <TPad.h>
49#include <TText.h>
50#include <TPaveText.h>
51
52#include "MLog.h"
53#include "MLogManip.h"
54
55#include "MCalibrationPix.h"
56#include "MCalibrationCam.h"
57
58#include "MHGausEvents.h"
59
60#include "MBadPixelsPix.h"
61#include "MBadPixelsCam.h"
62
63#include "MGeomCam.h"
64#include "MGeomPix.h"
65
66#include "MParList.h"
67
68#include "MRawRunHeader.h"
69
70ClassImp(MHCalibrationCam);
71
72using namespace std;
73
74const Int_t MHCalibrationCam::fgAverageNbins = 2000;
75const Int_t MHCalibrationCam::fgPulserFrequency = 500;
76// --------------------------------------------------------------------------
77//
78// Default Constructor.
79//
80// Sets:
81// - all pointers to NULL
82//
83// Initializes and sets owner of:
84// - fHiGainArray, fLoGainArray
85// - fAverageHiGainAreas, fAverageLoGainAreas
86// - fAverageHiGainSectors, fAverageLoGainSectors
87//
88// Initializes:
89// - fPulserFrequency to fgPulserFrequency
90//
91MHCalibrationCam::MHCalibrationCam(const char *name, const char *title)
92 : fColor(MCalibrationCam::kNONE),
93 fBadPixels(NULL), fCam(NULL), fGeom(NULL), fRunHeader(NULL), fDebug(kFALSE)
94{
95
96 fHiGainArray = new TObjArray;
97 fHiGainArray->SetOwner();
98
99 fLoGainArray = new TObjArray;
100 fLoGainArray->SetOwner();
101
102 fAverageHiGainAreas = new TObjArray;
103 fAverageHiGainAreas->SetOwner();
104
105 fAverageLoGainAreas = new TObjArray;
106 fAverageLoGainAreas->SetOwner();
107
108 fAverageHiGainSectors = new TObjArray;
109 fAverageHiGainSectors->SetOwner();
110
111 fAverageLoGainSectors = new TObjArray;
112 fAverageLoGainSectors->SetOwner();
113
114 SetAverageNbins();
115 SetPulserFrequency();
116
117 fHiGainOverFlow = 0;
118 fLoGainOverFlow = 0;
119}
120
121// --------------------------------------------------------------------------
122//
123// Deletes the TClonesArray of:
124// - fHiGainArray, fLoGainArray
125// - fAverageHiGainAreas, fAverageLoGainAreas
126// - fAverageHiGainSectors, fAverageLoGainSectors
127//
128MHCalibrationCam::~MHCalibrationCam()
129{
130
131 delete fHiGainArray;
132 delete fLoGainArray;
133
134 delete fAverageHiGainAreas;
135 delete fAverageLoGainAreas;
136
137 delete fAverageHiGainSectors;
138 delete fAverageLoGainSectors;
139}
140
141// --------------------------------------------------------------------------
142//
143// Get i-th High Gain pixel (pixel number)
144//
145MHGausEvents &MHCalibrationCam::operator[](UInt_t i)
146{
147 return *static_cast<MHGausEvents*>(fHiGainArray->UncheckedAt(i));
148}
149
150// --------------------------------------------------------------------------
151//
152// Get i-th High Gain pixel (pixel number)
153//
154const MHGausEvents &MHCalibrationCam::operator[](UInt_t i) const
155{
156 return *static_cast<MHGausEvents*>(fHiGainArray->UncheckedAt(i));
157}
158
159// --------------------------------------------------------------------------
160//
161// Get i-th Low Gain pixel (pixel number)
162//
163MHGausEvents &MHCalibrationCam::operator()(UInt_t i)
164{
165 return *static_cast<MHGausEvents*>(fLoGainArray->UncheckedAt(i));
166}
167
168// --------------------------------------------------------------------------
169//
170// Get i-th Low Gain pixel (pixel number)
171//
172const MHGausEvents &MHCalibrationCam::operator()(UInt_t i) const
173{
174 return *static_cast<MHGausEvents*>(fLoGainArray->UncheckedAt(i));
175}
176
177// --------------------------------------------------------------------------
178//
179// Returns the current size of the TObjArray fAverageHiGainAreas
180// independently if the MHGausEvents is filled with values or not.
181//
182const Int_t MHCalibrationCam::GetAverageAreas() const
183{
184 return fAverageHiGainAreas->GetEntries();
185}
186
187// --------------------------------------------------------------------------
188//
189// Get i-th High Gain pixel Area (area number)
190//
191MHGausEvents &MHCalibrationCam::GetAverageHiGainArea(UInt_t i)
192{
193 return *static_cast<MHGausEvents*>(fAverageHiGainAreas->UncheckedAt(i));
194}
195
196// --------------------------------------------------------------------------
197//
198// Get i-th High Gain pixel Area (area number)
199//
200const MHGausEvents &MHCalibrationCam::GetAverageHiGainArea(UInt_t i) const
201{
202 return *static_cast<MHGausEvents *>(fAverageHiGainAreas->UncheckedAt(i));
203}
204
205// --------------------------------------------------------------------------
206//
207// Get i-th Low Gain pixel Area (area number)
208//
209MHGausEvents &MHCalibrationCam::GetAverageLoGainArea(UInt_t i)
210{
211 return *static_cast<MHGausEvents*>(fAverageLoGainAreas->UncheckedAt(i));
212}
213
214// --------------------------------------------------------------------------
215//
216// Get i-th Low Gain pixel Area (area number)
217//
218const MHGausEvents &MHCalibrationCam::GetAverageLoGainArea(UInt_t i) const
219{
220 return *static_cast<MHGausEvents*>(fAverageLoGainAreas->UncheckedAt(i));
221}
222
223// --------------------------------------------------------------------------
224//
225// Returns the current size of the TObjArray fAverageHiGainSectors
226// independently if the MHGausEvents is filled with values or not.
227//
228const Int_t MHCalibrationCam::GetAverageSectors() const
229{
230 return fAverageHiGainSectors->GetEntries();
231}
232
233// --------------------------------------------------------------------------
234//
235// Get i-th High Gain Sector (sector number)
236//
237MHGausEvents &MHCalibrationCam::GetAverageHiGainSector(UInt_t i)
238{
239 return *static_cast<MHGausEvents*>(fAverageHiGainSectors->UncheckedAt(i));
240}
241
242// --------------------------------------------------------------------------
243//
244// Get i-th High Gain Sector (sector number)
245//
246const MHGausEvents &MHCalibrationCam::GetAverageHiGainSector(UInt_t i) const
247{
248 return *static_cast<MHGausEvents*>(fAverageHiGainSectors->UncheckedAt(i));
249}
250
251// --------------------------------------------------------------------------
252//
253// Get i-th Low Gain Sector (sector number)
254//
255MHGausEvents &MHCalibrationCam::GetAverageLoGainSector(UInt_t i)
256{
257 return *static_cast<MHGausEvents*>(fAverageLoGainSectors->UncheckedAt(i));
258}
259
260// --------------------------------------------------------------------------
261//
262// Get i-th Low Gain Sector (sector number)
263//
264const MHGausEvents &MHCalibrationCam::GetAverageLoGainSector(UInt_t i) const
265{
266 return *static_cast<MHGausEvents*>(fAverageLoGainSectors->UncheckedAt(i));
267}
268
269
270const TArrayI &MHCalibrationCam::GetRunNumbers() const
271{
272 return fRunNumbers;
273}
274
275// --------------------------------------------------------------------------
276//
277// Our own clone function is necessary since root 3.01/06 or Mars 0.4
278// I don't know the reason.
279//
280// Creates new MHCalibrationCam
281// Deletes the TObjArray's and Clones them individually
282// Copies the TArray's
283// Copies the fPulserFrequency
284//
285TObject *MHCalibrationCam::Clone(const char *) const
286{
287
288 // const Int_t nhi = fHiGainArray->GetEntries();
289 // const Int_t nlo = fLoGainArray->GetEntries();
290 const Int_t navhi = fAverageHiGainAreas->GetEntries();
291 const Int_t navlo = fAverageLoGainAreas->GetEntries();
292 const Int_t nsehi = fAverageHiGainSectors->GetEntries();
293 const Int_t nselo = fAverageLoGainSectors->GetEntries();
294
295 //
296 // FIXME, this might be done faster and more elegant, by direct copy.
297 //
298 MHCalibrationCam *cam = new MHCalibrationCam();
299
300 // cam->fHiGainArray->Expand(nhi);
301 // cam->fLoGainArray->Expand(nlo);
302 cam->fAverageHiGainAreas->Expand(navhi);
303 cam->fAverageLoGainAreas->Expand(navlo);
304 cam->fAverageHiGainSectors->Expand(nsehi);
305 cam->fAverageLoGainSectors->Expand(nselo);
306
307 /*
308 for (int i=0; i<nhi; i++)
309 {
310 delete (*cam->fHiGainArray)[i];
311 (*cam->fHiGainArray)[i] = (*fHiGainArray)[i]->Clone();
312 }
313 for (int i=0; i<nlo; i++)
314 {
315 delete (*cam->fLoGainArray)[i];
316 (*cam->fLoGainArray)[i] = (*fLoGainArray)[i]->Clone();
317 }
318 */
319
320 for (int i=0; i<navhi; i++)
321 {
322 // delete (*cam->fAverageHiGainAreas)[i];
323 (*cam->fAverageHiGainAreas)[i] = (*fAverageHiGainAreas)[i]->Clone();
324 }
325 for (int i=0; i<navlo; i++)
326 {
327 // delete (*cam->fAverageLoGainAreas)[i];
328 (*cam->fAverageLoGainAreas)[i] = (*fAverageLoGainAreas)[i]->Clone();
329 }
330 for (int i=0; i<nsehi; i++)
331 {
332 // delete (*cam->fAverageHiGainSectors)[i];
333 (*cam->fAverageHiGainSectors)[i] = (*fAverageHiGainSectors)[i]->Clone();
334 }
335 for (int i=0; i<nselo; i++)
336 {
337 // delete (*cam->fAverageLoGainSectors)[i];
338 (*cam->fAverageLoGainSectors)[i] = (*fAverageLoGainSectors)[i]->Clone();
339 }
340
341 cam->fAverageAreaNum = fAverageAreaNum;
342 cam->fAverageAreaSat = fAverageAreaSat;
343 cam->fAverageAreaSigma = fAverageAreaSigma;
344 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
345 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
346 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
347 cam->fAverageSectorNum = fAverageSectorNum;
348 cam->fRunNumbers = fRunNumbers;
349
350 cam->fPulserFrequency = fPulserFrequency;
351 cam->fAverageNbins = fAverageNbins;
352
353 return cam;
354}
355
356// --------------------------------------------------------------------------
357//
358// Gets the pointers to:
359// - MGeomCam
360//
361// Calls SetupHists(const MParList *pList)
362//
363// Calls Delete-Function of:
364// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
365// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
366// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
367//
368Bool_t MHCalibrationCam::SetupFill(const MParList *pList)
369{
370
371 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
372 if (!fGeom)
373 {
374 *fLog << err << GetDescriptor()
375 << ": MGeomCam not found... aborting." << endl;
376 return kFALSE;
377 }
378
379 fRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
380 if (!fRunHeader)
381 {
382 *fLog << warn << GetDescriptor()
383 << ": MRawRunHeader not found... will not store run numbers." << endl;
384 }
385
386 fHiGainArray->Delete();
387 fLoGainArray->Delete();
388
389 fAverageHiGainAreas->Delete();
390 fAverageLoGainAreas->Delete();
391
392 fAverageHiGainSectors->Delete();
393 fAverageLoGainSectors->Delete();
394
395 return SetupHists(pList);
396}
397
398
399Bool_t MHCalibrationCam::SetupHists(const MParList *pList)
400{
401 return kTRUE;
402}
403
404// --------------------------------------------------------------------------
405//
406// Gets or creates the pointers to:
407// - MBadPixelsCam
408//
409// Searches pointer to:
410// - MArrivalTimeCam
411//
412// Initializes, if empty to MArrivalTimeCam::GetSize() for:
413// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
414//
415// Initializes, if empty to MGeomCam::GetNumAreas() for:
416// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
417//
418// Initializes, if empty to MGeomCam::GetNumSectors() for:
419// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
420//
421// Initializes TArray's to MGeomCam::GetNumAreas and MGeomCam::GetNumSectors, respectively
422// Fills with number of valid pixels (if !MBadPixelsPix::IsBad()):
423// - MHCalibrationCam::fAverageAreaNum[area index]
424// - MHCalibrationCam::fAverageSectorNum[area index]
425//
426// Calls InitializeHists() for every entry in:
427// - MHCalibrationCam::fHiGainArray
428// - MHCalibrationCam::fAverageHiGainAreas
429// - MHCalibrationCam::fAverageHiGainSectors
430//
431// Sets Titles and Names for the Histograms
432// - MHCalibrationCam::fAverageHiGainAreas
433// - MHCalibrationCam::fAverageHiGainSectors
434//
435// Retrieves the run numbers from MRawRunHeader and stores them in fRunNumbers
436//
437Bool_t MHCalibrationCam::ReInit(MParList *pList)
438{
439
440 const Int_t npixels = fGeom->GetNumPixels();
441 const Int_t nsectors = fGeom->GetNumSectors();
442 const Int_t nareas = fGeom->GetNumAreas();
443
444 fBadPixels = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
445 if (!fBadPixels)
446 {
447
448 fBadPixels = (MBadPixelsCam*)pList->FindCreateObj(AddSerialNumber("MBadPixelsCam"));
449 if (!fBadPixels)
450 {
451 *fLog << err << "Cannot find nor create MBadPixelsCam ... abort." << endl;
452 return kFALSE;
453 }
454 else
455 fBadPixels->InitSize(npixels);
456 }
457
458 //
459 // The function TArrayF::Set() already sets all entries to 0.
460 //
461 fAverageAreaNum. Set(nareas);
462 fAverageAreaSat. Set(nareas);
463 fAverageAreaSigma. Set(nareas);
464 fAverageAreaSigmaVar. Set(nareas);
465 fAverageAreaRelSigma. Set(nareas);
466 fAverageAreaRelSigmaVar.Set(nareas);
467 fAverageSectorNum. Set(nsectors);
468 fRunNumbers. Set(fRunNumbers.GetSize()+1);
469
470 for (Int_t aidx=0; aidx<nareas; aidx++)
471 fAverageAreaNum[aidx] = 0;
472
473 for (Int_t sector=0; sector<nsectors; sector++)
474 fAverageSectorNum[sector] = 0;
475
476 for (Int_t i=0; i<npixels; i++)
477 {
478
479 if ((*fBadPixels)[i].IsBad())
480 continue;
481
482 fAverageAreaNum [(*fGeom)[i].GetAidx() ]++;
483 fAverageSectorNum[(*fGeom)[i].GetSector()]++;
484 }
485
486 if (fRunHeader)
487 fRunNumbers[fRunNumbers.GetSize()-1] = fRunHeader->GetRunNumber();
488
489 if (!ReInitHists(pList))
490 return kFALSE;
491
492 if (!fRunHeader)
493 return kTRUE;
494
495 for (Int_t i=0; i<fHiGainArray->GetEntries(); i++)
496 {
497 TH1F *h = (*this)[i].GetHGausHist();
498 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
499 }
500
501 for (Int_t i=0; i<fLoGainArray->GetEntries(); i++)
502 {
503 TH1F *h = (*this)(i).GetHGausHist();
504 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
505 }
506
507 for (Int_t j=0; j<nareas; j++)
508 {
509 TH1F *h = GetAverageHiGainArea(j).GetHGausHist();
510 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
511 }
512
513 for (Int_t j=0; j<nareas; j++)
514 {
515 TH1F *h = GetAverageLoGainArea(j).GetHGausHist();
516 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
517 }
518
519 for (Int_t j=0; j<nsectors; j++)
520 {
521 TH1F *h = GetAverageHiGainSector(j).GetHGausHist();
522 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
523 }
524
525 for (Int_t j=0; j<nsectors; j++)
526 {
527 TH1F *h = GetAverageLoGainSector(j).GetHGausHist();
528 h->SetTitle( Form("%s%i%s", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]," "));
529 }
530
531 return kTRUE;
532}
533
534
535
536Bool_t MHCalibrationCam::ReInitHists(MParList *pList)
537{
538 return kTRUE;
539}
540
541
542
543//--------------------------------------------------------------------------------
544//
545// Retrieves from MGeomCam:
546// - number of pixels
547// - number of pixel areas
548// - number of sectors
549//
550// For all TObjArray's (including the averaged ones), the following steps are performed:
551//
552// 1) Test size and return kFALSE if not matching
553// 2)
554//
555Bool_t MHCalibrationCam::Fill(const MParContainer *par, const Stat_t w)
556{
557
558 const Int_t npixels = fGeom->GetNumPixels();
559 const Int_t nareas = fGeom->GetNumAreas();
560 const Int_t nsectors = fGeom->GetNumSectors();
561
562 if (fHiGainArray->GetEntries() != npixels)
563 {
564 *fLog << err << "ERROR - Size mismatch... abort." << endl;
565 return kFALSE;
566 }
567
568 if (fLoGainArray->GetEntries() != npixels)
569 {
570 *fLog << err << "ERROR - Size mismatch... abort." << endl;
571 return kFALSE;
572 }
573
574 if (fAverageHiGainAreas->GetEntries() != nareas)
575 {
576 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
577 return kFALSE;
578 }
579
580 if (fAverageLoGainAreas->GetEntries() != nareas)
581 {
582 *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
583 return kFALSE;
584 }
585
586 if (fAverageHiGainSectors->GetEntries() != nsectors)
587 {
588 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
589 return kFALSE;
590 }
591
592 if (fAverageLoGainSectors->GetEntries() != nsectors)
593 {
594 *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
595 return kFALSE;
596 }
597
598 return FillHists(par, w);
599}
600
601Bool_t MHCalibrationCam::FillHists(const MParContainer *par, const Stat_t w)
602{
603 *fLog << warn << GetDescriptor() << "FillHists not overloaded! Can't be used!" << endl;
604 return kFALSE;
605}
606
607// --------------------------------------------------------------------------
608//
609// 0) Ask if fHiGainArray and fLoGainArray have been initialized,
610// otherwise return kFALSE.
611// 1) FinalizeHists()
612// 2) FinalizeBadPixels()
613// 3) CalcAverageSigma()
614//
615Bool_t MHCalibrationCam::Finalize()
616{
617
618 if (fHiGainArray->GetEntries() == 0 && fLoGainArray->GetEntries() == 0)
619 {
620 *fLog << err << GetDescriptor()
621 << ": ERROR: Both (HiGain and LoGain) histogram arrays have not been initialized... abort." << endl;
622 return kFALSE;
623 }
624
625 if (!FinalizeHists())
626 return kFALSE;
627
628 FinalizeBadPixels();
629 CalcAverageSigma();
630
631 return kTRUE;
632}
633
634Bool_t MHCalibrationCam::FinalizeHists()
635{
636 return kTRUE;
637}
638
639void MHCalibrationCam::FinalizeBadPixels()
640{
641}
642
643
644// -------------------------------------------------------------
645//
646// If MBadPixelsPix::IsBad():
647// - calls MHGausEvents::SetExcluded()
648//
649// Calls:
650// - MHGausEvents::InitBins()
651// - MHGausEvents::ChangeHistId(i)
652// - MHGausEvents::SetEventFrequency(fPulserFrequency)
653//
654void MHCalibrationCam::InitHists(MHGausEvents &hist, MBadPixelsPix &bad, const Int_t i)
655{
656
657 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
658 hist.SetExcluded();
659
660 hist.InitBins();
661 hist.ChangeHistId(i);
662 hist.SetEventFrequency(fPulserFrequency);
663
664 TH1F *h = hist.GetHGausHist();
665 h->SetTitle( Form("%s%s", h->GetTitle()," Runs: "));
666}
667
668void MHCalibrationCam::FitHiGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
669 MBadPixelsPix::UncalibratedType_t fittyp,
670 MBadPixelsPix::UncalibratedType_t osctyp)
671{
672
673 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
674 {
675
676 MHGausEvents &hist = (*this)[i];
677
678 if (hist.IsExcluded())
679 continue;
680
681 MCalibrationPix &pix = calcam[i];
682 MBadPixelsPix &bad = badcam[i];
683
684 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
685
686 }
687
688 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
689 {
690
691 MHGausEvents &hist = GetAverageHiGainArea(j);
692 MCalibrationPix &pix = calcam.GetAverageArea(j);
693 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
694
695 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
696 }
697
698
699 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
700 {
701
702 MHGausEvents &hist = GetAverageHiGainSector(j);
703 MCalibrationPix &pix = calcam.GetAverageSector(j);
704 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
705
706 FitHiGainHists(hist,pix,bad,fittyp,osctyp);
707 }
708
709}
710
711void MHCalibrationCam::FitLoGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
712 MBadPixelsPix::UncalibratedType_t fittyp,
713 MBadPixelsPix::UncalibratedType_t osctyp)
714{
715
716 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
717 {
718
719 MHGausEvents &hist = (*this)(i);
720
721 if (hist.IsExcluded())
722 continue;
723
724 MCalibrationPix &pix = calcam[i];
725 MBadPixelsPix &bad = badcam[i];
726
727 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
728
729 }
730
731 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
732 {
733
734 MHGausEvents &hist = GetAverageLoGainArea(j);
735 MCalibrationPix &pix = calcam.GetAverageArea(j);
736 MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
737
738 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
739 }
740
741
742 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
743 {
744
745 MHGausEvents &hist = GetAverageLoGainSector(j);
746 MCalibrationPix &pix = calcam.GetAverageSector(j);
747 MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
748
749 FitLoGainHists(hist,pix,bad,fittyp,osctyp);
750 }
751}
752
753//------------------------------------------------------------
754//
755// For all averaged areas, the fitted sigma is multiplied with the square root of
756// the number involved pixels
757//
758void MHCalibrationCam::CalcAverageSigma()
759{
760
761 for (UInt_t j=0; j<fGeom->GetNumAreas(); j++)
762 {
763
764 MCalibrationPix &pix = fCam->GetAverageArea(j);
765
766 const Float_t numsqr = TMath::Sqrt((Float_t)fAverageAreaNum[j]);
767 fAverageAreaSigma[j] = pix.GetSigma () * numsqr;
768 fAverageAreaSigmaVar[j] = pix.GetSigmaErr () * pix.GetSigmaErr() * numsqr;
769
770 pix.SetSigma (fAverageAreaSigma[j]);
771 pix.SetSigmaVar(fAverageAreaSigmaVar[j]);
772
773 fAverageAreaRelSigma [j] = fAverageAreaSigma[j] / pix.GetMean();
774 fAverageAreaRelSigmaVar[j] = fAverageAreaSigmaVar[j] / (fAverageAreaSigma[j]*fAverageAreaSigma[j]);
775 fAverageAreaRelSigmaVar[j] += pix.GetMeanRelVar();
776 fAverageAreaRelSigmaVar[j] *= fAverageAreaRelSigma[j];
777 }
778}
779
780// ---------------------------------------------------------------------------
781//
782// Returns if the histogram is empty and sets the following flag:
783// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
784//
785// Fits the histograms with a Gaussian, in case of failure
786// calls MHGausEvents::RepeatFit(), in case of repeated failure
787// calls MHGausEvents::BypassFit() and sets the following flags:
788// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
789// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
790//
791// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
792// In case no, sets the following flags:
793// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
794// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
795//
796// Retrieves the results and store them in MCalibrationPix
797//
798void MHCalibrationCam::FitHiGainHists(MHGausEvents &hist,
799 MCalibrationPix &pix,
800 MBadPixelsPix &bad,
801 MBadPixelsPix::UncalibratedType_t fittyp,
802 MBadPixelsPix::UncalibratedType_t osctyp)
803{
804
805
806 if (hist.IsEmpty())
807 return;
808
809 //
810 // 2) Fit the Hi Gain histograms with a Gaussian
811 //
812 if (!hist.FitGaus())
813 //
814 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
815 //
816 if (!hist.RepeatFit())
817 {
818 hist.BypassFit();
819 bad.SetUncalibrated( fittyp );
820 }
821
822 hist.Renorm();
823 //
824 // 4) Check for oscillations
825 //
826 hist.CreateFourierSpectrum();
827
828
829 if (!hist.IsFourierSpectrumOK())
830 bad.SetUncalibrated( osctyp );
831
832 //
833 // 5) Retrieve the results and store them in this class
834 //
835 pix.SetHiGainMean ( hist.GetMean() );
836 pix.SetHiGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
837 pix.SetHiGainSigma ( hist.GetSigma() );
838 pix.SetHiGainSigmaVar ( hist.GetSigmaErr()* hist.GetSigmaErr() );
839 pix.SetHiGainProb ( hist.GetProb() );
840 pix.SetHiGainNumBlackout( hist.GetBlackout() );
841 pix.SetHiGainNumPickup ( hist.GetPickup() );
842
843 if (IsDebug())
844 {
845 *fLog << dbginf << GetDescriptor() << ": ID " << hist.GetPixId()
846 << " HiGainSaturation: " << pix.IsHiGainSaturation()
847 << " HiGainMean: " << hist.GetMean()
848 << " HiGainMeanErr: " << hist.GetMeanErr()
849 << " HiGainMeanSigma: " << hist.GetSigma()
850 << " HiGainMeanSigmaErr: " << hist.GetSigmaErr()
851 << " HiGainMeanProb: " << hist.GetProb()
852 << " HiGainNumBlackout: " << hist.GetBlackout()
853 << " HiGainNumPickup : " << hist.GetPickup ()
854 << endl;
855 }
856
857}
858
859
860// ---------------------------------------------------------------------------
861//
862// Returns if the histogram is empty and sets the following flag:
863// - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
864//
865// Fits the histograms with a Gaussian, in case of failure
866// calls MHGausEvents::RepeatFit(), in case of repeated failure
867// calls MHGausEvents::BypassFit() and sets the following flags:
868// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
869// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
870//
871// Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
872// In case no, sets the following flags:
873// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
874// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
875//
876// Retrieves the results and store them in MCalibrationPix
877//
878void MHCalibrationCam::FitLoGainHists(MHGausEvents &hist,
879 MCalibrationPix &pix,
880 MBadPixelsPix &bad,
881 MBadPixelsPix::UncalibratedType_t fittyp,
882 MBadPixelsPix::UncalibratedType_t osctyp)
883{
884
885 if (hist.IsEmpty())
886 return;
887
888
889 //
890 // 2) Fit the Hi Gain histograms with a Gaussian
891 //
892 if (!hist.FitGaus())
893 //
894 // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
895 //
896 if (!hist.RepeatFit())
897 {
898 hist.BypassFit();
899 bad.SetUncalibrated( fittyp );
900 }
901
902 //
903 // 4) Check for oscillations
904 //
905 hist.CreateFourierSpectrum();
906
907 if (!hist.IsFourierSpectrumOK())
908 bad.SetUncalibrated( osctyp );
909
910 //
911 // 5) Retrieve the results and store them in this class
912 //
913 pix.SetLoGainMean ( hist.GetMean() );
914 pix.SetLoGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
915 pix.SetLoGainSigma ( hist.GetSigma() );
916 pix.SetLoGainSigmaVar ( hist.GetSigmaErr() * hist.GetSigmaErr() );
917 pix.SetLoGainProb ( hist.GetProb() );
918 pix.SetLoGainNumBlackout( hist.GetBlackout() );
919 pix.SetLoGainNumPickup ( hist.GetPickup() );
920
921 if (IsDebug())
922 {
923 *fLog << dbginf << GetDescriptor() << "ID: " << hist.GetPixId()
924 << " HiGainSaturation: " << pix.IsHiGainSaturation()
925 << " LoGainMean: " << hist.GetMean()
926 << " LoGainMeanErr: " << hist.GetMeanErr()
927 << " LoGainMeanSigma: " << hist.GetSigma()
928 << " LoGainMeanSigmaErr: " << hist.GetSigmaErr()
929 << " LoGainMeanProb: " << hist.GetProb()
930 << " LoGainNumBlackout: " << hist.GetBlackout()
931 << " LoGainNumPickup : " << hist.GetPickup ()
932 << endl;
933 }
934
935}
936
937
938
939// --------------------------------------------------------------------------
940//
941// Dummy, needed by MCamEvent
942//
943Bool_t MHCalibrationCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
944{
945 return kTRUE;
946}
947
948// --------------------------------------------------------------------------
949//
950// What MHCamera needs in order to draw an individual pixel in the camera
951//
952void MHCalibrationCam::DrawPixelContent(Int_t idx) const
953{
954}
955
956// -----------------------------------------------------------------------------
957//
958// Default draw:
959//
960// Displays the averaged areas, both High Gain and Low Gain
961//
962// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
963//
964void MHCalibrationCam::Draw(const Option_t *opt)
965{
966
967 const Int_t nareas = fAverageHiGainAreas->GetEntries();
968 if (nareas == 0)
969 return;
970
971 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
972 pad->SetBorderMode(0);
973
974 pad->Divide(2,nareas);
975
976 for (Int_t i=0; i<nareas;i++)
977 {
978
979 pad->cd(2*(i+1)-1);
980 GetAverageHiGainArea(i).Draw(opt);
981
982 if (!fAverageAreaSat[i])
983 DrawAverageSigma(fAverageAreaSat[i], i,
984 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
985 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
986
987 pad->cd(2*(i+1));
988 GetAverageLoGainArea(i).Draw(opt);
989
990 if (fAverageAreaSat[i])
991 DrawAverageSigma(fAverageAreaSat[i], i,
992 fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
993 fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
994 }
995
996}
997
998// -----------------------------------------------------------------------------
999//
1000// Default draw:
1001//
1002// Displays a TPaveText with the re-normalized sigmas of the average area
1003//
1004void MHCalibrationCam::DrawAverageSigma(Bool_t sat, Bool_t inner,
1005 Float_t sigma, Float_t sigmavar,
1006 Float_t relsigma, Float_t relsigmavar) const
1007{
1008
1009 if (sigma != 0 && sigmavar >= 0 && relsigmavar >= 0.)
1010 {
1011
1012 TPad *newpad = new TPad("newpad","transparent",0,0,1,1);
1013 newpad->SetFillStyle(4000);
1014 newpad->Draw();
1015 newpad->cd();
1016
1017 TPaveText *text = new TPaveText(sat? 0.1 : 0.35,0.7,sat ? 0.4 : 0.7,1.0);
1018 text->SetTextSize(0.07);
1019 const TString line1 = Form("%s%s%s",inner ? "Outer" : "Inner",
1020 " Pixels ", sat ? "Low Gain" : "High Gain");
1021 TText *txt1 = text->AddText(line1.Data());
1022 const TString line2 = Form("#sigma per pix: %2.2f #pm %2.2f",sigma,TMath::Sqrt(sigmavar));
1023 TText *txt2 = text->AddText(line2.Data());
1024 const TString line3 = Form("Rel. #sigma per pix: %2.2f #pm %2.2f",relsigma,TMath::Sqrt(relsigmavar));
1025 TText *txt3 = text->AddText(line3.Data());
1026 text->Draw("");
1027
1028 text->SetBit(kCanDelete);
1029 txt1->SetBit(kCanDelete);
1030 txt2->SetBit(kCanDelete);
1031 txt3->SetBit(kCanDelete);
1032 newpad->SetBit(kCanDelete);
1033 }
1034}
1035
1036Int_t MHCalibrationCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1037{
1038 Bool_t rc = kFALSE;
1039 if (IsEnvDefined(env, prefix, "Debug", print))
1040 {
1041 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
1042 rc = kTRUE;
1043 }
1044
1045 return rc;
1046}
Note: See TracBrowser for help on using the repository browser.