source: trunk/Mars/mcalib/MCalibrationChargeCam.cc@ 9854

Last change on this file since 9854 was 9374, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 37.3 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MCalibrationChargeCam.cc,v 1.73 2009-03-02 14:32:49 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without express
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Markus Gaug 11/2003 <mailto:markus@ifae.es>
21!
22! Copyright: MAGIC Software Development, 2000-2008
23!
24!
25\* ======================================================================== */
26/////////////////////////////////////////////////////////////////////////////
27//
28// MCalibrationChargeCam
29//
30// Storage container for charge calibration results from the signal distribution
31// fits (see MHCalibrationChargeCam and MHCalibrationChargePix), the calculation
32// of reduced sigmas and number of photo-electrons (this class) and conversion
33// factors sum FADC slices to photo-electrons (see MCalibrationChargeCalc)
34//
35// Individual pixels have to be cast when retrieved e.g.:
36// MCalibrationChargePix &avpix = (MCalibrationChargePix&)(*fChargeCam)[i]
37//
38// Averaged values over one whole area index (e.g. inner or outer pixels for
39// the MAGIC camera), can be retrieved via:
40// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageArea(i)
41//
42// Averaged values over one whole camera sector can be retrieved via:
43// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageSector(i)
44//
45// Note the averageing has been done on an event-by-event basis. Resulting
46// Sigma's of the Gauss fit have been multiplied with the square root of the number
47// of involved pixels in order to make a direct comparison possible with the mean of
48// sigmas.
49//
50// Final numbers of uncalibrated or unreliable pixels can be retrieved via the commands:
51// GetNumUncalibrated(aidx) and GetNumUnreliable(aidx) where aidx is the area index (0 for
52// inner and 1 for outer pixels in the MAGIC camera).
53//
54// The following "calibration" constants are used for the calibration of each pixel
55// (see MCalibrate):
56//
57// - MCalibrationQEPix::GetMeanConvFADC2Phe(): The mean conversion factor from the
58// summed FADC slices to the number of photo-electrons (in first order independent
59// on colour and intensity)
60// - MCalibrationQEPix::GetMeanFFactorFADC2Phot(): The mean F-Factor of the total
61// readout chain dividing the signal to noise of the incoming number of photons
62// (= sqrt(number photons)) by the signal to noise of the outgoing summed FADC slices
63// signal (= MCalibrationPix::GetMean() / MCalibrationChargePix::GetRSigma() )
64//
65// The following calibration constants can be retrieved directly from this class:
66//
67// - GetConversionFactorFFactor ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
68//
69// where:
70// - idx is the pixel software ID
71// - "mean" is the mean conversion constant, to be multiplied with the retrieved signal
72// in order to get a calibrated number of PHOTONS.
73// - "err" is the pure statistical uncertainty about the MEAN
74// - "sigma", if mulitplied with the square root of signal, gives the approximate sigma of the
75// retrieved mean number of incident Cherenkov photons.
76//
77// Note, Conversion is ALWAYS (included the F-Factor method!) from summed FADC slices to PHOTONS.
78//
79// See also: MCalibrationChargePix, MCalibrationChargeCalc, MCalibrationQECam
80// MHCalibrationChargePix, MHCalibrationChargeCam
81// MCalibrationBlindPix, MCalibrationChargePINDiode
82//
83/////////////////////////////////////////////////////////////////////////////
84#include "MCalibrationChargeCam.h"
85#include "MCalibrationChargePix.h"
86
87#include <TMath.h>
88
89#include <TH1.h>
90#include <TF1.h>
91#include <TOrdCollection.h>
92
93#include "MLog.h"
94#include "MLogManip.h"
95
96#include "MGeomCam.h"
97#include "MGeom.h"
98
99#include "MBadPixelsCam.h"
100#include "MBadPixelsPix.h"
101
102#include "MCalibrationQECam.h"
103#include "MCalibrationQEPix.h"
104
105#include "MHCamera.h"
106
107ClassImp(MCalibrationChargeCam);
108
109using namespace std;
110// --------------------------------------------------------------------------
111//
112// Default constructor.
113//
114// Calls:
115// - Clear()
116//
117MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
118{
119 fName = name ? name : "MCalibrationChargeCam";
120 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
121
122 Clear();
123}
124
125/*
126// --------------------------------------------------------------------------
127//
128// Creates new MCalibrationCam only for the averaged areas:
129// the rest has to be retrieved directly, e.g. via:
130//
131TObject *MCalibrationChargeCam::Clone(const char *) const
132{
133
134 //
135 // FIXME, this might be done faster and more elegant, by direct copy.
136 //
137 MCalibrationChargeCam *cam = new MCalibrationChargeCam(fName,fTitle);
138
139 cam->fNumUnsuitable = fNumUnsuitable;
140 cam->fNumUnreliable = fNumUnreliable;
141 cam->fNumHiGainFADCSlices = fNumHiGainFADCSlices;
142 cam->fNumLoGainFADCSlices = fNumLoGainFADCSlices;
143 cam->fPulserColor = fPulserColor;
144
145 cam->fFlags = fFlags;
146
147 cam->fNumPhotonsBlindPixelMethod = fNumPhotonsBlindPixelMethod;
148 cam->fNumPhotonsFFactorMethod = fNumPhotonsFFactorMethod;
149 cam->fNumPhotonsPINDiodeMethod = fNumPhotonsPINDiodeMethod;
150 cam->fNumPhotonsBlindPixelMethodErr = fNumPhotonsBlindPixelMethodErr;
151 cam->fNumPhotonsFFactorMethodErr = fNumPhotonsFFactorMethodErr;
152 cam->fNumPhotonsPINDiodeMethodErr = fNumPhotonsPINDiodeMethodErr;
153
154 for (Int_t i=0; i<GetSize(); i++)
155 cam->fPixels->AddAt((*this)[i].Clone(),i);
156
157 for (Int_t i=0; i<GetAverageAreas(); i++)
158 {
159 cam->fAverageAreas->AddAt(GetAverageArea(i).Clone(),i);
160 cam->fAverageBadAreas->AddAt(GetAverageBadArea(i).Clone(),i);
161 }
162 for (Int_t i=0; i<GetAverageSectors(); i++)
163 {
164 cam->fAverageSectors->AddAt(GetAverageSector(i).Clone(),i);
165 cam->fAverageBadSectors->AddAt(GetAverageBadSector(i).Clone(),i);
166 }
167
168 return cam;
169}
170*/
171
172// -------------------------------------------------------------------
173//
174// Add MCalibrationChargePix's in the ranges from - to to fPixels
175//
176void MCalibrationChargeCam::Add(const UInt_t a, const UInt_t b)
177{
178 for (UInt_t i=a; i<b; i++)
179 {
180 fPixels->AddAt(new MCalibrationChargePix,i);
181 (*this)[i].SetPixId(i);
182 }
183}
184
185// -------------------------------------------------------------------
186//
187// Add MCalibrationChargePix's in the ranges from - to to fAverageAreas
188//
189void MCalibrationChargeCam::AddArea(const UInt_t a, const UInt_t b)
190{
191 for (UInt_t i=a; i<b; i++)
192 {
193 fAverageAreas->AddAt(new MCalibrationChargePix,i);
194 GetAverageArea(i).SetPixId(i);
195 }
196}
197
198// -------------------------------------------------------------------
199//
200// Add MCalibrationChargePix's in the ranges from - to to fAverageSectors
201//
202void MCalibrationChargeCam::AddSector(const UInt_t a, const UInt_t b)
203{
204 for (UInt_t i=a; i<b; i++)
205 {
206 fAverageSectors->AddAt(new MCalibrationChargePix,i);
207 GetAverageSector(i).SetPixId(i);
208 }
209}
210
211
212// --------------------------------------
213//
214// Sets all variable to 0.
215// Sets all flags to kFALSE
216// Calls MCalibrationCam::Clear()
217//
218void MCalibrationChargeCam::Clear(Option_t *o)
219{
220 // DO NOT RESET THE HI-/LO-GAIN CONSTANTS
221 SetFFactorMethodValid ( kFALSE );
222
223 fNumPhotonsBlindPixelMethod = 0.;
224 fNumPhotonsFFactorMethod = 0.;
225 fNumPhotonsPINDiodeMethod = 0.;
226 fNumPhotonsBlindPixelMethodErr = 0.;
227 fNumPhotonsFFactorMethodErr = 0.;
228 fNumPhotonsPINDiodeMethodErr = 0.;
229
230 MCalibrationCam::Clear();
231}
232
233// -----------------------------------------------
234//
235// Sets the kFFactorMethodValid bit from outside
236//
237void MCalibrationChargeCam::SetFFactorMethodValid(const Bool_t b)
238{
239 b ? SETBIT(fFlags, kFFactorMethodValid) : CLRBIT(fFlags, kFFactorMethodValid);
240}
241
242
243// --------------------------------------------------------------------------
244//
245// Test bit kFFactorMethodValid
246//
247Bool_t MCalibrationChargeCam::IsFFactorMethodValid() const
248{
249 return TESTBIT(fFlags,kFFactorMethodValid);
250}
251
252// --------------------------------------------------------------------------
253//
254// Copy High-gain vs. low-gain conversion factors from cam to this.
255//
256Bool_t MCalibrationChargeCam::MergeHiLoConversionFactors(const MCalibrationChargeCam &cam) const
257{
258
259 if (GetSize() != cam.GetSize())
260 {
261 *fLog << warn << "Sizes mismatch! Cannot merge high-gain vs. low-gain convertion factors" << endl;
262 return kFALSE;
263 }
264
265 for (Int_t i=0; i<GetSize(); i++)
266 {
267 ((MCalibrationChargePix&)(*this)[i]).SetConversionHiLo (((MCalibrationChargePix&)cam[i]).GetConversionHiLo());
268 ((MCalibrationChargePix&)(*this)[i]).SetConversionHiLoErr(((MCalibrationChargePix&)cam[i]).GetConversionHiLoErr());
269 }
270
271 return kTRUE;
272}
273
274// --------------------------------------------------------------------------
275//
276// Print first the well fitted pixels
277// and then the ones which are not FitValid
278//
279void MCalibrationChargeCam::Print(Option_t *o) const
280{
281
282 *fLog << all << GetDescriptor() << ":" << endl;
283 int id = 0;
284
285 *fLog << all << "Calibrated pixels:" << endl;
286 *fLog << all << endl;
287
288 TIter Next(fPixels);
289 MCalibrationChargePix *pix;
290 while ((pix=(MCalibrationChargePix*)Next()))
291 {
292
293 if (pix->IsFFactorMethodValid())
294 {
295
296 *fLog << all
297 << Form("%s%3i","Pixel: ",pix->GetPixId())
298 << Form("%s%4.2f%s%4.2f"," Ped.Rms: ",pix->GetPedRms(),"+-",pix->GetPedRmsErr())
299 << Form("%s%4.2f%s%4.2f"," Charge: " ,pix->GetConvertedMean(),"+-",pix->GetConvertedSigma())
300 << Form("%s%4.2f%s%4.2f"," Red.Sigma: ",pix->GetConvertedRSigma(),"+-",pix->GetConvertedRSigmaErr())
301 << Form("%s%4.2f%s%4.2f"," Num.Phes: ",pix->GetPheFFactorMethod(),"+-",pix->GetPheFFactorMethodErr())
302 << Form("%s%4.2f%s%4.2f"," Conv.FADC2Phe: ",pix->GetMeanConvFADC2Phe(),"+-",pix->GetMeanConvFADC2PheErr())
303 << " Saturated? :" << pix->IsHiGainSaturation()
304 << Form("%s%4.2f%s%4.2f"," Conv.HiLo: ",pix->GetConversionHiLo(),"+-",pix->GetConversionHiLoErr())
305 << endl;
306 id++;
307 }
308 }
309
310 *fLog << all << id << " pixels" << endl;
311 id = 0;
312
313
314 *fLog << all << endl;
315 *fLog << all << "Previously Excluded pixels:" << endl;
316 *fLog << all << endl;
317
318 id = 0;
319
320 TIter Next4(fPixels);
321 while ((pix=(MCalibrationChargePix*)Next4()))
322 {
323 if (pix->IsExcluded())
324 {
325 *fLog << all << pix->GetPixId() << " ";
326 id++;
327
328 if (!(id % 25))
329 *fLog << endl;
330 }
331 }
332
333 *fLog << all << endl;
334 *fLog << all << "New Excluded pixels:" << endl;
335 *fLog << all << endl;
336
337 TIter Next5(fPixels);
338 while ((pix=(MCalibrationChargePix*)Next5()))
339 {
340 if (!pix->IsFFactorMethodValid() && !pix->IsExcluded())
341 {
342 *fLog << all << pix->GetPixId() << " ";
343 id++;
344
345 if (!(id % 25))
346 *fLog << endl;
347 }
348 }
349
350 *fLog << endl;
351 *fLog << all << id << " Excluded pixels " << endl;
352 *fLog << endl;
353
354 *fLog << all << endl;
355 *fLog << all << "Averaged Areas:" << endl;
356 *fLog << all << endl;
357
358 TIter Next6(fAverageAreas);
359 while ((pix=(MCalibrationChargePix*)Next6()))
360 {
361 *fLog << all << Form("%s%3i","Area Idx: ",pix->GetPixId())
362 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
363 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
364 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
365 << " Reduced Sigma: " << pix->GetRSigma()
366 << " Nr Phe's: " << pix->GetPheFFactorMethod()
367 << endl;
368 }
369
370 *fLog << all << endl;
371 *fLog << all << "Averaged Sectors:" << endl;
372 *fLog << all << endl;
373
374 TIter Next7(fAverageSectors);
375 while ((pix=(MCalibrationChargePix*)Next7()))
376 {
377 *fLog << all << Form("%s%3i","Sector: ",pix->GetPixId())
378 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
379 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
380 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
381 << " Reduced Sigma: " << pix->GetRSigma()
382 << " Nr Phe's: " << pix->GetPheFFactorMethod()
383 << endl;
384 }
385 *fLog << all << endl;
386}
387
388
389// --------------------------------------------------------------------------
390//
391// The types are as follows:
392//
393// Fitted values:
394// ==============
395//
396// 0: Fitted Charge (see MCalibrationPix::GetMean())
397// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
398// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
399// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
400//
401// Useful variables derived from the fit results:
402// =============================================
403//
404// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
405// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
406// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
407// 7: Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerCharge())
408// 8: Error of Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerChargeErr())
409//
410// Results of the F-Factor calibration Method:
411// ===========================================
412//
413// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
414// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
415// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2Phe()
416// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2PheErr()
417// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2Phot()
418// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2PhotErr()
419// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
420//
421// Results of the Low-Gain vs. High-Gain Conversion:
422// =================================================
423//
424// 16: Mean Signal Hi-Gain / Mean Signal Lo-Gain (see MCalibrationPix::GetHiLoMeansDivided())
425// 17: Error Signal High-Gain / Signal Low Gain (see MCalibrationPix::GetHiLoMeansDividedErr())
426// 18: Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDivided())
427// 19: Error Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDividedErr())
428//
429// Localized defects:
430// ==================
431//
432// 20: Excluded Pixels
433// 21: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
434// 22: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
435// 23: Number of blackout events in the Hi Gain (see MCalibrationPix::GetHiGainNumBlackout())
436// 24: Number of blackout events in the Lo Gain (see MCalibrationPix::GetLoGainNumBlackout())
437//
438// Other classifications of pixels:
439// ================================
440//
441// 25: Pixels with saturated High-Gain (see MCalibrationPix::IsHiGainSaturation())
442//
443// Calculated absolute arrival times (very low precision!):
444// ========================================================
445//
446// 26: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
447// 27: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
448//
449// Used Pedestals:
450// ===============
451//
452// 28: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
453// 29: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
454// 30: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
455// 31: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
456//
457// Special variables (for data check):
458// ====================================
459//
460// 32: HiGain RMS divided by Mean for every pixel (with inclusion of the excluded pixels)
461//
462Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
463{
464 if (idx > GetSize())
465 return kFALSE;
466
467 const Float_t area = cam[idx].GetA();
468 if (area == 0)
469 return kFALSE;
470
471 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
472 switch (type)
473 {
474 case 0:
475 if (pix.IsExcluded() || pix.GetConvertedMeanErr()<0)
476 return kFALSE;
477 val = pix.GetConvertedMean();
478 return kTRUE;
479 case 1:
480 if (pix.IsExcluded() || pix.GetConvertedMeanErr()<0)
481 return kFALSE;
482 val = pix.GetConvertedMeanErr();
483 return kTRUE;
484 case 2:
485 if (pix.IsExcluded() || pix.GetConvertedSigmaErr()<0)
486 return kFALSE;
487 val = pix.GetConvertedSigma();
488 return kTRUE;
489 case 3:
490 if (pix.IsExcluded() || pix.GetConvertedSigmaErr()<0)
491 return kFALSE;
492 val = pix.GetConvertedSigmaErr();
493 return kTRUE;
494 case 4:
495 if (pix.IsExcluded())
496 return kFALSE;
497 val = pix.GetProb();
498 return val>=0;
499 case 5:
500 if (!pix.IsFFactorMethodValid())
501 return kFALSE;
502 if (pix.GetRSigma() == -1.)
503 return kFALSE;
504 val = pix.GetConvertedRSigma();
505 break;
506 case 6:
507 if (!pix.IsFFactorMethodValid())
508 return kFALSE;
509 if (pix.GetRSigma() == -1.)
510 return kFALSE;
511 val = pix.GetConvertedRSigmaErr();
512 break;
513 case 7:
514 if (!pix.IsFFactorMethodValid())
515 return kFALSE;
516 val = pix.GetRSigmaPerCharge();
517 break;
518 case 8:
519 if (!pix.IsFFactorMethodValid())
520 return kFALSE;
521 val = pix.GetRSigmaPerChargeErr();
522 break;
523 case 9:
524 // if (!pix.IsFFactorMethodValid())
525 // return kFALSE;
526 if (pix.IsExcluded() || pix.GetPheFFactorMethodErr()<0)
527 return kFALSE;
528 val = pix.GetPheFFactorMethod();
529 return kTRUE;
530 case 10:
531 // if (!pix.IsFFactorMethodValid())
532 // return kFALSE;
533 if (pix.IsExcluded() || pix.GetPheFFactorMethodErr()<=0)
534 return kFALSE;
535 val = pix.GetPheFFactorMethodErr();
536 return kTRUE;
537 case 11:
538 // if (!pix.IsFFactorMethodValid())
539 // return kFALSE;
540 if (pix.IsExcluded() || pix.GetMeanConvFADC2PheErr()<0)
541 return kFALSE;
542 val = pix.GetMeanConvFADC2Phe();
543 return kTRUE;
544 case 12:
545 // if (!pix.IsFFactorMethodValid())
546 // return kFALSE;
547 if (pix.IsExcluded() || pix.GetMeanConvFADC2PheErr()<=0)
548 return kFALSE;
549 val = pix.GetMeanConvFADC2PheErr();
550 return kTRUE;
551 case 13:
552 // if (!pix.IsFFactorMethodValid())
553 // return kFALSE;
554 if (pix.IsExcluded())
555 return kFALSE;
556 val = pix.GetMeanFFactorFADC2Phot();
557 break;
558 case 14:
559 // if (!pix.IsFFactorMethodValid())
560 // return kFALSE;
561 if (pix.IsExcluded())
562 return kFALSE;
563 val = pix.GetMeanFFactorFADC2PhotErr();
564 if (val <= 0.)
565 val = 0.00001;
566 break;
567 case 15:
568 if (pix.IsExcluded())
569 return kFALSE;
570 if (!pix.IsFFactorMethodValid())
571 return kFALSE;
572 val = 1;
573 return kTRUE;
574 case 16:
575 if (pix.IsExcluded())
576 return kFALSE;
577 val = pix.GetHiLoMeansDivided();
578 break;
579 case 17:
580 if (pix.IsExcluded())
581 return kFALSE;
582 val = pix.GetHiLoMeansDividedErr();
583 break;
584 case 18:
585 if (pix.IsExcluded())
586 return kFALSE;
587 val = pix.GetHiLoSigmasDivided();
588 break;
589 case 19:
590 if (pix.IsExcluded())
591 return kFALSE;
592 val = pix.GetHiLoSigmasDividedErr();
593 break;
594 case 20:
595 if (!pix.IsExcluded())
596 return kFALSE;
597 val = 1.;
598 return kTRUE;
599 case 21:
600 if (pix.IsExcluded())
601 return kFALSE;
602 val = pix.GetHiGainNumPickup();
603 break;
604 case 22:
605 if (pix.IsExcluded())
606 return kFALSE;
607 val = pix.GetLoGainNumPickup();
608 break;
609 case 23:
610 if (pix.IsExcluded())
611 return kFALSE;
612 val = pix.GetHiGainNumBlackout();
613 break;
614 case 24:
615 if (pix.IsExcluded())
616 return kFALSE;
617 val = pix.GetLoGainNumBlackout();
618 break;
619 case 25:
620 if (pix.IsExcluded())
621 return kFALSE;
622 val = pix.IsHiGainSaturation();
623 break;
624 case 26:
625 if (pix.IsExcluded())
626 return kFALSE;
627 val = pix.GetAbsTimeMean();
628 break;
629 case 27:
630 if (pix.IsExcluded())
631 return kFALSE;
632 val = pix.GetAbsTimeRms();
633 break;
634 case 28:
635 if (pix.IsExcluded())
636 return kFALSE;
637 val = pix.GetPed();
638 break;
639 case 29:
640 if (pix.IsExcluded())
641 return kFALSE;
642 val = pix.GetPedErr();
643 break;
644 case 30:
645 if (pix.IsExcluded())
646 return kFALSE;
647 val = pix.GetPedRms();
648 break;
649 case 31:
650 if (pix.IsExcluded())
651 return kFALSE;
652 val = pix.GetPedErr()/2.;
653 break;
654 case 32:
655 if (pix.IsExcluded() || pix.GetRms()<0)
656 return kFALSE;
657 val = pix.GetMean() == 0. ? 0. : pix.GetRms()/pix.GetMean();
658 return kTRUE;
659 case 33:
660 if (pix.IsExcluded() || pix.GetRms()<0)
661 return kFALSE;
662 if (pix.GetMean() == 0.)
663 val = 0.;
664 else
665 val = pix.GetSigmaErr()/pix.GetMean() + pix.GetRms()*pix.GetMeanErr()/pix.GetMean()/pix.GetMean();
666 return kTRUE;
667 default:
668 return kFALSE;
669 }
670
671 return val!=-1.;
672}
673
674Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &ferr, Float_t &ffactor)
675{
676
677 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[ipx];
678
679 Float_t conv = pix.GetMeanConvFADC2Phe();
680
681 if (conv < 0.)
682 return kFALSE;
683
684 mean = conv;
685 ferr = pix.GetMeanConvFADC2PheErr();
686 ffactor = pix.GetMeanFFactorFADC2Phot();
687
688 return kTRUE;
689}
690
691
692// --------------------------------------------------------------------------
693//
694// Calculates the average conversion factor FADC counts to photons for pixel sizes.
695// The geometry container is used to get the necessary
696// geometry information (area index). The MCalibrationQECam container is necessary for
697// the quantum efficiency information.
698// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
699// in the calculation of the size average.
700//
701// Returns a TArrayF of dimension 2:
702// arr[0]: averaged conversion factors (default: -1.)
703// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
704//
705TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
706 const UInt_t ai, MBadPixelsCam *bad)
707{
708
709 const Int_t np = GetSize();
710
711 Double_t mean = 0.;
712 Double_t mean2 = 0.;
713 Int_t nr = 0;
714
715 MHCamera convcam(geom,"ConvFactors","Conversion Factors;Conv Factor [phot/FADC cnts];channels");
716
717 for (int i=0; i<np; i++)
718 {
719 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
720 continue;
721
722 if (bad && (*bad)[i].IsUncalibrated(MBadPixelsPix::kDeviatingNumPhes))
723 continue;
724
725 const UInt_t aidx = geom[i].GetAidx();
726
727 if (ai != aidx)
728 continue;
729
730 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
731 if (!qepix.IsAverageQEFFactorAvailable())
732 continue;
733
734 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
735 const Float_t conv = pix.GetMeanConvFADC2Phe()/qepix.GetQECascadesFFactor();
736
737 mean += conv;
738 mean2 += conv*conv;
739 nr ++;
740
741 convcam.Fill(i,conv);
742 convcam.SetUsed(i);
743 }
744
745 Float_t mn = nr ? mean/nr : -1.;
746 Float_t sg = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0.;
747
748 const Int_t aidx = (Int_t)ai;
749
750 TH1D *h = convcam.ProjectionS(TArrayI(),TArrayI(1,&aidx),"_py",750);
751 h->SetDirectory(NULL);
752
753 TF1 *fit = NULL;
754
755 if (geom.InheritsFrom("MGeomCamMagic"))
756 {
757
758 fit = new TF1("fit","gaus",0.4,5.);
759
760 // Fix the ranges, as found by Nadia
761 if(aidx == 0)
762 h->Fit(fit, "REQ0", "",0.4,1.5);
763 else
764 h->Fit(fit ,"REQ0", "",1.,5.);
765 }
766 else
767 {
768 h->Fit("gaus","Q0");
769 fit = h->GetFunction("gaus");
770 }
771
772 Float_t ci2 = fit->GetChisquare();
773 Float_t sigma = fit->GetParameter(2);
774
775 if (ci2 > 500. || sigma > sg)
776 {
777 if (geom.InheritsFrom("MGeomCamMagic"))
778 {
779 // Fix the ranges, as found by Nadia
780 if(aidx == 0)
781 h->Fit(fit, "REQ0", "",0.4,1.5);
782 else
783 h->Fit(fit, "REQ0", "",1.,5.);
784 }
785 else
786 {
787 h->Fit("gaus","MREQ0");
788 fit = h->GetFunction("gaus");
789 }
790
791 ci2 = fit->GetChisquare();
792 sigma = fit->GetParameter(2);
793 }
794
795 const Int_t ndf = fit->GetNDF();
796
797 if (ci2 < 500. && sigma < sg && ndf > 2)
798 {
799 mn = fit->GetParameter(1);
800 sg = sigma;
801 }
802
803 *fLog << inf << "Conversion Factors to photons area idx: " << aidx << ":" << endl;
804 *fLog << inf << "Mean: " << Form("%4.3f",mn)
805 << "+-" << Form("%4.3f",fit->GetParError(1))
806 << " Sigma: " << Form("%4.3f",sg) << "+-" << Form("%4.3f",fit->GetParError(2))
807 << " Chisquare: " << Form("%4.3f",fit->GetChisquare()) << " NDF : " << ndf << endl;
808
809 delete fit;
810 delete h;
811 gROOT->GetListOfFunctions()->Remove(fit);
812
813 TArrayF arr(2);
814 arr[0] = mn;
815 arr[1] = sg;
816
817 return arr;
818}
819
820// --------------------------------------------------------------------------
821//
822// Calculates the average conversion factor FADC counts to equiv. photo-electrons for pixel sizes.
823// The geometry container is used to get the necessary
824// geometry information (area index). The MCalibrationQECam container is necessary for
825// the quantum efficiency information.
826// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
827// in the calculation of the size average.
828//
829// Returns a TArrayF of dimension 2:
830// arr[0]: averaged conversion factors (default: -1.)
831// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
832//
833TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhePerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
834 const UInt_t ai, MBadPixelsCam *bad)
835{
836
837 const Int_t np = GetSize();
838
839 Double_t mean = 0.;
840 Double_t mean2 = 0.;
841 Int_t nr = 0;
842
843 MHCamera convcam(geom,"ConvFactors","Conversion Factors;Conv Factor [phe/FADC cnts];channels");
844
845 for (int i=0; i<np; i++)
846 {
847 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
848 continue;
849
850 if (bad && (*bad)[i].IsUncalibrated(MBadPixelsPix::kDeviatingNumPhes))
851 continue;
852
853 const UInt_t aidx = geom[i].GetAidx();
854
855 if (ai != aidx)
856 continue;
857
858 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
859 if (!qepix.IsAverageQEFFactorAvailable())
860 continue;
861
862 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
863 const Float_t conv = pix.GetMeanConvFADC2Phe()/qepix.GetQECascadesFFactor()*MCalibrationQEPix::gkDefaultAverageQE;
864
865 mean += conv;
866 mean2 += conv*conv;
867 nr ++;
868
869 convcam.Fill(i,conv);
870 convcam.SetUsed(i);
871 }
872
873 Float_t mn = nr ? mean/nr : -1.;
874 Float_t sg = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0.;
875
876 const Int_t aidx = (Int_t)ai;
877
878 TH1D *h = convcam.ProjectionS(TArrayI(),TArrayI(1,&aidx),"_py",750);
879 h->SetDirectory(NULL);
880
881 TF1 *fit = NULL;
882
883 if (geom.InheritsFrom("MGeomCamMagic"))
884 {
885
886 fit = new TF1("fit","gaus",0.01,1.);
887
888 // Fix the ranges, as found by Nadia
889 if(aidx == 0)
890 {h->Fit("fit","REQ0", "",0.07,0.3);}
891 else
892 {h->Fit("fit","REQ0", "",0.15,1.0);}
893 }
894 else
895 {
896 h->Fit("gaus","Q0");
897 fit = h->GetFunction("gaus");
898 }
899
900 Float_t ci2 = fit->GetChisquare();
901 Float_t sigma = fit->GetParameter(2);
902
903 if (ci2 > 500. || sigma > sg)
904 {
905 if (geom.InheritsFrom("MGeomCamMagic"))
906 {
907 // Fix the ranges, as found by Nadia
908 if(aidx == 0)
909 {h->Fit("fit","REQ0", "",0.07,0.3);}
910 else
911 {h->Fit("fit","REQ0", "",0.15,1.0);}
912 }
913 else
914 {
915 h->Fit("gaus","MREQ0");
916 fit = h->GetFunction("gaus");
917 }
918
919 ci2 = fit->GetChisquare();
920 sigma = fit->GetParameter(2);
921 }
922
923 const Int_t ndf = fit->GetNDF();
924
925 if (ci2 < 500. && sigma < sg && ndf > 2)
926 {
927 mn = fit->GetParameter(1);
928 sg = sigma;
929 }
930
931 *fLog << inf << "Conversion Factors to equiv. photo-electrons area idx: " << aidx << ":" << endl;
932 *fLog << inf << "Mean: " << Form("%4.3f",mn)
933 << "+-" << Form("%4.3f",fit->GetParError(1))
934 << " Sigma: " << Form("%4.3f",sg) << "+-" << Form("%4.3f",fit->GetParError(2))
935 << " Chisquare: " << Form("%4.3f",fit->GetChisquare()) << " NDF : " << ndf << endl;
936
937 delete fit;
938 delete h;
939 gROOT->GetListOfFunctions()->Remove(fit);
940
941 TArrayF arr(2);
942 arr[0] = mn;
943 arr[1] = sg;
944
945 return arr;
946}
947
948// --------------------------------------------------------------------------
949//
950// Calculates the average conversion factor FADC counts to photons for camera sectors.
951// The geometry container is used to get the necessary
952// geometry information (area index). The MCalibrationQECam container is necessary for
953// the quantum efficiency information.
954// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
955// in the calculation of the size average.
956//
957// Returns a TArrayF of dimension 2:
958// arr[0]: averaged conversion factors (default: -1.)
959// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
960//
961TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam,
962 const UInt_t sec, MBadPixelsCam *bad)
963{
964 const Int_t np = GetSize();
965
966 Double_t mean = 0.;
967 Double_t mean2 = 0.;
968 Int_t nr = 0;
969
970 for (int i=0; i<np; i++)
971 {
972 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
973 continue;
974
975 const UInt_t sector = geom[i].GetSector();
976
977 if (sector != sec)
978 continue;
979
980 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
981 if (!qepix.IsAverageQEFFactorAvailable())
982 continue;
983
984 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
985 const Float_t conv = pix.GetMeanConvFADC2Phe();
986 const Float_t qe = qepix.GetQECascadesFFactor();
987
988 mean += conv/qe;
989 mean2 += conv*conv/qe/qe;
990 nr ++;
991
992 }
993
994 TArrayF arr(2);
995 arr[0] = nr ? mean/nr : -1;
996 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
997 return arr;
998}
999
1000// --------------------------------------------------------------------------
1001//
1002// Calculates the average mean arrival times for pixel sizes.
1003// The geometry container is used to get the necessary
1004// geometry information (area index).
1005// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1006// in the calculation of the size average.
1007//
1008// Returns a TArrayF of dimension 2:
1009// arr[0]: averaged mean arrival times (default: -1.)
1010// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
1011//
1012TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerArea(const MGeomCam &geom,
1013 const UInt_t ai, MBadPixelsCam *bad)
1014{
1015
1016 const Int_t np = GetSize();
1017
1018 Double_t mean = 0.;
1019 Double_t mean2 = 0.;
1020 Int_t nr = 0;
1021
1022 for (int i=0; i<np; i++)
1023 {
1024 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1025 continue;
1026
1027 const UInt_t aidx = geom[i].GetAidx();
1028
1029 if (ai != aidx)
1030 continue;
1031
1032 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1033 const Float_t time = pix.GetAbsTimeMean();
1034
1035 mean += time ;
1036 mean2 += time*time;
1037 nr ++;
1038
1039 }
1040
1041 TArrayF arr(2);
1042 arr[0] = nr ? mean/nr : -1;
1043 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1044 return arr;
1045}
1046
1047// --------------------------------------------------------------------------
1048//
1049// Calculates the average mean arrival times for camera sectors.
1050// The geometry container is used to get the necessary
1051// geometry information (area index).
1052// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1053// in the calculation of the size average.
1054//
1055// Returns a TArrayF of dimension 2:
1056// arr[0]: averaged mean arrival times (default: -1.)
1057// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
1058//
1059TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerSector(const MGeomCam &geom,
1060 const UInt_t sec, MBadPixelsCam *bad)
1061{
1062 const Int_t np = GetSize();
1063
1064 Double_t mean = 0.;
1065 Double_t mean2 = 0.;
1066 Int_t nr = 0;
1067
1068 for (int i=0; i<np; i++)
1069 {
1070 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1071 continue;
1072
1073 const UInt_t sector = geom[i].GetSector();
1074
1075 if (sector != sec)
1076 continue;
1077
1078 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1079 const Float_t time = pix.GetAbsTimeMean();
1080
1081 mean += time;
1082 mean2 += time*time;
1083 nr ++;
1084
1085 }
1086
1087 TArrayF arr(2);
1088 arr[0] = nr ? mean/nr : -1;
1089 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1090 return arr;
1091}
1092
1093// --------------------------------------------------------------------------
1094//
1095// Calculates the average arrival time RMSs for pixel sizes.
1096// The geometry container is used to get the necessary
1097// geometry information (area index).
1098// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1099// in the calculation of the size average.
1100//
1101// Returns a TArrayF of dimension 2:
1102// arr[0]: averaged arrival time RMSs (default: -1.)
1103// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
1104//
1105TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerArea ( const MGeomCam &geom,
1106 const UInt_t ai, MBadPixelsCam *bad)
1107{
1108
1109 const Int_t np = GetSize();
1110
1111 Double_t mean = 0.;
1112 Double_t mean2 = 0.;
1113 Int_t nr = 0;
1114
1115 for (int i=0; i<np; i++)
1116 {
1117 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1118 continue;
1119
1120 const UInt_t aidx = geom[i].GetAidx();
1121
1122 if (ai != aidx)
1123 continue;
1124
1125 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1126 const Float_t rms = pix.GetAbsTimeRms();
1127
1128 mean += rms;
1129 mean2 += rms*rms;
1130 nr ++;
1131
1132 }
1133
1134 TArrayF arr(2);
1135 arr[0] = nr ? mean/nr : -1;
1136 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1137 return arr;
1138}
1139
1140// --------------------------------------------------------------------------
1141//
1142// Calculates the average arrival time RMSs for camera sectors.
1143// The geometry container is used to get the necessary
1144// geometry information (area index).
1145// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1146// in the calculation of the size average.
1147//
1148// Returns a TArrayF of dimension 2:
1149// arr[0]: averaged arrival time RMSs (default: -1.)
1150// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
1151//
1152TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerSector( const MGeomCam &geom, const UInt_t sec, MBadPixelsCam *bad)
1153{
1154 const Int_t np = GetSize();
1155
1156 Double_t mean = 0.;
1157 Double_t mean2 = 0.;
1158 Int_t nr = 0;
1159
1160 for (int i=0; i<np; i++)
1161 {
1162 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1163 continue;
1164
1165 const UInt_t sector = geom[i].GetSector();
1166
1167 if (sector != sec)
1168 continue;
1169
1170 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1171 const Float_t rms = pix.GetAbsTimeRms();
1172
1173 mean += rms;
1174 mean2 += rms*rms;
1175 nr ++;
1176 }
1177
1178 TArrayF arr(2);
1179 arr[0] = nr ? mean/nr : -1;
1180 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1181
1182 return arr;
1183}
Note: See TracBrowser for help on using the repository browser.