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

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