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

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