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

Last change on this file since 3455 was 3455, checked in by gaug, 21 years ago
*** empty log message ***
File size: 39.6 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!
19! Author(s): Markus Gaug 11/2003 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MCalibrationChargeCam
29//
30// Hold the whole Calibration results of the camera:
31//
32// 1) MCalibrationChargeCam initializes a TClonesArray whose elements are
33// pointers to MCalibrationChargePix Containers
34// 2) It initializes a pointer to an MCalibrationBlindPix container
35// 3) It initializes a pointer to an MCalibrationPINDiode container
36//
37//
38// The calculated values (types of GetPixelContent) are:
39//
40// --------------------------------------------------------------------------
41//
42// The types are as follows:
43//
44// Fitted values:
45// ==============
46//
47// 0: Fitted Charge
48// 1: Error of fitted Charge
49// 2: Sigma of fitted Charge
50// 3: Error of Sigma of fitted Charge
51//
52// Useful variables derived from the fit results:
53// =============================================
54//
55// 4: Returned probability of Gauss fit to Charge distribution
56// 5: Reduced Sigma of fitted Charge --> sqrt(sigma_Q^2 - PedRMS^2)
57// 6: Error Reduced Sigma of fitted Charge
58// 7: Reduced Sigma per Charge
59// 8: Error of Reduced Sigma per Charge
60//
61// Results of the different calibration methods:
62// =============================================
63//
64// 9: Number of Photo-electrons obtained with the F-Factor method
65// 10: Error on Number of Photo-electrons obtained with the F-Factor method
66// 11: Mean conversion factor obtained with the F-Factor method
67// 12: Error on the mean conversion factor obtained with the F-Factor method
68// 13: Overall F-Factor of the readout obtained with the F-Factor method
69// 14: Error on Overall F-Factor of the readout obtained with the F-Factor method
70// 15: Number of Photons inside Plexiglass obtained with the Blind Pixel method
71// 16: Error on Number of Photons inside Plexiglass obtained with the Blind Pixel method
72// 17: Mean conversion factor obtained with the Blind Pixel method
73// 18: Error on the mean conversion factor obtained with the Blind Pixel method
74// 19: Overall F-Factor of the readout obtained with the Blind Pixel method
75// 20: Error on Overall F-Factor of the readout obtained with the Blind Pixel method
76// 21: Number of Photons outside Plexiglass obtained with the PIN Diode method
77// 22: Error on Number of Photons outside Plexiglass obtained with the PIN Diode method
78// 23: Mean conversion factor obtained with the PIN Diode method
79// 24: Error on the mean conversion factor obtained with the PIN Diode method
80// 25: Overall F-Factor of the readout obtained with the PIN Diode method
81// 26: Error on Overall F-Factor of the readout obtained with the PIN Diode method
82//
83// Localized defects:
84// ==================
85//
86// 27: Excluded Pixels
87// 28: Pixels where the fit did not succeed --> results obtained only from the histograms
88// 29: Pixels with apparently wrong results
89// 30: Pixels with un-expected behavior in the Hi Gain fourier spectrum (e.g. oscillations)
90// 31: Pixels with un-expected behavior in the Lo Gain fourier spectrum (e.g. oscillations)a
91// 32: Number of probable pickup events in the Hi Gain
92// 33: Number of probable pickup events in the Lo Gain
93//
94// Other classifications of pixels:
95// ================================
96//
97// 34: Pixels with saturated Hi-Gain
98//
99// Classification of validity of the calibrations:
100// ===============================================
101//
102// 35: Pixels with valid calibration by the F-Factor-Method
103// 36: Pixels with valid calibration by the Blind Pixel-Method
104// 37: Pixels with valid calibration by the PIN Diode-Method
105//
106// Used Pedestals:
107// ===============
108//
109// 38: Mean Pedestal over the entire range of signal extraction
110// 39: Error on the Mean Pedestal over the entire range of signal extraction
111// 40: Pedestal RMS over the entire range of signal extraction
112// 41: Error on the Pedestal RMS over the entire range of signal extraction
113//
114// Calculated absolute arrival times (very low precision!):
115// ========================================================
116//
117// 42: Absolute Arrival time of the signal
118// 43: RMS of the Absolute Arrival time of the signal
119//
120/////////////////////////////////////////////////////////////////////////////
121#include "MCalibrationChargeCam.h"
122
123#include <TH2.h>
124#include <TCanvas.h>
125#include <TClonesArray.h>
126
127#include "MLog.h"
128#include "MLogManip.h"
129
130#include "MGeomCam.h"
131#include "MGeomPix.h"
132
133#include "MBadPixelsCam.h"
134#include "MBadPixelsPix.h"
135
136#include "MCalibrationChargePix.h"
137#include "MCalibrationChargeBlindPix.h"
138#include "MCalibrationChargePINDiode.h"
139
140ClassImp(MCalibrationChargeCam);
141
142using namespace std;
143
144const Float_t MCalibrationChargeCam::gkAverageQE = 0.18;
145const Float_t MCalibrationChargeCam::gkAverageQEErr = 0.02;
146const Float_t MCalibrationChargeCam::fgConvFFactorRelErrLimit = 0.25;
147// --------------------------------------------------------------------------
148//
149// Default constructor.
150//
151// Creates a TClonesArray of MCalibrationPix containers, initialized to 1 entry
152// Later, a call to MCalibrationChargeCam::InitSize(Int_t size) has to be performed
153//
154// Creates an MCalibrationBlindPix container
155//
156MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
157 : fBlindPixel(NULL),
158 fPINDiode(NULL),
159 fGeomCam(NULL),
160 fBadPixels(NULL),
161 fOffsets(NULL),
162 fSlopes(NULL),
163 fOffvsSlope(NULL)
164{
165 fName = name ? name : "MCalibrationChargeCam";
166 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
167
168 fPixels = new TClonesArray("MCalibrationChargePix",1);
169 fAverageInnerPix = new MCalibrationChargePix("AverageInnerPix","Container of the fit results of the camera average inner pixels");
170 fAverageOuterPix = new MCalibrationChargePix("AverageOuterPix","Container of the fit results of the camera average outer pixels");
171 fAverageInnerBadPix = new MBadPixelsPix("AverageInnerBadPix","Bad Pixel Container of the camera average inner pixels");
172 fAverageOuterBadPix = new MBadPixelsPix("AverageOuterBadPix","Bad Pixel Container of the camera average outer pixels");
173
174 Clear();
175
176 SetAverageQE();
177 SetConvFFactorRelErrLimit();
178}
179
180// --------------------------------------------------------------------------
181//
182// Delete the TClonesArray of MCalibrationPix containers
183// Delete the MCalibrationPINDiode and the MCalibrationBlindPix
184//
185// Delete the histograms if they exist
186//
187MCalibrationChargeCam::~MCalibrationChargeCam()
188{
189
190 //
191 // delete fPixels should delete all Objects stored inside
192 //
193 delete fPixels;
194 delete fAverageInnerPix;
195 delete fAverageOuterPix;
196
197 delete fAverageInnerBadPix;
198 delete fAverageOuterBadPix;
199
200 if (fOffsets)
201 delete fOffsets;
202 if (fSlopes)
203 delete fSlopes;
204 if (fOffvsSlope)
205 delete fOffvsSlope;
206
207}
208
209// -------------------------------------------------------------------
210//
211// This function simply allocates memory via the ROOT command:
212// (TObject**) TStorage::ReAlloc(fCont, newSize * sizeof(TObject*),
213// fSize * sizeof(TObject*));
214// newSize corresponds to size in our case
215// fSize is the old size (in most cases: 1)
216//
217void MCalibrationChargeCam::InitSize(const UInt_t i)
218{
219
220 fPixels->ExpandCreate(i);
221
222}
223
224// --------------------------------------------------------------------------
225//
226// This function returns the current size of the TClonesArray
227// independently if the MCalibrationPix is filled with values or not.
228//
229// It is the size of the array fPixels.
230//
231Int_t MCalibrationChargeCam::GetSize() const
232{
233 return fPixels->GetEntriesFast();
234}
235
236
237// --------------------------------------------------------------------------
238//
239// Get i-th pixel (pixel number)
240//
241MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i)
242{
243 return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
244}
245
246// --------------------------------------------------------------------------
247//
248// Get i-th pixel (pixel number)
249//
250const MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i) const
251{
252 return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
253}
254
255
256// --------------------------------------
257//
258void MCalibrationChargeCam::Clear(Option_t *o)
259{
260
261 fPixels->ForEach(TObject, Clear)();
262 fAverageInnerPix->Clear();
263 fAverageOuterPix->Clear();
264
265 fAverageInnerBadPix->Clear();
266 fAverageOuterBadPix->Clear();
267
268 fNumExcludedPixels = 0;
269 fMeanFluxPhesInnerPixel = 0.;
270 fMeanFluxPhesInnerPixelErr = 0.;
271 fMeanFluxPhesOuterPixel = 0.;
272 fMeanFluxPhesOuterPixelErr = 0.;
273
274 CLRBIT(fFlags,kBlindPixelMethodValid);
275 CLRBIT(fFlags,kFFactorMethodValid);
276 CLRBIT(fFlags,kPINDiodeMethodValid);
277
278 return;
279}
280
281void MCalibrationChargeCam::SetFFactorMethodValid(const Bool_t b)
282{
283 b ? SETBIT(fFlags, kFFactorMethodValid) : CLRBIT(fFlags, kFFactorMethodValid);
284}
285
286void MCalibrationChargeCam::SetBlindPixelMethodValid(const Bool_t b)
287{
288 b ? SETBIT(fFlags, kBlindPixelMethodValid) : CLRBIT(fFlags, kBlindPixelMethodValid);
289}
290
291void MCalibrationChargeCam::SetPINDiodeMethodValid(const Bool_t b)
292{
293 b ? SETBIT(fFlags, kPINDiodeMethodValid) : CLRBIT(fFlags, kPINDiodeMethodValid);
294}
295
296Bool_t MCalibrationChargeCam::IsBlindPixelMethodValid() const
297{
298 return TESTBIT(fFlags,kBlindPixelMethodValid);
299}
300
301Bool_t MCalibrationChargeCam::IsPINDiodeMethodValid() const
302{
303 return TESTBIT(fFlags,kPINDiodeMethodValid);
304}
305
306
307// --------------------------------------------------------------------------
308//
309// Print first the well fitted pixels
310// and then the ones which are not FitValid
311//
312void MCalibrationChargeCam::Print(Option_t *o) const
313{
314
315 *fLog << all << GetDescriptor() << ":" << endl;
316 int id = 0;
317
318 *fLog << all << "Succesfully calibrated pixels:" << endl;
319 *fLog << all << endl;
320
321 TIter Next(fPixels);
322 MCalibrationChargePix *pix;
323 while ((pix=(MCalibrationChargePix*)Next()))
324 {
325
326 if ((*fBadPixels)[pix->GetPixId()].IsCalibrationResultOK() && !pix->IsExcluded())
327 {
328
329 *fLog << all << "Pix " << pix->GetPixId()
330 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
331 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
332 << " Reduced Sigma: " << pix->GetRSigmaCharge()
333 << " Nr Phe's: " << pix->GetPheFFactorMethod()
334 << " Saturated? :" << pix->IsHiGainSaturation()
335 << endl;
336 id++;
337 }
338 }
339
340 *fLog << all << id << " succesful pixels :-))" << endl;
341 id = 0;
342
343 *fLog << all << endl;
344 *fLog << all << "Pixels with errors:" << endl;
345 *fLog << all << endl;
346
347 TIter Next2(fPixels);
348 while ((pix=(MCalibrationChargePix*)Next2()))
349 {
350
351 if (!pix->IsExcluded() && !(*fBadPixels)[pix->GetPixId()].IsCalibrationResultOK())
352 {
353
354
355 *fLog << all << "Pix " << pix->GetPixId()
356 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
357 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
358 << " Reduced Sigma: " << pix->GetRSigmaCharge()
359 << " Nr Phe's: " << pix->GetPheFFactorMethod()
360 << " Saturated? :" << pix->IsHiGainSaturation()
361 << endl;
362 id++;
363 }
364 }
365 *fLog << all << id << " pixels with errors :-((" << endl;
366
367 *fLog << all << endl;
368 *fLog << all << "Pixels with oscillations:" << endl;
369 *fLog << all << endl;
370
371 id = 0;
372
373 TIter Next3(fPixels);
374 while ((pix=(MCalibrationChargePix*)Next3()))
375 {
376
377 if ( (*fBadPixels)[pix->GetPixId()].IsCalibrationOscillating() && !pix->IsExcluded())
378 {
379
380 *fLog << all << "Pix " << pix->GetPixId()
381 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
382 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
383 << " Reduced Sigma: " << pix->GetRSigmaCharge()
384 << " Nr Phe's: " << pix->GetPheFFactorMethod()
385 << " Saturated? :" << pix->IsHiGainSaturation()
386 << endl;
387 id++;
388 }
389 }
390 *fLog << all << id << " Oscillating pixels :-((" << endl;
391
392
393 *fLog << all << endl;
394 *fLog << all << "Excluded pixels:" << endl;
395 *fLog << all << endl;
396
397 id = 0;
398
399 TIter Next4(fPixels);
400 while ((pix=(MCalibrationChargePix*)Next4()))
401 {
402 if (pix->IsExcluded())
403 {
404 *fLog << all << pix->GetPixId() << endl;
405 id++;
406 }
407 }
408 *fLog << all << id << " Excluded pixels " << endl;
409 *fLog << endl;
410 *fLog << all << "Average Inner Pix:"
411 << " Ped. Rms: " << fAverageInnerPix->GetPedRms() << " +- " << fAverageInnerPix->GetPedRmsErr()
412 << " Mean signal: " << fAverageInnerPix->GetMeanCharge() << " +- " << fAverageInnerPix->GetMeanChargeErr()
413 << " Sigma signal: " << fAverageInnerPix->GetSigmaCharge() << " +- "<< fAverageInnerPix->GetSigmaChargeErr()
414 << " Reduced Sigma: " << fAverageInnerPix->GetRSigmaCharge()
415 << " Nr Phe's: " << fAverageInnerPix->GetPheFFactorMethod()
416 << endl;
417 *fLog << all << "Average Outer Pix:"
418 << " Ped. Rms: " << fAverageOuterPix->GetPedRms() << " +- " << fAverageOuterPix->GetPedRmsErr()
419 << " Mean signal: " << fAverageOuterPix->GetMeanCharge() << " +- " << fAverageOuterPix->GetMeanChargeErr()
420 << " Sigma signal: " << fAverageOuterPix->GetSigmaCharge() << " +- "<< fAverageOuterPix->GetSigmaChargeErr()
421 << " Reduced Sigma: " << fAverageOuterPix->GetRSigmaCharge()
422 << " Nr Phe's: " << fAverageOuterPix->GetPheFFactorMethod()
423 << endl;
424
425}
426
427
428// --------------------------------------------------------------------------
429//
430// The types are as follows:
431//
432// Fitted values:
433// ==============
434//
435// 0: Fitted Charge
436// 1: Error of fitted Charge
437// 2: Sigma of fitted Charge
438// 3: Error of Sigma of fitted Charge
439//
440// Useful variables derived from the fit results:
441// =============================================
442//
443// 4: Returned probability of Gauss fit to Charge distribution
444// 5: Reduced Sigma of fitted Charge --> sqrt(sigma_Q^2 - PedRMS^2)
445// 6: Error Reduced Sigma of fitted Charge
446// 7: Reduced Sigma per Charge
447// 8: Error of Reduced Sigma per Charge
448//
449// Results of the different calibration methods:
450// =============================================
451//
452// 9: Number of Photo-electrons obtained with the F-Factor method
453// 10: Error on Number of Photo-electrons obtained with the F-Factor method
454// 11: Mean conversion factor obtained with the F-Factor method
455// 12: Error on the mean conversion factor obtained with the F-Factor method
456// 13: Overall F-Factor of the readout obtained with the F-Factor method
457// 14: Error on Overall F-Factor of the readout obtained with the F-Factor method
458// 15: Number of Photons inside Plexiglass obtained with the Blind Pixel method
459// 16: Error on Number of Photons inside Plexiglass obtained with the Blind Pixel method
460// 17: Mean conversion factor obtained with the Blind Pixel method
461// 18: Error on the mean conversion factor obtained with the Blind Pixel method
462// 19: Overall F-Factor of the readout obtained with the Blind Pixel method
463// 20: Error on Overall F-Factor of the readout obtained with the Blind Pixel method
464// 21: Number of Photons outside Plexiglass obtained with the PIN Diode method
465// 22: Error on Number of Photons outside Plexiglass obtained with the PIN Diode method
466// 23: Mean conversion factor obtained with the PIN Diode method
467// 24: Error on the mean conversion factor obtained with the PIN Diode method
468// 25: Overall F-Factor of the readout obtained with the PIN Diode method
469// 26: Error on Overall F-Factor of the readout obtained with the PIN Diode method
470//
471// Localized defects:
472// ==================
473//
474// 27: Excluded Pixels
475// 28: Pixels where the fit did not succeed --> results obtained only from the histograms
476// 29: Pixels with apparently wrong results
477// 30: Pixels with un-expected behavior in the Hi Gain fourier spectrum (e.g. oscillations)
478// 31: Pixels with un-expected behavior in the Lo Gain fourier spectrum (e.g. oscillations)a
479// 32: Number of probable pickup events in the Hi Gain
480// 33: Number of probable pickup events in the Lo Gain
481//
482// Other classifications of pixels:
483// ================================
484//
485// 34: Pixels with saturated Hi-Gain
486//
487// Classification of validity of the calibrations:
488// ===============================================
489//
490// 35: Pixels with valid calibration by the F-Factor-Method
491// 36: Pixels with valid calibration by the Blind Pixel-Method
492// 37: Pixels with valid calibration by the PIN Diode-Method
493//
494// Used Pedestals:
495// ===============
496//
497// 38: Mean Pedestal over the entire range of signal extraction
498// 39: Error on the Mean Pedestal over the entire range of signal extraction
499// 40: Pedestal RMS over the entire range of signal extraction
500// 41: Error on the Pedestal RMS over the entire range of signal extraction
501//
502// Calculated absolute arrival times (very low precision!):
503// ========================================================
504//
505// 42: Absolute Arrival time of the signal
506// 43: RMS of the Absolute Arrival time of the signal
507//
508Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
509{
510
511 if (idx > GetSize())
512 return kFALSE;
513
514 Float_t area = cam[idx].GetA();
515
516 if (area == 0)
517 return kFALSE;
518
519 switch (type)
520 {
521 case 0:
522 if ((*this)[idx].IsExcluded())
523 return kFALSE;
524 val = (*this)[idx].GetMeanCharge();
525 break;
526 case 1:
527 if ((*this)[idx].IsExcluded())
528 return kFALSE;
529 val = (*this)[idx].GetMeanChargeErr();
530 break;
531 case 2:
532 if ((*this)[idx].IsExcluded())
533 return kFALSE;
534 val = (*this)[idx].GetSigmaCharge();
535 break;
536 case 3:
537 if ((*this)[idx].IsExcluded())
538 return kFALSE;
539 val = (*this)[idx].GetSigmaChargeErr();
540 break;
541 case 4:
542 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
543 return kFALSE;
544 val = (*this)[idx].GetChargeProb();
545 break;
546 case 5:
547 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
548 return kFALSE;
549 if ((*this)[idx].GetRSigmaCharge() == -1.)
550 return kFALSE;
551 val = (*this)[idx].GetRSigmaCharge();
552 break;
553 case 6:
554 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
555 return kFALSE;
556 if ((*this)[idx].GetRSigmaCharge() == -1.)
557 return kFALSE;
558 val = (*this)[idx].GetRSigmaChargeErr();
559 break;
560 case 7:
561 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
562 return kFALSE;
563 if ((*this)[idx].GetRSigmaCharge() == -1.)
564 return kFALSE;
565 val = (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
566 break;
567 case 8:
568 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
569 return kFALSE;
570 if ((*this)[idx].GetRSigmaCharge() == -1.)
571 return kFALSE;
572 // relative error RsigmaCharge square
573 val = (*this)[idx].GetRSigmaChargeErr()* (*this)[idx].GetRSigmaChargeErr()
574 / ((*this)[idx].GetRSigmaCharge() * (*this)[idx].GetRSigmaCharge() );
575 // relative error Charge square
576 val += (*this)[idx].GetMeanChargeErr() * (*this)[idx].GetMeanChargeErr()
577 / ((*this)[idx].GetMeanCharge() * (*this)[idx].GetMeanCharge() );
578 // calculate relative error out of squares
579 val = TMath::Sqrt(val) ;
580 // multiply with value to get absolute error
581 val *= (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
582 break;
583 case 9:
584 if ((*this)[idx].IsExcluded()
585 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
586 || !(*this)[idx].IsFFactorMethodValid())
587 return kFALSE;
588 val = (*this)[idx].GetPheFFactorMethod();
589 break;
590 case 10:
591 if ((*this)[idx].IsExcluded()
592 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
593 || !(*this)[idx].IsFFactorMethodValid())
594 return kFALSE;
595 val = (*this)[idx].GetPheFFactorMethodErr();
596 break;
597 case 11:
598 if ((*this)[idx].IsExcluded()
599 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
600 || !(*this)[idx].IsFFactorMethodValid())
601 return kFALSE;
602 val = (*this)[idx].GetMeanConversionFFactorMethod();
603 break;
604 case 12:
605 if ((*this)[idx].IsExcluded()
606 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
607 || !(*this)[idx].IsFFactorMethodValid())
608 return kFALSE;
609 val = (*this)[idx].GetConversionFFactorMethodErr();
610 break;
611 case 13:
612 if ((*this)[idx].IsExcluded()
613 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
614 || !(*this)[idx].IsFFactorMethodValid())
615 return kFALSE;
616 val = (*this)[idx].GetTotalFFactorFFactorMethod();
617 break;
618 case 14:
619 if ((*this)[idx].IsExcluded()
620 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
621 || !(*this)[idx].IsFFactorMethodValid())
622 return kFALSE;
623 val = (*this)[idx].GetTotalFFactorErrFFactorMethod();
624 break;
625 case 15:
626 if ((*this)[idx].IsExcluded()
627 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
628 || !(*this)[idx].IsBlindPixelMethodValid())
629 return kFALSE;
630 val = fBlindPixel->GetMeanFluxInsidePlexiglass()*area;
631 break;
632 case 16:
633 if ((*this)[idx].IsExcluded()
634 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
635 || !(*this)[idx].IsBlindPixelMethodValid())
636 return kFALSE;
637 val = fBlindPixel->GetMeanFluxErrInsidePlexiglass()*area;
638 break;
639 case 17:
640 if ((*this)[idx].IsExcluded()
641 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
642 || !(*this)[idx].IsBlindPixelMethodValid())
643 return kFALSE;
644 val = (*this)[idx].GetMeanConversionBlindPixelMethod();
645 break;
646 case 18:
647 if ((*this)[idx].IsExcluded()
648 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
649 || !(*this)[idx].IsBlindPixelMethodValid())
650 return kFALSE;
651 val = (*this)[idx].GetConversionBlindPixelMethodErr();
652 break;
653 case 19:
654 if ((*this)[idx].IsExcluded()
655 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
656 || !(*this)[idx].IsBlindPixelMethodValid())
657 return kFALSE;
658 val = (*this)[idx].GetTotalFFactorBlindPixelMethod();
659 break;
660 case 20:
661 if ((*this)[idx].IsExcluded()
662 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
663 || !(*this)[idx].IsBlindPixelMethodValid())
664 return kFALSE;
665 val = (*this)[idx].GetTotalFFactorErrBlindPixelMethod();
666 break;
667 case 21:
668 if ((*this)[idx].IsExcluded()
669 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
670 || !(*this)[idx].IsPINDiodeMethodValid())
671 return kFALSE;
672 val = fPINDiode->GetMeanFluxOutsidePlexiglass()*area;
673 break;
674 case 22:
675 if ((*this)[idx].IsExcluded()
676 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
677 || !(*this)[idx].IsPINDiodeMethodValid())
678 return kFALSE;
679 val = fPINDiode->GetMeanFluxErrOutsidePlexiglass()*area;
680 break;
681 case 23:
682 if ((*this)[idx].IsExcluded()
683 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
684 || !(*this)[idx].IsPINDiodeMethodValid())
685 return kFALSE;
686 val = (*this)[idx].GetMeanConversionPINDiodeMethod();
687 break;
688 case 24:
689 if ((*this)[idx].IsExcluded()
690 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
691 || !(*this)[idx].IsPINDiodeMethodValid())
692 return kFALSE;
693 val = (*this)[idx].GetConversionPINDiodeMethodErr();
694 break;
695 case 25:
696 if ((*this)[idx].IsExcluded()
697 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
698 || !(*this)[idx].IsPINDiodeMethodValid())
699 return kFALSE;
700 val = (*this)[idx].GetTotalFFactorPINDiodeMethod();
701 break;
702 case 26:
703 if ((*this)[idx].IsExcluded()
704 || !(*fBadPixels)[idx].IsCalibrationSignalOK()
705 || !(*this)[idx].IsPINDiodeMethodValid())
706 return kFALSE;
707 val = (*this)[idx].GetTotalFFactorErrPINDiodeMethod();
708 break;
709 case 27:
710 if ((*this)[idx].IsExcluded())
711 val = 1.;
712 else
713 return kFALSE;
714 break;
715 case 28:
716 if ((*this)[idx].IsExcluded())
717 return kFALSE;
718 if (!(*this)[idx].IsFitted())
719 val = 1;
720 else
721 return kFALSE;
722 break;
723 case 29:
724 if ((*this)[idx].IsExcluded())
725 return kFALSE;
726 if (!(*fBadPixels)[idx].IsCalibrationSignalOK())
727 val = 1;
728 else
729 return kFALSE;
730 break;
731 case 30:
732 if ((*this)[idx].IsExcluded())
733 return kFALSE;
734 if ((*fBadPixels)[idx].IsCalibrationOscillating())
735 val = 1;
736 else
737 return kFALSE;
738 break;
739 case 31:
740 if ((*this)[idx].IsExcluded())
741 return kFALSE;
742 if ((*fBadPixels)[idx].IsCalibrationOscillating())
743 val = 1;
744 else
745 return kFALSE;
746 break;
747 case 32:
748 if ((*this)[idx].IsExcluded())
749 return kFALSE;
750 val = (*this)[idx].GetHiGainNumPickup();
751 break;
752 case 33:
753 if ((*this)[idx].IsExcluded())
754 return kFALSE;
755 val = (*this)[idx].GetLoGainNumPickup();
756 break;
757 case 34:
758 if ((*this)[idx].IsExcluded())
759 return kFALSE;
760 val = (*this)[idx].IsHiGainSaturation();
761 break;
762 case 35:
763 if ((*this)[idx].IsExcluded())
764 return kFALSE;
765 if ((*this)[idx].IsFFactorMethodValid())
766 val = 1;
767 else
768 return kFALSE;
769 break;
770 case 36:
771 if ((*this)[idx].IsExcluded())
772 return kFALSE;
773 if ((*this)[idx].IsBlindPixelMethodValid())
774 val = 1;
775 else
776 return kFALSE;
777 break;
778 case 37:
779 if ((*this)[idx].IsExcluded())
780 return kFALSE;
781 if ((*this)[idx].IsPINDiodeMethodValid())
782 val = 1;
783 else
784 return kFALSE;
785 break;
786 case 38:
787 if ((*this)[idx].IsExcluded())
788 return kFALSE;
789 val = (*this)[idx].GetPed();
790 break;
791 case 39:
792 if ((*this)[idx].IsExcluded())
793 return kFALSE;
794 val = (*this)[idx].GetPedErr();
795 break;
796 case 40:
797 if ((*this)[idx].IsExcluded())
798 return kFALSE;
799 val = (*this)[idx].GetPedRms();
800 break;
801 case 41:
802 if ((*this)[idx].IsExcluded())
803 return kFALSE;
804 val = (*this)[idx].GetPedErr()/2.;
805 break;
806 case 42:
807 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
808 return kFALSE;
809 val = (*this)[idx].GetAbsTimeMean();
810 break;
811 case 43:
812 if ((*this)[idx].IsExcluded() || !(*fBadPixels)[idx].IsCalibrationSignalOK())
813 return kFALSE;
814 val = (*this)[idx].GetAbsTimeRms();
815 break;
816 default:
817 return kFALSE;
818 }
819 return val!=-1.;
820}
821
822// --------------------------------------------------------------------------
823//
824// What MHCamera needs in order to draw an individual pixel in the camera
825//
826void MCalibrationChargeCam::DrawPixelContent(Int_t idx) const
827{
828 if (idx == -1)
829 fAverageInnerPix->DrawClone();
830 if (idx == -2)
831 fAverageOuterPix->DrawClone();
832
833 (*this)[idx].DrawClone();
834}
835
836//
837// Calculate the weighted mean of the phe's of all inner and outer pixels, respectively.
838// Bad pixels are excluded from the calculation.
839//
840Bool_t MCalibrationChargeCam::CalcMeanFluxPhotonsFFactorMethod()
841{
842
843 const Float_t avQERelErrSquare = fAverageQEErr * fAverageQEErr
844 / (fAverageQE * fAverageQE );
845
846 Float_t sumweightsinner = 0.;
847 Float_t sumphesinner = 0.;
848 Float_t sumweightsouter = 0.;
849 Float_t sumphesouter = 0.;
850
851 TIter Next(fPixels);
852 MCalibrationChargePix *pix;
853 while ((pix=(MCalibrationChargePix*)Next()))
854 {
855
856 if (!pix->IsFFactorMethodValid())
857 continue;
858
859 const Int_t idx = pix->GetPixId();
860
861 if(!(*fBadPixels)[idx].IsCalibrationResultOK())
862 continue;
863
864 const Float_t nphe = pix->GetPheFFactorMethod();
865 const Float_t npheerr = pix->GetPheFFactorMethodErr();
866 const Float_t ratio = fGeomCam->GetPixRatio(idx);
867
868 if (npheerr > 0.)
869 {
870 //
871 // first the inner pixels:
872 //
873 if (ratio == 1.)
874 {
875 const Float_t weight = 1./npheerr/npheerr;
876 sumweightsinner += weight;
877 sumphesinner += weight*nphe;
878 }
879 else
880 {
881 //
882 // now the outers
883 //
884 const Float_t weight = 1./npheerr/npheerr;
885 sumweightsouter += weight;
886 sumphesouter += weight*nphe;
887 }
888 } /* if npheerr != 0 */
889 } /* while ((pix=(MCalibrationChargePix*)Next())) */
890
891 if (sumweightsinner <= 0. || sumphesinner <= 0.)
892 {
893 *fLog << warn << " Mean number of phe's from inner pixels cannot be calculated: "
894 << " Sum of weights: " << sumweightsinner
895 << " Sum of weighted phes: " << sumphesinner << endl;
896 return kFALSE;
897 }
898 else
899 {
900 fMeanFluxPhesInnerPixel = sumphesinner/sumweightsinner;
901 fMeanFluxPhesInnerPixelErr = TMath::Sqrt(1./sumweightsinner);
902
903 }
904
905 if (sumweightsouter <= 0. || sumphesouter <= 0.)
906 {
907 *fLog << warn << " Mean number of phe's from outer pixels cannot be calculated: "
908 << " Sum of weights or sum of weighted phes is 0. " << endl;
909 }
910 else
911 {
912 fMeanFluxPhesOuterPixel = sumphesouter/sumweightsouter;
913 fMeanFluxPhesOuterPixelErr = TMath::Sqrt(1./sumweightsouter);
914 }
915
916
917 const Float_t meanFluxPhotonsRelErrSquare = fMeanFluxPhesInnerPixelErr * fMeanFluxPhesInnerPixelErr
918 / (fMeanFluxPhesInnerPixel * fMeanFluxPhesInnerPixel);
919
920 fMeanFluxPhotonsInnerPixel = fMeanFluxPhesInnerPixel/fAverageQE;
921 fMeanFluxPhotonsInnerPixelErr = TMath::Sqrt(meanFluxPhotonsRelErrSquare + avQERelErrSquare)
922 * fMeanFluxPhotonsInnerPixel;
923
924 fMeanFluxPhotonsOuterPixel = 4.*fMeanFluxPhotonsInnerPixel;
925 fMeanFluxPhotonsOuterPixelErr = 4.*fMeanFluxPhotonsInnerPixelErr;
926
927 *fLog << inf << " Mean number of photo-electrons from inner pixels (F-Factor Method): "
928 << fMeanFluxPhesInnerPixel << " +- " << fMeanFluxPhesInnerPixelErr << endl;
929
930 *fLog << inf << " Mean number of photons from inner pixels (F-Factor Method): "
931 << fMeanFluxPhotonsInnerPixel << " +- " << fMeanFluxPhotonsInnerPixelErr << endl;
932
933
934
935 return kTRUE;
936}
937
938void MCalibrationChargeCam::ApplyFFactorCalibration()
939{
940
941 const Float_t meanphotRelErrSquare = fMeanFluxPhotonsInnerPixelErr * fMeanFluxPhotonsInnerPixelErr
942 /( fMeanFluxPhotonsInnerPixel * fMeanFluxPhotonsInnerPixel );
943
944 TIter Next(fPixels);
945 MCalibrationChargePix *pix;
946 while ((pix=(MCalibrationChargePix*)Next()))
947 {
948
949 if (!pix->IsFFactorMethodValid())
950 continue;
951
952 const Int_t idx = pix->GetPixId();
953
954 if(!(*fBadPixels)[idx].IsCalibrationResultOK())
955 continue;
956
957 const Float_t ratio = fGeomCam->GetPixRatio(idx);
958 //
959 // Calculate the conversion factor between PHOTONS and FADC counts
960 //
961 // Nphot = Nphe / avQE
962 // conv = Nphot / FADC counts
963 //
964 Float_t conv;
965
966 if (ratio == 1.)
967 conv = fMeanFluxPhotonsInnerPixel / pix->GetMeanCharge();
968 else
969 conv = fMeanFluxPhotonsOuterPixel / pix->GetMeanCharge();
970
971 if (conv <= 0.)
972 {
973 pix->SetFFactorMethodValid(kFALSE);
974 continue;
975 }
976
977 const Float_t chargeRelErrSquare = pix->GetMeanChargeErr() * pix->GetMeanChargeErr()
978 / ( pix->GetMeanCharge() * pix->GetMeanCharge());
979 const Float_t rsigmaChargeRelErrSquare = pix->GetRSigmaChargeErr() * pix->GetRSigmaChargeErr()
980 / (pix->GetRSigmaCharge() * pix->GetRSigmaCharge()) ;
981
982 const Float_t convrelerr = TMath::Sqrt(meanphotRelErrSquare + chargeRelErrSquare);
983
984 if (convrelerr > fConvFFactorRelErrLimit)
985 {
986 *fLog << warn << GetDescriptor() << ": Conversion Factor F-Factor Method Rel. Error: "
987 << convrelerr << " above limit of: " << fConvFFactorRelErrLimit
988 << " in pixel: " << idx << endl;
989 pix->SetFFactorMethodValid(kFALSE);
990 continue;
991 }
992
993 //
994 // Calculate the Total F-Factor of the camera (in photons)
995 //
996 const Float_t totalFFactor = (pix->GetRSigmaCharge()/pix->GetMeanCharge())
997 *TMath::Sqrt(fMeanFluxPhotonsInnerPixel);
998
999 //
1000 // Calculate the error of the Total F-Factor of the camera ( in photons )
1001 //
1002 const Float_t totalFFactorErr = TMath::Sqrt( rsigmaChargeRelErrSquare
1003 + chargeRelErrSquare
1004 + meanphotRelErrSquare );
1005
1006 pix->SetConversionFFactorMethod(conv,
1007 convrelerr*conv,
1008 totalFFactor*TMath::Sqrt(conv));
1009
1010 pix->SetTotalFFactorFFactorMethod( totalFFactor );
1011 pix->SetTotalFFactorErrFFactorMethod(totalFFactorErr);
1012 pix->SetFFactorMethodValid();
1013 }
1014}
1015
1016
1017
1018void MCalibrationChargeCam::ApplyBlindPixelCalibration()
1019{
1020
1021 Float_t flux = fBlindPixel->GetMeanFluxInsidePlexiglass();
1022 Float_t fluxerr = fBlindPixel->GetMeanFluxErrInsidePlexiglass();
1023
1024 TIter Next(fPixels);
1025 MCalibrationChargePix *pix;
1026 while ((pix=(MCalibrationChargePix*)Next()))
1027 {
1028
1029 if((*fBadPixels)[pix->GetPixId()].IsCalibrationResultOK())
1030 {
1031
1032 const Int_t idx = pix->GetPixId();
1033
1034 const Float_t charge = pix->GetMeanCharge();
1035 const Float_t area = (*fGeomCam)[idx].GetA();
1036 const Float_t chargeerr = pix->GetMeanChargeErr();
1037
1038 const Float_t nphot = flux * area;
1039 const Float_t nphoterr = fluxerr * area;
1040 const Float_t conversion = nphot/charge;
1041 Float_t conversionerr;
1042
1043 conversionerr = nphoterr/charge
1044 * nphoterr/charge ;
1045 conversionerr += chargeerr/charge
1046 * chargeerr/charge
1047 * conversion*conversion;
1048 conversionerr = TMath::Sqrt(conversionerr);
1049
1050 const Float_t conversionsigma = 0.;
1051
1052 pix->SetConversionBlindPixelMethod(conversion, conversionerr, conversionsigma);
1053
1054 if (conversionerr/conversion < 0.1)
1055 pix->SetBlindPixelMethodValid();
1056 }
1057 }
1058}
1059
1060
1061void MCalibrationChargeCam::ApplyPINDiodeCalibration()
1062{
1063
1064 Float_t flux = fPINDiode->GetMeanFluxOutsidePlexiglass();
1065 Float_t fluxerr = fPINDiode->GetMeanFluxErrOutsidePlexiglass();
1066
1067 TIter Next(fPixels);
1068 MCalibrationChargePix *pix;
1069 while ((pix=(MCalibrationChargePix*)Next()))
1070 {
1071
1072 if((*fBadPixels)[pix->GetPixId()].IsCalibrationResultOK())
1073 {
1074
1075 const Int_t idx = pix->GetPixId();
1076
1077 const Float_t charge = pix->GetMeanCharge();
1078 const Float_t area = (*fGeomCam)[idx].GetA();
1079 const Float_t chargeerr = pix->GetMeanChargeErr();
1080
1081 const Float_t nphot = flux * area;
1082 const Float_t nphoterr = fluxerr * area;
1083 const Float_t conversion = nphot/charge;
1084
1085 Float_t conversionerr;
1086
1087 conversionerr = nphoterr/charge
1088 * nphoterr/charge ;
1089 conversionerr += chargeerr/charge
1090 * chargeerr/charge
1091 * conversion*conversion;
1092 if (conversionerr > 0.)
1093 conversionerr = TMath::Sqrt(conversionerr);
1094
1095 const Float_t conversionsigma = 0.;
1096
1097 pix->SetConversionPINDiodeMethod(conversion, conversionerr, conversionsigma);
1098
1099 if (conversionerr/conversion < 0.1)
1100 pix->SetPINDiodeMethodValid();
1101
1102 }
1103 }
1104}
1105
1106
1107
1108Bool_t MCalibrationChargeCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
1109{
1110
1111 mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
1112 err = (*this)[ipx].GetConversionBlindPixelMethodErr();
1113 sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
1114
1115 return kTRUE;
1116}
1117
1118
1119Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
1120{
1121
1122 Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
1123
1124 if (conv < 0.)
1125 return kFALSE;
1126
1127 mean = conv;
1128 err = (*this)[ipx].GetConversionFFactorMethodErr();
1129 sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
1130
1131 return kTRUE;
1132}
1133
1134
1135//-----------------------------------------------------------------------------------
1136//
1137// Calculates the conversion factor between the integral of FADCs slices
1138// (as defined in the signal extractor MExtractSignal.cc)
1139// and the number of photons reaching the plexiglass for one Inner Pixel
1140//
1141// FIXME: The PINDiode is still not working and so is the code
1142//
1143Bool_t MCalibrationChargeCam::GetConversionFactorPINDiode(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
1144{
1145
1146 mean = (*this)[ipx].GetMeanConversionPINDiodeMethod();
1147 err = (*this)[ipx].GetConversionPINDiodeMethodErr();
1148 sigma = (*this)[ipx].GetSigmaConversionPINDiodeMethod();
1149
1150 return kFALSE;
1151
1152}
1153
1154//-----------------------------------------------------------------------------------
1155//
1156// Calculates the best combination of the three used methods possible
1157// between the integral of FADCs slices
1158// (as defined in the signal extractor MExtractSignal.cc)
1159// and the number of photons reaching one Inner Pixel.
1160// The procedure is not yet defined.
1161//
1162// FIXME: The PINDiode is still not working and so is the code
1163//
1164Bool_t MCalibrationChargeCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
1165{
1166 return kFALSE;
1167
1168}
1169
1170/*
1171void MCalibrationChargeCam::DrawHiLoFits()
1172{
1173
1174 if (!fOffsets)
1175 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
1176 if (!fSlopes)
1177 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
1178 if (!fOffvsSlope)
1179 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
1180
1181 TIter Next(fPixels);
1182 MCalibrationPix *pix;
1183 MHCalibrationPixel *hist;
1184 while ((pix=(MCalibrationPix*)Next()))
1185 {
1186 hist = pix->GetHist();
1187 hist->FitHiGainvsLoGain();
1188 fOffsets->Fill(hist->GetOffset(),1.);
1189 fSlopes->Fill(hist->GetSlope(),1.);
1190 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
1191 }
1192
1193 TCanvas *c1 = new TCanvas();
1194
1195 c1->Divide(1,3);
1196 c1->cd(1);
1197 fOffsets->Draw();
1198 gPad->Modified();
1199 gPad->Update();
1200
1201 c1->cd(2);
1202 fSlopes->Draw();
1203 gPad->Modified();
1204 gPad->Update();
1205
1206 c1->cd(3);
1207 fOffvsSlope->Draw("col1");
1208 gPad->Modified();
1209 gPad->Update();
1210}
1211
1212*/
Note: See TracBrowser for help on using the repository browser.