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

Last change on this file since 7977 was 7059, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 36.4 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::MergeHiLoConversionFactors(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->IsFFactorMethodValid())
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 << "Previously 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 << all << endl;
332 *fLog << all << "New Excluded pixels:" << endl;
333 *fLog << all << endl;
334
335 TIter Next5(fPixels);
336 while ((pix=(MCalibrationChargePix*)Next5()))
337 {
338 if (!pix->IsFFactorMethodValid() && !pix->IsExcluded())
339 {
340 *fLog << all << pix->GetPixId() << " ";
341 id++;
342
343 if (!(id % 25))
344 *fLog << endl;
345 }
346 }
347
348 *fLog << endl;
349 *fLog << all << id << " Excluded pixels " << endl;
350 *fLog << endl;
351
352 *fLog << all << endl;
353 *fLog << all << "Averaged Areas:" << endl;
354 *fLog << all << endl;
355
356 TIter Next6(fAverageAreas);
357 while ((pix=(MCalibrationChargePix*)Next6()))
358 {
359 *fLog << all << Form("%s%3i","Area Idx: ",pix->GetPixId())
360 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
361 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
362 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
363 << " Reduced Sigma: " << pix->GetRSigma()
364 << " Nr Phe's: " << pix->GetPheFFactorMethod()
365 << endl;
366 }
367
368 *fLog << all << endl;
369 *fLog << all << "Averaged Sectors:" << endl;
370 *fLog << all << endl;
371
372 TIter Next7(fAverageSectors);
373 while ((pix=(MCalibrationChargePix*)Next7()))
374 {
375 *fLog << all << Form("%s%3i","Sector: ",pix->GetPixId())
376 << " Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
377 << " Mean signal: " << pix->GetMean() << " +- " << pix->GetMeanErr()
378 << " Sigma signal: " << pix->GetSigma() << " +- "<< pix->GetSigmaErr()
379 << " Reduced Sigma: " << pix->GetRSigma()
380 << " Nr Phe's: " << pix->GetPheFFactorMethod()
381 << endl;
382 }
383 *fLog << all << endl;
384}
385
386
387// --------------------------------------------------------------------------
388//
389// The types are as follows:
390//
391// Fitted values:
392// ==============
393//
394// 0: Fitted Charge (see MCalibrationPix::GetMean())
395// 1: Error of fitted Charge (see MCalibrationPix::GetMeanErr())
396// 2: Sigma of fitted Charge (see MCalibrationPix::GetSigma())
397// 3: Error of Sigma of fitted Charge (see MCalibrationPix::GetSigmaErr())
398//
399// Useful variables derived from the fit results:
400// =============================================
401//
402// 4: Probability Gauss fit Charge distribution (see MCalibrationPix::GetProb())
403// 5: Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigma())
404// 6: Error Reduced Sigma of fitted Charge (see MCalibrationChargePix::GetRSigmaErr())
405// 7: Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerCharge())
406// 8: Error of Reduced Sigma per Charge (see MCalibrationChargePix::GetRSigmaPerChargeErr())
407//
408// Results of the F-Factor calibration Method:
409// ===========================================
410//
411// 9: Nr. Photo-electrons from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethod())
412// 10: Error Nr. Photo-el. from F-Factor Method (see MCalibrationChargePix::GetPheFFactorMethodErr())
413// 11: Conversion factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2Phe()
414// 12: Error conv. factor from F-Factor Method (see MCalibrationChargePix::GetMeanConvFADC2PheErr()
415// 13: Overall F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2Phot()
416// 14: Error F-Factor from F-Factor Method (see MCalibrationChargePix::GetMeanFFactorFADC2PhotErr()
417// 15: Pixels valid calibration F-Factor-Method (see MCalibrationChargePix::IsFFactorMethodValid())
418//
419// Results of the Low-Gain vs. High-Gain Conversion:
420// =================================================
421//
422// 16: Mean Signal Hi-Gain / Mean Signal Lo-Gain (see MCalibrationPix::GetHiLoMeansDivided())
423// 17: Error Signal High-Gain / Signal Low Gain (see MCalibrationPix::GetHiLoMeansDividedErr())
424// 18: Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDivided())
425// 19: Error Sigma High-Gain / Sigma Low Gain (see MCalibrationPix::GetHiLoSigmasDividedErr())
426//
427// Localized defects:
428// ==================
429//
430// 20: Excluded Pixels
431// 21: Number of pickup events in the Hi Gain (see MCalibrationPix::GetHiGainNumPickup())
432// 22: Number of pickup events in the Lo Gain (see MCalibrationPix::GetLoGainNumPickup())
433// 23: Number of blackout events in the Hi Gain (see MCalibrationPix::GetHiGainNumBlackout())
434// 24: Number of blackout events in the Lo Gain (see MCalibrationPix::GetLoGainNumBlackout())
435//
436// Other classifications of pixels:
437// ================================
438//
439// 25: Pixels with saturated High-Gain (see MCalibrationPix::IsHiGainSaturation())
440//
441// Calculated absolute arrival times (very low precision!):
442// ========================================================
443//
444// 26: Absolute Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeMean())
445// 27: RMS Ab. Arrival time of the signal (see MCalibrationChargePix::GetAbsTimeRms())
446//
447// Used Pedestals:
448// ===============
449//
450// 28: Pedestal for entire signal extr. range (see MCalibrationChargePix::Ped())
451// 29: Error Pedestal entire signal extr. range (see MCalibrationChargePix::PedErr())
452// 30: Ped. RMS entire signal extraction range (see MCalibrationChargePix::PedRms())
453// 31: Error Ped. RMS entire signal extr. range (see MCalibrationChargePix::PedRmsErr())
454//
455// Special variables (for data check):
456// ====================================
457//
458// 32: HiGain RMS divided by Mean for every pixel (with inclusion of the excluded pixels)
459//
460Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
461{
462
463 if (idx > GetSize())
464 return kFALSE;
465
466 Float_t area = cam[idx].GetA();
467
468 if (area == 0)
469 return kFALSE;
470
471 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[idx];
472
473 switch (type)
474 {
475 case 0:
476 if (pix.IsExcluded())
477 return kFALSE;
478 val = pix.GetConvertedMean();
479 break;
480 case 1:
481 if (pix.IsExcluded())
482 return kFALSE;
483 val = pix.GetConvertedMeanErr();
484 break;
485 case 2:
486 if (pix.IsExcluded())
487 return kFALSE;
488 val = pix.GetConvertedSigma();
489 break;
490 case 3:
491 if (pix.IsExcluded())
492 return kFALSE;
493 val = pix.GetConvertedSigmaErr();
494 break;
495 case 4:
496 if (pix.IsExcluded())
497 return kFALSE;
498 val = pix.GetProb();
499 break;
500 case 5:
501 if (!pix.IsFFactorMethodValid())
502 return kFALSE;
503 if (pix.GetRSigma() == -1.)
504 return kFALSE;
505 val = pix.GetConvertedRSigma();
506 break;
507 case 6:
508 if (!pix.IsFFactorMethodValid())
509 return kFALSE;
510 if (pix.GetRSigma() == -1.)
511 return kFALSE;
512 val = pix.GetConvertedRSigmaErr();
513 break;
514 case 7:
515 if (!pix.IsFFactorMethodValid())
516 return kFALSE;
517 val = pix.GetRSigmaPerCharge();
518 break;
519 case 8:
520 if (!pix.IsFFactorMethodValid())
521 return kFALSE;
522 val = pix.GetRSigmaPerChargeErr();
523 break;
524 case 9:
525 // if (!pix.IsFFactorMethodValid())
526 // return kFALSE;
527 if (pix.IsExcluded())
528 return kFALSE;
529 val = pix.GetPheFFactorMethod();
530 break;
531 case 10:
532 // if (!pix.IsFFactorMethodValid())
533 // return kFALSE;
534 if (pix.IsExcluded())
535 return kFALSE;
536 val = pix.GetPheFFactorMethodErr();
537 if (val < 0.)
538 val = 0.00001;
539 break;
540 case 11:
541 // if (!pix.IsFFactorMethodValid())
542 // return kFALSE;
543 if (pix.IsExcluded())
544 return kFALSE;
545 val = pix.GetMeanConvFADC2Phe();
546 break;
547 case 12:
548 // if (!pix.IsFFactorMethodValid())
549 // return kFALSE;
550 if (pix.IsExcluded())
551 return kFALSE;
552 val = pix.GetMeanConvFADC2PheErr();
553 if (val < 0.)
554 val = 0.00001;
555 break;
556 case 13:
557 // if (!pix.IsFFactorMethodValid())
558 // return kFALSE;
559 if (pix.IsExcluded())
560 return kFALSE;
561 val = pix.GetMeanFFactorFADC2Phot();
562 break;
563 case 14:
564 // if (!pix.IsFFactorMethodValid())
565 // return kFALSE;
566 if (pix.IsExcluded())
567 return kFALSE;
568 val = pix.GetMeanFFactorFADC2PhotErr();
569 if (val < 0.)
570 val = 0.00001;
571 break;
572 case 15:
573 if (pix.IsExcluded())
574 return kFALSE;
575 if (pix.IsFFactorMethodValid())
576 val = 1;
577 else
578 return kFALSE;
579 break;
580 case 16:
581 if (pix.IsExcluded())
582 return kFALSE;
583 val = pix.GetHiLoMeansDivided();
584 break;
585 case 17:
586 if (pix.IsExcluded())
587 return kFALSE;
588 val = pix.GetHiLoMeansDividedErr();
589 break;
590 case 18:
591 if (pix.IsExcluded())
592 return kFALSE;
593 val = pix.GetHiLoSigmasDivided();
594 break;
595 case 19:
596 if (pix.IsExcluded())
597 return kFALSE;
598 val = pix.GetHiLoSigmasDividedErr();
599 break;
600 case 20:
601 if (pix.IsExcluded())
602 val = 1.;
603 else
604 return kFALSE;
605 break;
606 case 21:
607 if (pix.IsExcluded())
608 return kFALSE;
609 val = pix.GetHiGainNumPickup();
610 break;
611 case 22:
612 if (pix.IsExcluded())
613 return kFALSE;
614 val = pix.GetLoGainNumPickup();
615 break;
616 case 23:
617 if (pix.IsExcluded())
618 return kFALSE;
619 val = pix.GetHiGainNumBlackout();
620 break;
621 case 24:
622 if (pix.IsExcluded())
623 return kFALSE;
624 val = pix.GetLoGainNumBlackout();
625 break;
626 case 25:
627 if (pix.IsExcluded())
628 return kFALSE;
629 val = pix.IsHiGainSaturation();
630 break;
631 case 26:
632 if (pix.IsExcluded())
633 return kFALSE;
634 val = pix.GetAbsTimeMean();
635 break;
636 case 27:
637 if (pix.IsExcluded())
638 return kFALSE;
639 val = pix.GetAbsTimeRms();
640 break;
641 case 28:
642 if (pix.IsExcluded())
643 return kFALSE;
644 val = pix.GetPed();
645 break;
646 case 29:
647 if (pix.IsExcluded())
648 return kFALSE;
649 val = pix.GetPedErr();
650 break;
651 case 30:
652 if (pix.IsExcluded())
653 return kFALSE;
654 val = pix.GetPedRms();
655 break;
656 case 31:
657 if (pix.IsExcluded())
658 return kFALSE;
659 val = pix.GetPedErr()/2.;
660 break;
661 case 32:
662 if (pix.IsExcluded())
663 return kFALSE;
664 val = pix.GetMean() == 0. ? 0. : pix.GetRms()/pix.GetMean();
665 break;
666 case 33:
667 if (pix.IsExcluded())
668 return kFALSE;
669 if (pix.GetMean() == 0.)
670 val = 0.;
671 else
672 val = pix.GetSigmaErr()/pix.GetMean() + pix.GetRms()/pix.GetMean()/pix.GetMean()*pix.GetMeanErr();
673 break;
674 default:
675 return kFALSE;
676 }
677
678 return val!=-1.;
679}
680
681
682Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &ffactor)
683{
684
685 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[ipx];
686
687 Float_t conv = pix.GetMeanConvFADC2Phe();
688
689 if (conv < 0.)
690 return kFALSE;
691
692 mean = conv;
693 err = pix.GetMeanConvFADC2PheErr();
694 ffactor = pix.GetMeanFFactorFADC2Phot();
695
696 return kTRUE;
697}
698
699
700// --------------------------------------------------------------------------
701//
702// Calculates the average conversion factor FADC counts to photons for pixel sizes.
703// The geometry container is used to get the necessary
704// geometry information (area index). The MCalibrationQECam container is necessary for
705// the quantum efficiency information.
706// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
707// in the calculation of the size average.
708//
709// Returns a TArrayF of dimension 2:
710// arr[0]: averaged conversion factors (default: -1.)
711// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
712//
713TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
714 const UInt_t ai, MBadPixelsCam *bad)
715{
716
717 const Int_t np = GetSize();
718
719 Double_t mean = 0.;
720 Double_t mean2 = 0.;
721 Int_t nr = 0;
722
723 MHCamera convcam(geom,"ConvFactors","Conversion Factors;Conv Factor [phot/FADC cnts];channels");
724
725 for (int i=0; i<np; i++)
726 {
727 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
728 continue;
729
730 if (bad && (*bad)[i].IsUncalibrated(MBadPixelsPix::kDeviatingNumPhes))
731 continue;
732
733 const UInt_t aidx = geom[i].GetAidx();
734
735 if (ai != aidx)
736 continue;
737
738 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
739 if (!qepix.IsAverageQEFFactorAvailable())
740 continue;
741
742 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
743 const Float_t conv = pix.GetMeanConvFADC2Phe()/qepix.GetQECascadesFFactor();
744
745 mean += conv;
746 mean2 += conv*conv;
747 nr ++;
748
749 convcam.Fill(i,conv);
750 convcam.SetUsed(i);
751 }
752
753 Float_t mn = nr ? mean/nr : -1.;
754 Float_t sg = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0.;
755
756 const Int_t aidx = (Int_t)ai;
757
758 TH1D *h = convcam.ProjectionS(TArrayI(),TArrayI(1,&aidx),"_py",750);
759 h->SetDirectory(NULL);
760
761 TF1 *fit = NULL;
762
763 if (geom.InheritsFrom("MGeomCamMagic"))
764 {
765
766 fit = new TF1("fit","gaus",0.4,5.);
767
768 // Fix the ranges, as found by Nadia
769 if(aidx == 0)
770 h->Fit(fit, "REQ0", "",0.4,1.5);
771 else
772 h->Fit(fit ,"REQ0", "",1.,5.);
773 }
774 else
775 {
776 h->Fit("gaus","Q0");
777 fit = h->GetFunction("gaus");
778 }
779
780 Float_t ci2 = fit->GetChisquare();
781 Float_t sigma = fit->GetParameter(2);
782
783 if (ci2 > 500. || sigma > sg)
784 {
785 if (geom.InheritsFrom("MGeomCamMagic"))
786 {
787 // Fix the ranges, as found by Nadia
788 if(aidx == 0)
789 h->Fit(fit, "REQ0", "",0.4,1.5);
790 else
791 h->Fit(fit, "REQ0", "",1.,5.);
792 }
793 else
794 {
795 h->Fit("gaus","MREQ0");
796 fit = h->GetFunction("gaus");
797 }
798
799 ci2 = fit->GetChisquare();
800 sigma = fit->GetParameter(2);
801 }
802
803 const Int_t ndf = fit->GetNDF();
804
805 if (ci2 < 500. && sigma < sg && ndf > 2)
806 {
807 mn = fit->GetParameter(1);
808 sg = sigma;
809 }
810
811 *fLog << inf << "Conversion Factors to photons area idx: " << aidx << ":" << endl;
812 *fLog << inf << "Mean: " << Form("%4.3f",mn)
813 << "+-" << Form("%4.3f",fit->GetParError(1))
814 << " Sigma: " << Form("%4.3f",sg) << "+-" << Form("%4.3f",fit->GetParError(2))
815 << " Chisquare: " << Form("%4.3f",fit->GetChisquare()) << " NDF : " << ndf << endl;
816
817 delete fit;
818 delete h;
819 gROOT->GetListOfFunctions()->Remove(fit);
820
821 TArrayF arr(2);
822 arr[0] = mn;
823 arr[1] = sg;
824
825 return arr;
826}
827
828// --------------------------------------------------------------------------
829//
830// Calculates the average conversion factor FADC counts to equiv. photo-electrons for pixel sizes.
831// The geometry container is used to get the necessary
832// geometry information (area index). The MCalibrationQECam container is necessary for
833// the quantum efficiency information.
834// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
835// in the calculation of the size average.
836//
837// Returns a TArrayF of dimension 2:
838// arr[0]: averaged conversion factors (default: -1.)
839// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
840//
841TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhePerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam,
842 const UInt_t ai, MBadPixelsCam *bad)
843{
844
845 const Int_t np = GetSize();
846
847 Double_t mean = 0.;
848 Double_t mean2 = 0.;
849 Int_t nr = 0;
850
851 MHCamera convcam(geom,"ConvFactors","Conversion Factors;Conv Factor [phe/FADC cnts];channels");
852
853 for (int i=0; i<np; i++)
854 {
855 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
856 continue;
857
858 if (bad && (*bad)[i].IsUncalibrated(MBadPixelsPix::kDeviatingNumPhes))
859 continue;
860
861 const UInt_t aidx = geom[i].GetAidx();
862
863 if (ai != aidx)
864 continue;
865
866 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
867 if (!qepix.IsAverageQEFFactorAvailable())
868 continue;
869
870 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
871 const Float_t conv = pix.GetMeanConvFADC2Phe()/qepix.GetQECascadesFFactor()*MCalibrationQEPix::gkDefaultAverageQE;
872
873 mean += conv;
874 mean2 += conv*conv;
875 nr ++;
876
877 convcam.Fill(i,conv);
878 convcam.SetUsed(i);
879 }
880
881 Float_t mn = nr ? mean/nr : -1.;
882 Float_t sg = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0.;
883
884 const Int_t aidx = (Int_t)ai;
885
886 TH1D *h = convcam.ProjectionS(TArrayI(),TArrayI(1,&aidx),"_py",750);
887 h->SetDirectory(NULL);
888
889 TF1 *fit = NULL;
890
891 if (geom.InheritsFrom("MGeomCamMagic"))
892 {
893
894 fit = new TF1("fit","gaus",0.01,1.);
895
896 // Fix the ranges, as found by Nadia
897 if(aidx == 0)
898 {h->Fit("fit","REQ0", "",0.07,0.3);}
899 else
900 {h->Fit("fit","REQ0", "",0.15,1.0);}
901 }
902 else
903 {
904 h->Fit("gaus","Q0");
905 fit = h->GetFunction("gaus");
906 }
907
908 Float_t ci2 = fit->GetChisquare();
909 Float_t sigma = fit->GetParameter(2);
910
911 if (ci2 > 500. || sigma > sg)
912 {
913 if (geom.InheritsFrom("MGeomCamMagic"))
914 {
915 // Fix the ranges, as found by Nadia
916 if(aidx == 0)
917 {h->Fit("fit","REQ0", "",0.07,0.3);}
918 else
919 {h->Fit("fit","REQ0", "",0.15,1.0);}
920 }
921 else
922 {
923 h->Fit("gaus","MREQ0");
924 fit = h->GetFunction("gaus");
925 }
926
927 ci2 = fit->GetChisquare();
928 sigma = fit->GetParameter(2);
929 }
930
931 const Int_t ndf = fit->GetNDF();
932
933 if (ci2 < 500. && sigma < sg && ndf > 2)
934 {
935 mn = fit->GetParameter(1);
936 sg = sigma;
937 }
938
939 *fLog << inf << "Conversion Factors to equiv. photo-electrons area idx: " << aidx << ":" << endl;
940 *fLog << inf << "Mean: " << Form("%4.3f",mn)
941 << "+-" << Form("%4.3f",fit->GetParError(1))
942 << " Sigma: " << Form("%4.3f",sg) << "+-" << Form("%4.3f",fit->GetParError(2))
943 << " Chisquare: " << Form("%4.3f",fit->GetChisquare()) << " NDF : " << ndf << endl;
944
945 delete fit;
946 delete h;
947 gROOT->GetListOfFunctions()->Remove(fit);
948
949 TArrayF arr(2);
950 arr[0] = mn;
951 arr[1] = sg;
952
953 return arr;
954}
955
956// --------------------------------------------------------------------------
957//
958// Calculates the average conversion factor FADC counts to photons for camera sectors.
959// The geometry container is used to get the necessary
960// geometry information (area index). The MCalibrationQECam container is necessary for
961// the quantum efficiency information.
962// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
963// in the calculation of the size average.
964//
965// Returns a TArrayF of dimension 2:
966// arr[0]: averaged conversion factors (default: -1.)
967// arr[1]: Error (rms) of averaged conversion factors (default: 0.)
968//
969TArrayF MCalibrationChargeCam::GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam,
970 const UInt_t sec, MBadPixelsCam *bad)
971{
972 const Int_t np = GetSize();
973
974 Double_t mean = 0.;
975 Double_t mean2 = 0.;
976 Int_t nr = 0;
977
978 for (int i=0; i<np; i++)
979 {
980 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
981 continue;
982
983 const UInt_t sector = geom[i].GetSector();
984
985 if (sector != sec)
986 continue;
987
988 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i];
989 if (!qepix.IsAverageQEFFactorAvailable())
990 continue;
991
992 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
993 const Float_t conv = pix.GetMeanConvFADC2Phe();
994 const Float_t qe = qepix.GetQECascadesFFactor();
995
996 mean += conv/qe;
997 mean2 += conv*conv/qe/qe;
998 nr ++;
999
1000 }
1001
1002 TArrayF arr(2);
1003 arr[0] = nr ? mean/nr : -1;
1004 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1005 return arr;
1006}
1007
1008// --------------------------------------------------------------------------
1009//
1010// Calculates the average mean arrival times for pixel sizes.
1011// The geometry container is used to get the necessary
1012// geometry information (area index).
1013// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1014// in the calculation of the size average.
1015//
1016// Returns a TArrayF of dimension 2:
1017// arr[0]: averaged mean arrival times (default: -1.)
1018// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
1019//
1020TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerArea(const MGeomCam &geom,
1021 const UInt_t ai, MBadPixelsCam *bad)
1022{
1023
1024 const Int_t np = GetSize();
1025
1026 Double_t mean = 0.;
1027 Double_t mean2 = 0.;
1028 Int_t nr = 0;
1029
1030 for (int i=0; i<np; i++)
1031 {
1032 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1033 continue;
1034
1035 const UInt_t aidx = geom[i].GetAidx();
1036
1037 if (ai != aidx)
1038 continue;
1039
1040 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1041 const Float_t time = pix.GetAbsTimeMean();
1042
1043 mean += time ;
1044 mean2 += time*time;
1045 nr ++;
1046
1047 }
1048
1049 TArrayF arr(2);
1050 arr[0] = nr ? mean/nr : -1;
1051 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1052 return arr;
1053}
1054
1055// --------------------------------------------------------------------------
1056//
1057// Calculates the average mean arrival times for camera sectors.
1058// The geometry container is used to get the necessary
1059// geometry information (area index).
1060// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1061// in the calculation of the size average.
1062//
1063// Returns a TArrayF of dimension 2:
1064// arr[0]: averaged mean arrival times (default: -1.)
1065// arr[1]: Error (rms) of averaged mean arrival times (default: 0.)
1066//
1067TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeMeanPerSector(const MGeomCam &geom,
1068 const UInt_t sec, MBadPixelsCam *bad)
1069{
1070 const Int_t np = GetSize();
1071
1072 Double_t mean = 0.;
1073 Double_t mean2 = 0.;
1074 Int_t nr = 0;
1075
1076 for (int i=0; i<np; i++)
1077 {
1078 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1079 continue;
1080
1081 const UInt_t sector = geom[i].GetSector();
1082
1083 if (sector != sec)
1084 continue;
1085
1086 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1087 const Float_t time = pix.GetAbsTimeMean();
1088
1089 mean += time;
1090 mean2 += time*time;
1091 nr ++;
1092
1093 }
1094
1095 TArrayF arr(2);
1096 arr[0] = nr ? mean/nr : -1;
1097 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1098 return arr;
1099}
1100
1101// --------------------------------------------------------------------------
1102//
1103// Calculates the average arrival time RMSs for pixel sizes.
1104// The geometry container is used to get the necessary
1105// geometry information (area index).
1106// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1107// in the calculation of the size average.
1108//
1109// Returns a TArrayF of dimension 2:
1110// arr[0]: averaged arrival time RMSs (default: -1.)
1111// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
1112//
1113TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerArea ( const MGeomCam &geom,
1114 const UInt_t ai, MBadPixelsCam *bad)
1115{
1116
1117 const Int_t np = GetSize();
1118
1119 Double_t mean = 0.;
1120 Double_t mean2 = 0.;
1121 Int_t nr = 0;
1122
1123 for (int i=0; i<np; i++)
1124 {
1125 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1126 continue;
1127
1128 const UInt_t aidx = geom[i].GetAidx();
1129
1130 if (ai != aidx)
1131 continue;
1132
1133 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1134 const Float_t rms = pix.GetAbsTimeRms();
1135
1136 mean += rms;
1137 mean2 += rms*rms;
1138 nr ++;
1139
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 return arr;
1146}
1147
1148// --------------------------------------------------------------------------
1149//
1150// Calculates the average arrival time RMSs for camera sectors.
1151// The geometry container is used to get the necessary
1152// geometry information (area index).
1153// If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored
1154// in the calculation of the size average.
1155//
1156// Returns a TArrayF of dimension 2:
1157// arr[0]: averaged arrival time RMSs (default: -1.)
1158// arr[1]: Error (rms) of averaged arrival time RMSs (default: 0.)
1159//
1160TArrayF MCalibrationChargeCam::GetAveragedArrivalTimeRmsPerSector( const MGeomCam &geom, const UInt_t sec, MBadPixelsCam *bad)
1161{
1162 const Int_t np = GetSize();
1163
1164 Double_t mean = 0.;
1165 Double_t mean2 = 0.;
1166 Int_t nr = 0;
1167
1168 for (int i=0; i<np; i++)
1169 {
1170 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1171 continue;
1172
1173 const UInt_t sector = geom[i].GetSector();
1174
1175 if (sector != sec)
1176 continue;
1177
1178 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i];
1179 const Float_t rms = pix.GetAbsTimeRms();
1180
1181 mean += rms;
1182 mean2 += rms*rms;
1183 nr ++;
1184 }
1185
1186 TArrayF arr(2);
1187 arr[0] = nr ? mean/nr : -1;
1188 arr[1] = nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0;
1189
1190 return arr;
1191}
Note: See TracBrowser for help on using the repository browser.