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

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