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

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