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

Last change on this file since 5710 was 5661, checked in by gaug, 20 years ago
*** empty log message ***
File size: 28.9 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 {
172 fPixels->AddAt(new MCalibrationChargePix,i);
173 (*this)[i].SetPixId(i);
174 }
175}
176
177// -------------------------------------------------------------------
178//
179// Add MCalibrationChargePix's in the ranges from - to to fAverageAreas
180//
181void MCalibrationChargeCam::AddArea(const UInt_t a, const UInt_t b)
182{
183 for (UInt_t i=a; i<b; i++)
184 {
185 fAverageAreas->AddAt(new MCalibrationChargePix,i);
186 GetAverageArea(i).SetPixId(i);
187 }
188}
189
190// -------------------------------------------------------------------
191//
192// Add MCalibrationChargePix's in the ranges from - to to fAverageSectors
193//
194void MCalibrationChargeCam::AddSector(const UInt_t a, const UInt_t b)
195{
196 for (UInt_t i=a; i<b; i++)
197 {
198 fAverageSectors->AddAt(new MCalibrationChargePix,i);
199 GetAverageSector(i).SetPixId(i);
200 }
201}
202
203
204// --------------------------------------
205//
206// Sets all variable to 0.
207// Sets all flags to kFALSE
208// Calls MCalibrationCam::Clear()
209//
210void MCalibrationChargeCam::Clear(Option_t *o)
211{
212
213 SetFFactorMethodValid ( kFALSE );
214
215 fNumPhotonsBlindPixelMethod = 0.;
216 fNumPhotonsFFactorMethod = 0.;
217 fNumPhotonsPINDiodeMethod = 0.;
218 fNumPhotonsBlindPixelMethodErr = 0.;
219 fNumPhotonsFFactorMethodErr = 0.;
220 fNumPhotonsPINDiodeMethodErr = 0.;
221
222 MCalibrationCam::Clear();
223
224 return;
225}
226
227// -----------------------------------------------
228//
229// Sets the kFFactorMethodValid bit from outside
230//
231void MCalibrationChargeCam::SetFFactorMethodValid(const Bool_t b)
232{
233 b ? SETBIT(fFlags, kFFactorMethodValid) : CLRBIT(fFlags, kFFactorMethodValid);
234}
235
236
237// --------------------------------------------------------------------------
238//
239// Test bit kFFactorMethodValid
240//
241Bool_t MCalibrationChargeCam::IsFFactorMethodValid() const
242{
243 return TESTBIT(fFlags,kFFactorMethodValid);
244}
245
246// --------------------------------------------------------------------------
247//
248// Print first the well fitted pixels
249// and then the ones which are not FitValid
250//
251void MCalibrationChargeCam::Print(Option_t *o) const
252{
253
254 *fLog << all << GetDescriptor() << ":" << endl;
255 int id = 0;
256
257 *fLog << all << "Calibrated pixels:" << endl;
258 *fLog << all << endl;
259
260 TIter Next(fPixels);
261 MCalibrationChargePix *pix;
262 while ((pix=(MCalibrationChargePix*)Next()))
263 {
264
265 if (!pix->IsExcluded())
266 {
267
268 *fLog << all
269 << Form("%s%3i","Pixel: ",pix->GetPixId())
270 << Form("%s%4.2f%s%4.2f"," Ped.Rms: ",pix->GetPedRms(),"+-",pix->GetPedRmsErr())
271 << Form("%s%4.2f%s%4.2f"," Charge: " ,pix->GetConvertedMean(),"+-",pix->GetConvertedSigma())
272 << Form("%s%4.2f%s%4.2f"," Red.Sigma: ",pix->GetConvertedRSigma(),"+-",pix->GetConvertedRSigmaErr())
273 << Form("%s%4.2f%s%4.2f"," Num.Phes: ",pix->GetPheFFactorMethod(),"+-",pix->GetPheFFactorMethodErr())
274 << Form("%s%4.2f%s%4.2f"," Conv.FADC2Phe: ",pix->GetMeanConvFADC2Phe(),"+-",pix->GetMeanConvFADC2PheErr())
275 << " Saturated? :" << pix->IsHiGainSaturation()
276 << endl;
277 id++;
278 }
279 }
280
281 *fLog << all << id << " pixels" << endl;
282 id = 0;
283
284
285 *fLog << all << endl;
286 *fLog << all << "Excluded pixels:" << endl;
287 *fLog << all << endl;
288
289 id = 0;
290
291 TIter Next4(fPixels);
292 while ((pix=(MCalibrationChargePix*)Next4()))
293 {
294 if (pix->IsExcluded())
295 {
296 *fLog << all << pix->GetPixId() << " ";
297 id++;
298
299 if (!(id % 25))
300 *fLog << endl;
301 }
302 }
303
304 *fLog << endl;
305 *fLog << all << id << " Excluded pixels " << endl;
306 *fLog << endl;
307
308 *fLog << all << endl;
309 *fLog << all << "Averaged Areas:" << endl;
310 *fLog << all << endl;
311
312 TIter Next5(fAverageAreas);
313 while ((pix=(MCalibrationChargePix*)Next5()))
314 {
315 *fLog << all << Form("%s%3i","Area Idx: ",pix->GetPixId())
316 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
317 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
318 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
319 << " Reduced Sigma: " << pix->GetRSigma()
320 << " Nr Phe's: " << pix->GetPheFFactorMethod()
321 << endl;
322 }
323
324 *fLog << all << endl;
325 *fLog << all << "Averaged Sectors:" << endl;
326 *fLog << all << endl;
327
328 TIter Next6(fAverageSectors);
329 while ((pix=(MCalibrationChargePix*)Next6()))
330 {
331 *fLog << all << Form("%s%3i","Sector: ",pix->GetPixId())
332 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
333 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
334 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
335 << " Reduced Sigma: " << pix->GetRSigma()
336 << " Nr Phe's: " << pix->GetPheFFactorMethod()
337 << endl;
338 }
339 *fLog << all << endl;
340}
341
342
343// --------------------------------------------------------------------------
344//
345// The types are as follows:
346//
347// Fitted values:
348// ==============
349//
350// 0: Fitted Charge (see MCalibrationPix::GetMean())
351// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
352// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
353// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
354//
355// Useful variables derived from the fit results:
356// =============================================
357//
358// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
359// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
360// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
361// 7: Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerCharge())
362// 8: Error of Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerChargeErr())
363//
364// Results of the F-Factor calibration Method:
365// ===========================================
366//
367// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
368// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
369// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2Phe()
370// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2PheErr()
371// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2Phot()
372// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2PhotErr()
373// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
374//
375// Results of the Low-Gain vs. High-Gain Conversion:
376// =================================================
377//
378// 16: Mean Signal Hi-Gain / Mean Signal Lo-Gain (see MCalibrationPix::GetHiLoMeansDivided())
379// 17: Error Signal High-Gain / Signal Low Gain (see MCalibrationPix::GetHiLoMeansDividedErr())
380// 18: Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDivided())
381// 19: Error Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDividedErr())
382//
383// Localized defects:
384// ==================
385//
386// 20: Excluded Pixels
387// 21: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
388// 22: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
389// 23: Number of blackout events in the Hi Gain (see MCalibrationPix::GetHiGainNumBlackout())
390// 24: Number of blackout events in the Lo Gain (see MCalibrationPix::GetLoGainNumBlackout())
391//
392// Other classifications of pixels:
393// ================================
394//
395// 25: Pixels with saturated High-Gain (see MCalibrationPix::IsHiGainSaturation())
396//
397// Calculated absolute arrival times (very low precision!):
398// ========================================================
399//
400// 26: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
401// 27: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
402//
403// Used Pedestals:
404// ===============
405//
406// 28: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
407// 29: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
408// 30: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
409// 31: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
410//
411// Special variables (for data check):
412// ====================================
413//
414// 32: HiGain RMS divided by Mean for every pixel (with inclusion of the excluded pixels)
415//
416Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
417{
418
419 if (idx > GetSize())
420 return kFALSE;
421
422 Float_t area = cam[idx].GetA();
423
424 if (area == 0)
425 return kFALSE;
426
427 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
428
429 switch (type)
430 {
431 case 0:
432 if (pix.IsExcluded())
433 return kFALSE;
434 val = pix.GetMean();
435 break;
436 case 1:
437 if (pix.IsExcluded())
438 return kFALSE;
439 val = pix.GetMeanErr();
440 break;
441 case 2:
442 if (pix.IsExcluded())
443 return kFALSE;
444 val = pix.GetSigma();
445 break;
446 case 3:
447 if (pix.IsExcluded())
448 return kFALSE;
449 val = pix.GetSigmaErr();
450 break;
451 case 4:
452 if (pix.IsExcluded())
453 return kFALSE;
454 val = pix.GetProb();
455 break;
456 case 5:
457 if (pix.IsExcluded())
458 return kFALSE;
459 if (pix.GetRSigma() == -1.)
460 return kFALSE;
461 val = pix.GetRSigma();
462 break;
463 case 6:
464 if (pix.IsExcluded())
465 return kFALSE;
466 if (pix.GetRSigma() == -1.)
467 return kFALSE;
468 val = pix.GetRSigmaErr();
469 break;
470 case 7:
471 if (pix.IsExcluded())
472 return kFALSE;
473 val = pix.GetRSigmaPerCharge();
474 break;
475 case 8:
476 if (pix.IsExcluded())
477 return kFALSE;
478 val = pix.GetRSigmaPerChargeErr();
479 break;
480 case 9:
481 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
482 return kFALSE;
483 val = pix.GetPheFFactorMethod();
484 break;
485 case 10:
486 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
487 return kFALSE;
488 val = pix.GetPheFFactorMethodErr();
489 break;
490 case 11:
491 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
492 return kFALSE;
493 val = pix.GetMeanConvFADC2Phe();
494 break;
495 case 12:
496 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
497 return kFALSE;
498 val = pix.GetMeanConvFADC2PheErr();
499 break;
500 case 13:
501 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
502 return kFALSE;
503 val = pix.GetMeanFFactorFADC2Phot();
504 break;
505 case 14:
506 if (pix.IsExcluded() || !pix.IsFFactorMethodValid())
507 return kFALSE;
508 val = pix.GetMeanFFactorFADC2PhotErr();
509 break;
510 case 15:
511 if (pix.IsExcluded())
512 return kFALSE;
513 if (pix.IsFFactorMethodValid())
514 val = 1;
515 else
516 return kFALSE;
517 break;
518 case 16:
519 if (pix.IsExcluded())
520 return kFALSE;
521 val = pix.GetHiLoMeansDivided();
522 break;
523 case 17:
524 if (pix.IsExcluded())
525 return kFALSE;
526 val = pix.GetHiLoMeansDividedErr();
527 break;
528 case 18:
529 if (pix.IsExcluded())
530 return kFALSE;
531 val = pix.GetHiLoSigmasDivided();
532 break;
533 case 19:
534 if (pix.IsExcluded())
535 return kFALSE;
536 val = pix.GetHiLoSigmasDividedErr();
537 break;
538 case 20:
539 if (pix.IsExcluded())
540 val = 1.;
541 else
542 return kFALSE;
543 break;
544 case 21:
545 if (pix.IsExcluded())
546 return kFALSE;
547 val = pix.GetHiGainNumPickup();
548 break;
549 case 22:
550 if (pix.IsExcluded())
551 return kFALSE;
552 val = pix.GetLoGainNumPickup();
553 break;
554 case 23:
555 if (pix.IsExcluded())
556 return kFALSE;
557 val = pix.GetHiGainNumBlackout();
558 break;
559 case 24:
560 if (pix.IsExcluded())
561 return kFALSE;
562 val = pix.GetLoGainNumBlackout();
563 break;
564 case 25:
565 if (pix.IsExcluded())
566 return kFALSE;
567 val = pix.IsHiGainSaturation();
568 break;
569 case 26:
570 if (pix.IsExcluded())
571 return kFALSE;
572 val = pix.GetAbsTimeMean();
573 break;
574 case 27:
575 if (pix.IsExcluded())
576 return kFALSE;
577 val = pix.GetAbsTimeRms();
578 break;
579 case 28:
580 if (pix.IsExcluded())
581 return kFALSE;
582 val = pix.GetPed();
583 break;
584 case 29:
585 if (pix.IsExcluded())
586 return kFALSE;
587 val = pix.GetPedErr();
588 break;
589 case 30:
590 if (pix.IsExcluded())
591 return kFALSE;
592 val = pix.GetPedRms();
593 break;
594 case 31:
595 if (pix.IsExcluded())
596 return kFALSE;
597 val = pix.GetPedErr()/2.;
598 break;
599 case 32:
600 val = pix.GetMean() == 0. ? 0. : pix.GetRms()/pix.GetMean();
601 break;
602 default:
603 return kFALSE;
604 }
605
606 return val!=-1.;
607}
608
609
610
611Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &ffactor)
612{
613
614 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[ipx];
615
616 Float_t conv = pix.GetMeanConvFADC2Phe();
617
618 if (conv < 0.)
619 return kFALSE;
620
621 mean = conv;
622 err = pix.GetMeanConvFADC2PheErr();
623 ffactor = pix.GetMeanFFactorFADC2Phot();
624
625 return kTRUE;
626}
627
628
629// --------------------------------------------------------------------------
630//
631// Calculates the average conversion factor FADC counts to photons for pixel sizes.
632// The geometry container is used to get the necessary
633// geometry information (area index). The MCalibrationQECam container is necessary for
634// the quantum efficiency information.
635// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
636// in the calculation of the size average.
637//
638// Returns a TArrayF of dimension 2:
639// arr[0]: averaged conversion factors (default: -1.)
640// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
641//
642TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
643 const UInt_t ai, MBadPixelsCam *bad)
644{
645
646 const Int_t np = GetSize();
647
648 Double_t mean = 0.;
649 Double_t mean2 = 0.;
650 Int_t nr = 0;
651
652 for (int i=0; i<np; i++)
653 {
654 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
655 continue;
656
657 const UInt_t aidx = geom[i].GetAidx();
658
659 if (ai != aidx)
660 continue;
661
662 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
663 if (!qepix.IsAverageQECombinedAvailable())
664 continue;
665
666 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
667 const Float_t conv = pix.GetMeanConvFADC2Phe();
668 const Float_t qe = qepix.GetQECascadesCombined();
669
670 mean += conv/qe;
671 mean2 += conv*conv/qe/qe;
672 nr ++;
673
674 }
675
676 TArrayF arr(2);
677 arr[0] = nr ? mean/nr : -1.;
678 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0.;
679
680 return arr;
681}
682
683// --------------------------------------------------------------------------
684//
685// Calculates the average conversion factor FADC counts to photons for camera sectors.
686// The geometry container is used to get the necessary
687// geometry information (area index). The MCalibrationQECam container is necessary for
688// the quantum efficiency information.
689// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
690// in the calculation of the size average.
691//
692// Returns a TArrayF of dimension 2:
693// arr[0]: averaged conversion factors (default: -1.)
694// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
695//
696TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam,
697 const UInt_t sec, MBadPixelsCam *bad)
698{
699 const Int_t np = GetSize();
700
701 Double_t mean = 0.;
702 Double_t mean2 = 0.;
703 Int_t nr = 0;
704
705 for (int i=0; i<np; i++)
706 {
707 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
708 continue;
709
710 const UInt_t sector = geom[i].GetSector();
711
712 if (sector != sec)
713 continue;
714
715 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
716 if (!qepix.IsAverageQECombinedAvailable())
717 continue;
718
719 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
720 const Float_t conv = pix.GetMeanConvFADC2Phe();
721 const Float_t qe = qepix.GetQECascadesCombined();
722
723 mean += conv/qe;
724 mean2 += conv*conv/qe/qe;
725 nr ++;
726
727 }
728
729 TArrayF arr(2);
730 arr[0] = nr ? mean/nr : -1;
731 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
732 return arr;
733}
734
735// --------------------------------------------------------------------------
736//
737// Calculates the average mean arrival times for pixel sizes.
738// The geometry container is used to get the necessary
739// geometry information (area index).
740// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
741// in the calculation of the size average.
742//
743// Returns a TArrayF of dimension 2:
744// arr[0]: averaged mean arrival times (default: -1.)
745// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
746//
747TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerArea(const MGeomCam &geom,
748 const UInt_t ai, MBadPixelsCam *bad)
749{
750
751 const Int_t np = GetSize();
752
753 Double_t mean = 0.;
754 Double_t mean2 = 0.;
755 Int_t nr = 0;
756
757 for (int i=0; i<np; i++)
758 {
759 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
760 continue;
761
762 const UInt_t aidx = geom[i].GetAidx();
763
764 if (ai != aidx)
765 continue;
766
767 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
768 const Float_t time = pix.GetAbsTimeMean();
769
770 mean += time ;
771 mean2 += time*time;
772 nr ++;
773
774 }
775
776 TArrayF arr(2);
777 arr[0] = nr ? mean/nr : -1;
778 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
779 return arr;
780}
781
782// --------------------------------------------------------------------------
783//
784// Calculates the average mean arrival times for camera sectors.
785// The geometry container is used to get the necessary
786// geometry information (area index).
787// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
788// in the calculation of the size average.
789//
790// Returns a TArrayF of dimension 2:
791// arr[0]: averaged mean arrival times (default: -1.)
792// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
793//
794TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerSector(const MGeomCam &geom,
795 const UInt_t sec, MBadPixelsCam *bad)
796{
797 const Int_t np = GetSize();
798
799 Double_t mean = 0.;
800 Double_t mean2 = 0.;
801 Int_t nr = 0;
802
803 for (int i=0; i<np; i++)
804 {
805 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
806 continue;
807
808 const UInt_t sector = geom[i].GetSector();
809
810 if (sector != sec)
811 continue;
812
813 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
814 const Float_t time = pix.GetAbsTimeMean();
815
816 mean += time;
817 mean2 += time*time;
818 nr ++;
819
820 }
821
822 TArrayF arr(2);
823 arr[0] = nr ? mean/nr : -1;
824 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
825 return arr;
826}
827
828// --------------------------------------------------------------------------
829//
830// Calculates the average arrival time RMSs for pixel sizes.
831// The geometry container is used to get the necessary
832// geometry information (area index).
833// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
834// in the calculation of the size average.
835//
836// Returns a TArrayF of dimension 2:
837// arr[0]: averaged arrival time RMSs (default: -1.)
838// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
839//
840TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerArea ( const MGeomCam &geom,
841 const UInt_t ai, MBadPixelsCam *bad)
842{
843
844 const Int_t np = GetSize();
845
846 Double_t mean = 0.;
847 Double_t mean2 = 0.;
848 Int_t nr = 0;
849
850 for (int i=0; i<np; i++)
851 {
852 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
853 continue;
854
855 const UInt_t aidx = geom[i].GetAidx();
856
857 if (ai != aidx)
858 continue;
859
860 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
861 const Float_t rms = pix.GetAbsTimeRms();
862
863 mean += rms;
864 mean2 += rms*rms;
865 nr ++;
866
867 }
868
869 TArrayF arr(2);
870 arr[0] = nr ? mean/nr : -1;
871 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
872 return arr;
873}
874
875// --------------------------------------------------------------------------
876//
877// Calculates the average arrival time RMSs for camera sectors.
878// The geometry container is used to get the necessary
879// geometry information (area index).
880// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
881// in the calculation of the size average.
882//
883// Returns a TArrayF of dimension 2:
884// arr[0]: averaged arrival time RMSs (default: -1.)
885// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
886//
887TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerSector( const MGeomCam &geom, const UInt_t sec, MBadPixelsCam *bad)
888{
889 const Int_t np = GetSize();
890
891 Double_t mean = 0.;
892 Double_t mean2 = 0.;
893 Int_t nr = 0;
894
895 for (int i=0; i<np; i++)
896 {
897 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
898 continue;
899
900 const UInt_t sector = geom[i].GetSector();
901
902 if (sector != sec)
903 continue;
904
905 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
906 const Float_t rms = pix.GetAbsTimeRms();
907
908 mean += rms;
909 mean2 += rms*rms;
910 nr ++;
911 }
912
913 TArrayF arr(2);
914 arr[0] = nr ? mean/nr : -1;
915 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
916
917 return arr;
918}
Note: See TracBrowser for help on using the repository browser.