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

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