source: trunk/MagicSoft/Mars/mcalib/MCalibrationChargeCam.cc@ 5557

Last change on this file since 5557 was 5557, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 28.7 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 11/2003 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24/////////////////////////////////////////////////////////////////////////////
25//
26// MCalibrationChargeCam
27//
28// Storage container for charge calibration results from the signal distribution
29// fits (see MHCalibrationChargeCam and MHCalibrationChargePix), the calculation
30// of reduced sigmas and number of photo-electrons (this class) and conversion
31// factors sum FADC slices to photo-electrons (see MCalibrationChargeCalc)
32//
33// Individual pixels have to be cast when retrieved e.g.:
34// MCalibrationChargePix &avpix = (MCalibrationChargePix&)(*fChargeCam)[i]
35//
36// Averaged values over one whole area index (e.g. inner or outer pixels for
37// the MAGIC camera), can be retrieved via:
38// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageArea(i)
39//
40// Averaged values over one whole camera sector can be retrieved via:
41// MCalibrationChargePix &avpix = (MCalibrationChargePix&)fChargeCam->GetAverageSector(i)
42//
43// Note the averageing has been done on an event-by-event basis. Resulting
44// Sigma's of the Gauss fit have been multiplied with the square root of the number
45// of involved pixels in order to make a direct comparison possible with the mean of
46// sigmas.
47//
48// Final numbers of uncalibrated or unreliable pixels can be retrieved via the commands:
49// GetNumUncalibrated(aidx) and GetNumUnreliable(aidx) where aidx is the area index (0 for
50// inner and 1 for outer pixels in the MAGIC camera).
51//
52// The following "calibration" constants are used for the calibration of each pixel
53// (see MCalibrate):
54//
55// - MCalibrationQEPix::GetMeanConvFADC2Phe(): The mean conversion factor from the
56// summed FADC slices to the number of photo-electrons (in first order independent
57// on colour and intensity)
58// - MCalibrationQEPix::GetMeanFFactorFADC2Phot(): The mean F-Factor of the total
59// readout chain dividing the signal to noise of the incoming number of photons
60// (= sqrt(number photons)) by the signal to noise of the outgoing summed FADC slices
61// signal (= MCalibrationPix::GetMean() / MCalibrationChargePix::GetRSigma() )
62//
63// The following calibration constants can be retrieved directly from this class:
64//
65// - GetConversionFactorFFactor ( Int_t idx, Float_t &mean, Float_t &err, Float_t &sigma );
66//
67// where:
68// - idx is the pixel software ID
69// - "mean" is the mean conversion constant, to be multiplied with the retrieved signal
70// in order to get a calibrated number of PHOTONS.
71// - "err" is the pure statistical uncertainty about the MEAN
72// - "sigma", if mulitplied with the square root of signal, gives the approximate sigma of the
73// retrieved mean number of incident Cherenkov photons.
74//
75// Note, Conversion is ALWAYS (included the F-Factor method!) from summed FADC slices to PHOTONS.
76//
77// See also: MCalibrationChargePix, MCalibrationChargeCalc, MCalibrationQECam
78// MHCalibrationChargePix, MHCalibrationChargeCam
79// MCalibrationBlindPix, MCalibrationChargePINDiode
80//
81/////////////////////////////////////////////////////////////////////////////
82#include "MCalibrationChargeCam.h"
83#include "MCalibrationChargePix.h"
84
85#include <TOrdCollection.h>
86
87#include "MLog.h"
88#include "MLogManip.h"
89
90#include "MGeomCam.h"
91#include "MGeomPix.h"
92
93#include "MBadPixelsCam.h"
94#include "MBadPixelsPix.h"
95
96#include "MCalibrationQECam.h"
97#include "MCalibrationQEPix.h"
98
99ClassImp(MCalibrationChargeCam);
100
101using namespace std;
102// --------------------------------------------------------------------------
103//
104// Default constructor.
105//
106// Calls:
107// - Clear()
108//
109MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
110{
111 fName = name ? name : "MCalibrationChargeCam";
112 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
113
114 Clear();
115}
116
117/*
118// --------------------------------------------------------------------------
119//
120// Creates new MCalibrationCam only for the averaged areas:
121// the rest has to be retrieved directly, e.g. via:
122//
123TObject *MCalibrationChargeCam::Clone(const char *) const
124{
125
126 //
127 // FIXME, this might be done faster and more elegant, by direct copy.
128 //
129 MCalibrationChargeCam *cam = new MCalibrationChargeCam(fName,fTitle);
130
131 cam->fNumUnsuitable = fNumUnsuitable;
132 cam->fNumUnreliable = fNumUnreliable;
133 cam->fNumHiGainFADCSlices = fNumHiGainFADCSlices;
134 cam->fNumLoGainFADCSlices = fNumLoGainFADCSlices;
135 cam->fPulserColor = fPulserColor;
136
137 cam->fFlags = fFlags;
138
139 cam->fNumPhotonsBlindPixelMethod = fNumPhotonsBlindPixelMethod;
140 cam->fNumPhotonsFFactorMethod = fNumPhotonsFFactorMethod;
141 cam->fNumPhotonsPINDiodeMethod = fNumPhotonsPINDiodeMethod;
142 cam->fNumPhotonsBlindPixelMethodErr = fNumPhotonsBlindPixelMethodErr;
143 cam->fNumPhotonsFFactorMethodErr = fNumPhotonsFFactorMethodErr;
144 cam->fNumPhotonsPINDiodeMethodErr = fNumPhotonsPINDiodeMethodErr;
145
146 for (Int_t i=0; i<GetSize(); i++)
147 cam->fPixels->AddAt((*this)[i].Clone(),i);
148
149 for (Int_t i=0; i<GetAverageAreas(); i++)
150 {
151 cam->fAverageAreas->AddAt(GetAverageArea(i).Clone(),i);
152 cam->fAverageBadAreas->AddAt(GetAverageBadArea(i).Clone(),i);
153 }
154 for (Int_t i=0; i<GetAverageSectors(); i++)
155 {
156 cam->fAverageSectors->AddAt(GetAverageSector(i).Clone(),i);
157 cam->fAverageBadSectors->AddAt(GetAverageBadSector(i).Clone(),i);
158 }
159
160 return cam;
161}
162*/
163
164// -------------------------------------------------------------------
165//
166// Add MCalibrationChargePix's in the ranges from - to to fPixels
167//
168void MCalibrationChargeCam::Add(const UInt_t a, const UInt_t b)
169{
170 for (UInt_t i=a; i<b; i++)
171 fPixels->AddAt(new MCalibrationChargePix,i);
172}
173
174// -------------------------------------------------------------------
175//
176// Add MCalibrationChargePix's in the ranges from - to to fAverageAreas
177//
178void MCalibrationChargeCam::AddArea(const UInt_t a, const UInt_t b)
179{
180 for (UInt_t i=a; i<b; i++)
181 fAverageAreas->AddAt(new MCalibrationChargePix,i);
182}
183
184// -------------------------------------------------------------------
185//
186// Add MCalibrationChargePix's in the ranges from - to to fAverageSectors
187//
188void MCalibrationChargeCam::AddSector(const UInt_t a, const UInt_t b)
189{
190 for (UInt_t i=a; i<b; i++)
191 fAverageSectors->AddAt(new MCalibrationChargePix,i);
192}
193
194
195// --------------------------------------
196//
197// Sets all variable to 0.
198// Sets all flags to kFALSE
199// Calls MCalibrationCam::Clear()
200//
201void MCalibrationChargeCam::Clear(Option_t *o)
202{
203
204 SetFFactorMethodValid ( kFALSE );
205
206 fNumPhotonsBlindPixelMethod = 0.;
207 fNumPhotonsFFactorMethod = 0.;
208 fNumPhotonsPINDiodeMethod = 0.;
209 fNumPhotonsBlindPixelMethodErr = 0.;
210 fNumPhotonsFFactorMethodErr = 0.;
211 fNumPhotonsPINDiodeMethodErr = 0.;
212
213 MCalibrationCam::Clear();
214
215 return;
216}
217
218// -----------------------------------------------
219//
220// Sets the kFFactorMethodValid bit from outside
221//
222void MCalibrationChargeCam::SetFFactorMethodValid(const Bool_t b)
223{
224 b ? SETBIT(fFlags, kFFactorMethodValid) : CLRBIT(fFlags, kFFactorMethodValid);
225}
226
227
228// --------------------------------------------------------------------------
229//
230// Test bit kFFactorMethodValid
231//
232Bool_t MCalibrationChargeCam::IsFFactorMethodValid() const
233{
234 return TESTBIT(fFlags,kFFactorMethodValid);
235}
236
237// --------------------------------------------------------------------------
238//
239// Print first the well fitted pixels
240// and then the ones which are not FitValid
241//
242void MCalibrationChargeCam::Print(Option_t *o) const
243{
244
245 *fLog << all << GetDescriptor() << ":" << endl;
246 int id = 0;
247
248 *fLog << all << "Calibrated pixels:" << endl;
249 *fLog << all << endl;
250
251 TIter Next(fPixels);
252 MCalibrationChargePix *pix;
253 while ((pix=(MCalibrationChargePix*)Next()))
254 {
255
256 if (!pix->IsExcluded())
257 {
258
259 *fLog << all
260 << Form("%s%3i","Pixel: ",pix->GetPixId())
261 << Form("%s%4.2f%s%4.2f"," Ped.Rms: ",pix->GetPedRms(),"+-",pix->GetPedRmsErr())
262 << Form("%s%4.2f%s%4.2f"," Charge: " ,pix->GetConvertedMean(),"+-",pix->GetConvertedSigma())
263 << Form("%s%4.2f%s%4.2f"," Red.Sigma: ",pix->GetConvertedRSigma(),"+-",pix->GetConvertedRSigmaErr())
264 << Form("%s%4.2f%s%4.2f"," Num.Phes: ",pix->GetPheFFactorMethod(),"+-",pix->GetPheFFactorMethodErr())
265 << Form("%s%4.2f%s%4.2f"," Conv.FADC2Phe: ",pix->GetMeanConvFADC2Phe(),"+-",pix->GetMeanConvFADC2PheErr())
266 << " Saturated? :" << pix->IsHiGainSaturation()
267 << endl;
268 id++;
269 }
270 }
271
272 *fLog << all << id << " pixels" << endl;
273 id = 0;
274
275
276 *fLog << all << endl;
277 *fLog << all << "Excluded pixels:" << endl;
278 *fLog << all << endl;
279
280 id = 0;
281
282 TIter Next4(fPixels);
283 while ((pix=(MCalibrationChargePix*)Next4()))
284 {
285 if (pix->IsExcluded())
286 {
287 *fLog << all << pix->GetPixId() << " ";
288 id++;
289
290 if (!(id % 25))
291 *fLog << endl;
292 }
293 }
294
295 *fLog << endl;
296 *fLog << all << id << " Excluded pixels " << endl;
297 *fLog << endl;
298
299 *fLog << all << endl;
300 *fLog << all << "Averaged Areas:" << endl;
301 *fLog << all << endl;
302
303 TIter Next5(fAverageAreas);
304 while ((pix=(MCalibrationChargePix*)Next5()))
305 {
306 *fLog << all << Form("%s%3i","Area Idx: ",pix->GetPixId())
307 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
308 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
309 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
310 << " Reduced Sigma: " << pix->GetRSigma()
311 << " Nr Phe's: " << pix->GetPheFFactorMethod()
312 << endl;
313 }
314
315 *fLog << all << endl;
316 *fLog << all << "Averaged Sectors:" << endl;
317 *fLog << all << endl;
318
319 TIter Next6(fAverageSectors);
320 while ((pix=(MCalibrationChargePix*)Next6()))
321 {
322 *fLog << all << Form("%s%3i","Sector: ",pix->GetPixId())
323 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
324 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
325 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
326 << " Reduced Sigma: " << pix->GetRSigma()
327 << " Nr Phe's: " << pix->GetPheFFactorMethod()
328 << endl;
329 }
330 *fLog << all << endl;
331}
332
333
334// --------------------------------------------------------------------------
335//
336// The types are as follows:
337//
338// Fitted values:
339// ==============
340//
341// 0: Fitted Charge (see MCalibrationPix::GetMean())
342// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
343// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
344// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
345//
346// Useful variables derived from the fit results:
347// =============================================
348//
349// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
350// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
351// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
352// 7: Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerCharge())
353// 8: Error of Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerChargeErr())
354//
355// Results of the F-Factor calibration Method:
356// ===========================================
357//
358// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
359// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
360// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2Phe()
361// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2PheErr()
362// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2Phot()
363// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2PhotErr()
364// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
365//
366// Results of the Low-Gain vs. High-Gain Conversion:
367// =================================================
368//
369// 16: Mean Signal Hi-Gain / Mean Signal Lo-Gain (see MCalibrationPix::GetHiLoMeansDivided())
370// 17: Error Signal High-Gain / Signal Low Gain (see MCalibrationPix::GetHiLoMeansDividedErr())
371// 18: Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDivided())
372// 19: Error Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDividedErr())
373//
374// Localized defects:
375// ==================
376//
377// 20: Excluded Pixels
378// 21: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
379// 22: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
380// 23: Number of blackout events in the Hi Gain (see MCalibrationPix::GetHiGainNumBlackout())
381// 24: Number of blackout events in the Lo Gain (see MCalibrationPix::GetLoGainNumBlackout())
382//
383// Other classifications of pixels:
384// ================================
385//
386// 25: Pixels with saturated High-Gain (see MCalibrationPix::IsHiGainSaturation())
387//
388// Calculated absolute arrival times (very low precision!):
389// ========================================================
390//
391// 26: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
392// 27: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
393//
394// Used Pedestals:
395// ===============
396//
397// 28: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
398// 29: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
399// 30: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
400// 31: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
401//
402// Special variables (for data check):
403// ====================================
404//
405// 32: HiGain RMS divided by Mean for every pixel (with inclusion of the excluded pixels)
406//
407Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
408{
409
410 if (idx > GetSize())
411 return kFALSE;
412
413 Float_t area = cam[idx].GetA();
414
415 if (area == 0)
416 return kFALSE;
417
418 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
419
420 switch (type)
421 {
422 case 0:
423 if (pix.IsExcluded())
424 return kFALSE;
425 val = pix.GetMean();
426 break;
427 case 1:
428 if (pix.IsExcluded())
429 return kFALSE;
430 val = pix.GetMeanErr();
431 break;
432 case 2:
433 if (pix.IsExcluded())
434 return kFALSE;
435 val = pix.GetSigma();
436 break;
437 case 3:
438 if (pix.IsExcluded())
439 return kFALSE;
440 val = pix.GetSigmaErr();
441 break;
442 case 4:
443 if (pix.IsExcluded())
444 return kFALSE;
445 val = pix.GetProb();
446 break;
447 case 5:
448 if (pix.IsExcluded())
449 return kFALSE;
450 if (pix.GetRSigma() == -1.)
451 return kFALSE;
452 val = pix.GetRSigma();
453 break;
454 case 6:
455 if (pix.IsExcluded())
456 return kFALSE;
457 if (pix.GetRSigma() == -1.)
458 return kFALSE;
459 val = pix.GetRSigmaErr();
460 break;
461 case 7:
462 if (pix.IsExcluded())
463 return kFALSE;
464 val = pix.GetRSigmaPerCharge();
465 break;
466 case 8:
467 if (pix.IsExcluded())
468 return kFALSE;
469 val = pix.GetRSigmaPerChargeErr();
470 break;
471 case 9:
472 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
473 return kFALSE;
474 val = pix.GetPheFFactorMethod();
475 break;
476 case 10:
477 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
478 return kFALSE;
479 val = pix.GetPheFFactorMethodErr();
480 break;
481 case 11:
482 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
483 return kFALSE;
484 val = pix.GetMeanConvFADC2Phe();
485 break;
486 case 12:
487 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
488 return kFALSE;
489 val = pix.GetMeanConvFADC2PheErr();
490 break;
491 case 13:
492 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
493 return kFALSE;
494 val = pix.GetMeanFFactorFADC2Phot();
495 break;
496 case 14:
497 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
498 return kFALSE;
499 val = pix.GetMeanFFactorFADC2PhotErr();
500 break;
501 case 15:
502 if (pix.IsExcluded())
503 return kFALSE;
504 if (pix.IsFFactorMethodValid())
505 val = 1;
506 else
507 return kFALSE;
508 break;
509 case 16:
510 if (pix.IsExcluded())
511 return kFALSE;
512 val = pix.GetHiLoMeansDivided();
513 break;
514 case 17:
515 if (pix.IsExcluded())
516 return kFALSE;
517 val = pix.GetHiLoMeansDividedErr();
518 break;
519 case 18:
520 if (pix.IsExcluded())
521 return kFALSE;
522 val = pix.GetHiLoSigmasDivided();
523 break;
524 case 19:
525 if (pix.IsExcluded())
526 return kFALSE;
527 val = pix.GetHiLoSigmasDividedErr();
528 break;
529 case 20:
530 if (pix.IsExcluded())
531 val = 1.;
532 else
533 return kFALSE;
534 break;
535 case 21:
536 if (pix.IsExcluded())
537 return kFALSE;
538 val = pix.GetHiGainNumPickup();
539 break;
540 case 22:
541 if (pix.IsExcluded())
542 return kFALSE;
543 val = pix.GetLoGainNumPickup();
544 break;
545 case 23:
546 if (pix.IsExcluded())
547 return kFALSE;
548 val = pix.GetHiGainNumBlackout();
549 break;
550 case 24:
551 if (pix.IsExcluded())
552 return kFALSE;
553 val = pix.GetLoGainNumBlackout();
554 break;
555 case 25:
556 if (pix.IsExcluded())
557 return kFALSE;
558 val = pix.IsHiGainSaturation();
559 break;
560 case 26:
561 if (pix.IsExcluded())
562 return kFALSE;
563 val = pix.GetAbsTimeMean();
564 break;
565 case 27:
566 if (pix.IsExcluded())
567 return kFALSE;
568 val = pix.GetAbsTimeRms();
569 break;
570 case 28:
571 if (pix.IsExcluded())
572 return kFALSE;
573 val = pix.GetPed();
574 break;
575 case 29:
576 if (pix.IsExcluded())
577 return kFALSE;
578 val = pix.GetPedErr();
579 break;
580 case 30:
581 if (pix.IsExcluded())
582 return kFALSE;
583 val = pix.GetPedRms();
584 break;
585 case 31:
586 if (pix.IsExcluded())
587 return kFALSE;
588 val = pix.GetPedErr()/2.;
589 break;
590 case 32:
591 val = pix.GetMean() == 0. ? 0. : pix.GetHiGainRms()/pix.GetMean();
592 break;
593 default:
594 return kFALSE;
595 }
596
597 return val!=-1.;
598}
599
600
601
602Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &ffactor)
603{
604
605 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[ipx];
606
607 Float_t conv = pix.GetMeanConvFADC2Phe();
608
609 if (conv < 0.)
610 return kFALSE;
611
612 mean = conv;
613 err = pix.GetMeanConvFADC2PheErr();
614 ffactor = pix.GetMeanFFactorFADC2Phot();
615
616 return kTRUE;
617}
618
619
620// --------------------------------------------------------------------------
621//
622// Calculates the average conversion factor FADC counts to photons for pixel sizes.
623// The geometry container is used to get the necessary
624// geometry information (area index). The MCalibrationQECam container is necessary for
625// the quantum efficiency information.
626// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
627// in the calculation of the size average.
628//
629// Returns a TArrayF of dimension 2:
630// arr[0]: averaged conversion factors (default: -1.)
631// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
632//
633TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
634 const UInt_t ai, MBadPixelsCam *bad)
635{
636
637 const Int_t np = GetSize();
638
639 Double_t mean = 0.;
640 Double_t mean2 = 0.;
641 Int_t nr = 0;
642
643 for (int i=0; i<np; i++)
644 {
645 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
646 continue;
647
648 const UInt_t aidx = geom[i].GetAidx();
649
650 if (ai != aidx)
651 continue;
652
653 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
654 if (!qepix.IsAverageQECombinedAvailable())
655 continue;
656
657 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
658 const Float_t conv = pix.GetMeanConvFADC2Phe();
659 const Float_t qe = qepix.GetQECascadesCombined();
660
661 mean += conv/qe;
662 mean2 += conv*conv/qe/qe;
663 nr ++;
664
665 }
666
667 TArrayF arr(2);
668 arr[0] = nr ? mean/nr : -1;
669 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
670
671 return arr;
672}
673
674// --------------------------------------------------------------------------
675//
676// Calculates the average conversion factor FADC counts to photons for camera sectors.
677// The geometry container is used to get the necessary
678// geometry information (area index). The MCalibrationQECam container is necessary for
679// the quantum efficiency information.
680// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
681// in the calculation of the size average.
682//
683// Returns a TArrayF of dimension 2:
684// arr[0]: averaged conversion factors (default: -1.)
685// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
686//
687TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam,
688 const UInt_t sec, MBadPixelsCam *bad)
689{
690 const Int_t np = GetSize();
691
692 Double_t mean = 0.;
693 Double_t mean2 = 0.;
694 Int_t nr = 0;
695
696 for (int i=0; i<np; i++)
697 {
698 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
699 continue;
700
701 const UInt_t sector = geom[i].GetSector();
702
703 if (sector != sec)
704 continue;
705
706 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
707 if (!qepix.IsAverageQECombinedAvailable())
708 continue;
709
710 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
711 const Float_t conv = pix.GetMeanConvFADC2Phe();
712 const Float_t qe = qepix.GetQECascadesCombined();
713
714 mean += conv/qe;
715 mean2 += conv*conv/qe/qe;
716 nr ++;
717
718 }
719
720 TArrayF arr(2);
721 arr[0] = nr ? mean/nr : -1;
722 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
723 return arr;
724}
725
726// --------------------------------------------------------------------------
727//
728// Calculates the average mean arrival times for pixel sizes.
729// The geometry container is used to get the necessary
730// geometry information (area index).
731// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
732// in the calculation of the size average.
733//
734// Returns a TArrayF of dimension 2:
735// arr[0]: averaged mean arrival times (default: -1.)
736// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
737//
738TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerArea(const MGeomCam &geom,
739 const UInt_t ai, MBadPixelsCam *bad)
740{
741
742 const Int_t np = GetSize();
743
744 Double_t mean = 0.;
745 Double_t mean2 = 0.;
746 Int_t nr = 0;
747
748 for (int i=0; i<np; i++)
749 {
750 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
751 continue;
752
753 const UInt_t aidx = geom[i].GetAidx();
754
755 if (ai != aidx)
756 continue;
757
758 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
759 const Float_t time = pix.GetAbsTimeMean();
760
761 mean += time ;
762 mean2 += time*time;
763 nr ++;
764
765 }
766
767 TArrayF arr(2);
768 arr[0] = nr ? mean/nr : -1;
769 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
770 return arr;
771}
772
773// --------------------------------------------------------------------------
774//
775// Calculates the average mean arrival times for camera sectors.
776// The geometry container is used to get the necessary
777// geometry information (area index).
778// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
779// in the calculation of the size average.
780//
781// Returns a TArrayF of dimension 2:
782// arr[0]: averaged mean arrival times (default: -1.)
783// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
784//
785TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerSector(const MGeomCam &geom,
786 const UInt_t sec, MBadPixelsCam *bad)
787{
788 const Int_t np = GetSize();
789
790 Double_t mean = 0.;
791 Double_t mean2 = 0.;
792 Int_t nr = 0;
793
794 for (int i=0; i<np; i++)
795 {
796 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
797 continue;
798
799 const UInt_t sector = geom[i].GetSector();
800
801 if (sector != sec)
802 continue;
803
804 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
805 const Float_t time = pix.GetAbsTimeMean();
806
807 mean += time;
808 mean2 += time*time;
809 nr ++;
810
811 }
812
813 TArrayF arr(2);
814 arr[0] = nr ? mean/nr : -1;
815 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
816 return arr;
817}
818
819// --------------------------------------------------------------------------
820//
821// Calculates the average arrival time RMSs for pixel sizes.
822// The geometry container is used to get the necessary
823// geometry information (area index).
824// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
825// in the calculation of the size average.
826//
827// Returns a TArrayF of dimension 2:
828// arr[0]: averaged arrival time RMSs (default: -1.)
829// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
830//
831TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerArea ( const MGeomCam &geom,
832 const UInt_t ai, MBadPixelsCam *bad)
833{
834
835 const Int_t np = GetSize();
836
837 Double_t mean = 0.;
838 Double_t mean2 = 0.;
839 Int_t nr = 0;
840
841 for (int i=0; i<np; i++)
842 {
843 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
844 continue;
845
846 const UInt_t aidx = geom[i].GetAidx();
847
848 if (ai != aidx)
849 continue;
850
851 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
852 const Float_t rms = pix.GetAbsTimeRms();
853
854 mean += rms;
855 mean2 += rms*rms;
856 nr ++;
857
858 }
859
860 TArrayF arr(2);
861 arr[0] = nr ? mean/nr : -1;
862 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
863 return arr;
864}
865
866// --------------------------------------------------------------------------
867//
868// Calculates the average arrival time RMSs for camera sectors.
869// The geometry container is used to get the necessary
870// geometry information (area index).
871// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
872// in the calculation of the size average.
873//
874// Returns a TArrayF of dimension 2:
875// arr[0]: averaged arrival time RMSs (default: -1.)
876// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
877//
878TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerSector( const MGeomCam &geom, const UInt_t sec, MBadPixelsCam *bad)
879{
880 const Int_t np = GetSize();
881
882 Double_t mean = 0.;
883 Double_t mean2 = 0.;
884 Int_t nr = 0;
885
886 for (int i=0; i<np; i++)
887 {
888 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
889 continue;
890
891 const UInt_t sector = geom[i].GetSector();
892
893 if (sector != sec)
894 continue;
895
896 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
897 const Float_t rms = pix.GetAbsTimeRms();
898
899 mean += rms;
900 mean2 += rms*rms;
901 nr ++;
902 }
903
904 TArrayF arr(2);
905 arr[0] = nr ? mean/nr : -1;
906 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
907
908 return arr;
909}
Note: See TracBrowser for help on using the repository browser.