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

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