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

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