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

Last change on this file since 3361 was 3354, checked in by gaug, 21 years ago
*** empty log message ***
File size: 31.4 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15
16! *
17!
18!
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 "MCalibrationChargePix.h"
134#include "MCalibrationChargeBlindPix.h"
135#include "MCalibrationChargePINDiode.h"
136
137
138ClassImp(MCalibrationChargeCam);
139
140using namespace std;
141
142// --------------------------------------------------------------------------
143//
144// Default constructor.
145//
146// Creates a TClonesArray of MCalibrationPix containers, initialized to 1 entry
147// Later, a call to MCalibrationChargeCam::InitSize(Int_t size) has to be performed
148//
149// Creates an MCalibrationBlindPix container
150//
151MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
152 : fBlindPixel(NULL),
153 fPINDiode(NULL),
154 fOffsets(NULL),
155 fSlopes(NULL),
156 fOffvsSlope(NULL)
157{
158 fName = name ? name : "MCalibrationChargeCam";
159 fTitle = title ? title : "Storage container for the Calibration Information in the camera";
160
161 fPixels = new TClonesArray("MCalibrationChargePix",1);
162 fAverageInnerPix = new MCalibrationChargePix("AverageInnerPix","Container of the fit results of the camera average inner pixels");
163 fAverageOuterPix = new MCalibrationChargePix("AverageOuterPix","Container of the fit results of the camera average outer pixels");
164
165 Clear();
166}
167
168// --------------------------------------------------------------------------
169//
170// Delete the TClonesArray of MCalibrationPix containers
171// Delete the MCalibrationPINDiode and the MCalibrationBlindPix
172//
173// Delete the histograms if they exist
174//
175MCalibrationChargeCam::~MCalibrationChargeCam()
176{
177
178 //
179 // delete fPixels should delete all Objects stored inside
180 //
181 delete fPixels;
182 delete fAverageInnerPix;
183 delete fAverageOuterPix;
184
185 if (fOffsets)
186 delete fOffsets;
187 if (fSlopes)
188 delete fSlopes;
189 if (fOffvsSlope)
190 delete fOffvsSlope;
191
192}
193
194// -------------------------------------------------------------------
195//
196// This function simply allocates memory via the ROOT command:
197// (TObject**) TStorage::ReAlloc(fCont, newSize * sizeof(TObject*),
198// fSize * sizeof(TObject*));
199// newSize corresponds to size in our case
200// fSize is the old size (in most cases: 1)
201//
202void MCalibrationChargeCam::InitSize(const UInt_t i)
203{
204
205 fPixels->ExpandCreate(i);
206
207}
208
209// --------------------------------------------------------------------------
210//
211// This function returns the current size of the TClonesArray
212// independently if the MCalibrationPix is filled with values or not.
213//
214// It is the size of the array fPixels.
215//
216Int_t MCalibrationChargeCam::GetSize() const
217{
218 return fPixels->GetEntriesFast();
219}
220
221
222// --------------------------------------------------------------------------
223//
224// Get i-th pixel (pixel number)
225//
226MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i)
227{
228 return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
229}
230
231// --------------------------------------------------------------------------
232//
233// Get i-th pixel (pixel number)
234//
235const MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i) const
236{
237 return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
238}
239
240
241// --------------------------------------
242//
243void MCalibrationChargeCam::Clear(Option_t *o)
244{
245
246 fPixels->ForEach(TObject, Clear)();
247 fAverageInnerPix->Clear();
248 fAverageOuterPix->Clear();
249
250 fNumExcludedPixels = 0;
251
252 CLRBIT(fFlags,kBlindPixelMethodValid);
253 CLRBIT(fFlags,kPINDiodeMethodValid);
254
255 return;
256}
257
258void MCalibrationChargeCam::SetBlindPixelMethodValid(const Bool_t b)
259{
260 b ? SETBIT(fFlags, kBlindPixelMethodValid) : CLRBIT(fFlags, kBlindPixelMethodValid);
261}
262
263void MCalibrationChargeCam::SetPINDiodeMethodValid(const Bool_t b)
264{
265 b ? SETBIT(fFlags, kPINDiodeMethodValid) : CLRBIT(fFlags, kPINDiodeMethodValid);
266}
267
268Bool_t MCalibrationChargeCam::IsBlindPixelMethodValid() const
269{
270 return TESTBIT(fFlags,kBlindPixelMethodValid);
271}
272
273Bool_t MCalibrationChargeCam::IsPINDiodeMethodValid() const
274{
275 return TESTBIT(fFlags,kPINDiodeMethodValid);
276}
277
278
279// --------------------------------------------------------------------------
280//
281// Print first the well fitted pixels
282// and then the ones which are not FitValid
283//
284void MCalibrationChargeCam::Print(Option_t *o) const
285{
286
287 *fLog << all << GetDescriptor() << ":" << endl;
288 int id = 0;
289
290 *fLog << all << "Succesfully calibrated pixels:" << endl;
291 *fLog << all << endl;
292
293 TIter Next(fPixels);
294 MCalibrationChargePix *pix;
295 while ((pix=(MCalibrationChargePix*)Next()))
296 {
297
298 if (pix->IsChargeValid() && !pix->IsExcluded() && !pix->IsOscillating())
299 {
300
301 *fLog << all << "Pix " << pix->GetPixId()
302 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
303 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
304 << " Reduced Sigma: " << pix->GetRSigmaCharge()
305 << " Nr Phe's: " << pix->GetPheFFactorMethod()
306 << endl;
307 id++;
308 }
309 }
310
311 *fLog << all << id << " succesful pixels :-))" << endl;
312 id = 0;
313
314 *fLog << all << endl;
315 *fLog << all << "Pixels with errors:" << endl;
316 *fLog << all << endl;
317
318 TIter Next2(fPixels);
319 while ((pix=(MCalibrationChargePix*)Next2()))
320 {
321
322 if (!pix->IsExcluded() && !pix->IsChargeValid())
323 {
324
325
326 *fLog << all << "Pix " << pix->GetPixId()
327 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
328 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
329 << " Reduced Sigma: " << pix->GetRSigmaCharge()
330 << " Nr Phe's: " << pix->GetPheFFactorMethod()
331 << endl;
332 id++;
333 }
334 }
335 *fLog << all << id << " pixels with errors :-((" << endl;
336
337 *fLog << all << endl;
338 *fLog << all << "Pixels with oscillations:" << endl;
339 *fLog << all << endl;
340
341 id = 0;
342
343 TIter Next3(fPixels);
344 while ((pix=(MCalibrationChargePix*)Next3()))
345 {
346
347 if ( pix->IsOscillating() && !pix->IsExcluded())
348 {
349
350 *fLog << all << "Pix " << pix->GetPixId()
351 << ": Ped. Rms: " << pix->GetPedRms() << " +- " << pix->GetPedRmsErr()
352 << " Mean signal: " << pix->GetMeanCharge() << " +- " << pix->GetSigmaCharge()
353 << " Reduced Sigma: " << pix->GetRSigmaCharge()
354 << " Nr Phe's: " << pix->GetPheFFactorMethod()
355 << endl;
356 id++;
357 }
358 }
359 *fLog << all << id << " Oscillating pixels :-((" << endl;
360
361
362 *fLog << all << endl;
363 *fLog << all << "Excluded pixels:" << endl;
364 *fLog << all << endl;
365
366 id = 0;
367
368 TIter Next4(fPixels);
369 while ((pix=(MCalibrationChargePix*)Next4()))
370 {
371 if (pix->IsExcluded())
372 {
373 *fLog << all << pix->GetPixId() << endl;
374 id++;
375 }
376 }
377 *fLog << all << id << " Excluded pixels " << endl;
378 *fLog << endl;
379 *fLog << all << "Average Inner Pix:"
380 << " Ped. Rms: " << fAverageInnerPix->GetPedRms() << " +- " << fAverageInnerPix->GetPedRmsErr()
381 << " Mean signal: " << fAverageInnerPix->GetMeanCharge() << " +- " << fAverageInnerPix->GetMeanChargeErr()
382 << " Sigma signal: " << fAverageInnerPix->GetSigmaCharge() << " +- "<< fAverageInnerPix->GetSigmaChargeErr()
383 << " Reduced Sigma: " << fAverageInnerPix->GetRSigmaCharge()
384 << " Nr Phe's: " << fAverageInnerPix->GetPheFFactorMethod()
385 << endl;
386 *fLog << all << "Average Outer Pix:"
387 << " Ped. Rms: " << fAverageOuterPix->GetPedRms() << " +- " << fAverageOuterPix->GetPedRmsErr()
388 << " Mean signal: " << fAverageOuterPix->GetMeanCharge() << " +- " << fAverageOuterPix->GetMeanChargeErr()
389 << " Sigma signal: " << fAverageOuterPix->GetSigmaCharge() << " +- "<< fAverageOuterPix->GetSigmaChargeErr()
390 << " Reduced Sigma: " << fAverageOuterPix->GetRSigmaCharge()
391 << " Nr Phe's: " << fAverageOuterPix->GetPheFFactorMethod()
392 << endl;
393
394}
395
396
397// --------------------------------------------------------------------------
398//
399// The types are as follows:
400//
401// Fitted values:
402// ==============
403//
404// 0: Fitted Charge
405// 1: Error of fitted Charge
406// 2: Sigma of fitted Charge
407// 3: Error of Sigma of fitted Charge
408//
409// Useful variables derived from the fit results:
410// =============================================
411//
412// 4: Returned probability of Gauss fit to Charge distribution
413// 5: Reduced Sigma of fitted Charge --> sqrt(sigma_Q^2 - PedRMS^2)
414// 6: Error Reduced Sigma of fitted Charge
415// 7: Reduced Sigma per Charge
416// 8: Error of Reduced Sigma per Charge
417//
418// Results of the different calibration methods:
419// =============================================
420//
421// 9: Number of Photo-electrons obtained with the F-Factor method
422// 10: Error on Number of Photo-electrons obtained with the F-Factor method
423// 11: Mean conversion factor obtained with the F-Factor method
424// 12: Error on the mean conversion factor obtained with the F-Factor method
425// 13: Overall F-Factor of the readout obtained with the F-Factor method
426// 14: Error on Overall F-Factor of the readout obtained with the F-Factor method
427// 15: Number of Photons inside Plexiglass obtained with the Blind Pixel method
428// 16: Error on Number of Photons inside Plexiglass obtained with the Blind Pixel method
429// 17: Mean conversion factor obtained with the Blind Pixel method
430// 18: Error on the mean conversion factor obtained with the Blind Pixel method
431// 19: Overall F-Factor of the readout obtained with the Blind Pixel method
432// 20: Error on Overall F-Factor of the readout obtained with the Blind Pixel method
433// 21: Number of Photons outside Plexiglass obtained with the PIN Diode method
434// 22: Error on Number of Photons outside Plexiglass obtained with the PIN Diode method
435// 23: Mean conversion factor obtained with the PIN Diode method
436// 24: Error on the mean conversion factor obtained with the PIN Diode method
437// 25: Overall F-Factor of the readout obtained with the PIN Diode method
438// 26: Error on Overall F-Factor of the readout obtained with the PIN Diode method
439//
440// Localized defects:
441// ==================
442//
443// 27: Excluded Pixels
444// 28: Pixels where the fit did not succeed --> results obtained only from the histograms
445// 29: Pixels with apparently wrong results
446// 30: Pixels with un-expected behavior in the Hi Gain fourier spectrum (e.g. oscillations)
447// 31: Pixels with un-expected behavior in the Lo Gain fourier spectrum (e.g. oscillations)a
448// 32: Number of probable pickup events in the Hi Gain
449// 33: Number of probable pickup events in the Lo Gain
450//
451// Other classifications of pixels:
452// ================================
453//
454// 34: Pixels with saturated Hi-Gain
455//
456// Classification of validity of the calibrations:
457// ===============================================
458//
459// 35: Pixels with valid calibration by the F-Factor-Method
460// 36: Pixels with valid calibration by the Blind Pixel-Method
461// 37: Pixels with valid calibration by the PIN Diode-Method
462//
463// Used Pedestals:
464// ===============
465//
466// 38: Mean Pedestal over the entire range of signal extraction
467// 39: Error on the Mean Pedestal over the entire range of signal extraction
468// 40: Pedestal RMS over the entire range of signal extraction
469// 41: Error on the Pedestal RMS over the entire range of signal extraction
470//
471// Calculated absolute arrival times (very low precision!):
472// ========================================================
473//
474// 42: Absolute Arrival time of the signal
475// 43: RMS of the Absolute Arrival time of the signal
476//
477Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
478{
479
480 if (idx > GetSize())
481 return kFALSE;
482
483 Float_t area = cam[idx].GetA();
484
485 if (area == 0)
486 return kFALSE;
487
488 switch (type)
489 {
490 case 0:
491 if ((*this)[idx].IsExcluded())
492 return kFALSE;
493 val = (*this)[idx].GetMeanCharge();
494 break;
495 case 1:
496 if ((*this)[idx].IsExcluded())
497 return kFALSE;
498 val = (*this)[idx].GetMeanChargeErr();
499 break;
500 case 2:
501 if ((*this)[idx].IsExcluded())
502 return kFALSE;
503 val = (*this)[idx].GetSigmaCharge();
504 break;
505 case 3:
506 if ((*this)[idx].IsExcluded())
507 return kFALSE;
508 val = (*this)[idx].GetSigmaChargeErr();
509 break;
510 case 4:
511 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
512 return kFALSE;
513 val = (*this)[idx].GetChargeProb();
514 break;
515 case 5:
516 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
517 return kFALSE;
518 if ((*this)[idx].GetRSigmaCharge() == -1.)
519 return kFALSE;
520 val = (*this)[idx].GetRSigmaCharge();
521 break;
522 case 6:
523 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
524 return kFALSE;
525 if ((*this)[idx].GetRSigmaCharge() == -1.)
526 return kFALSE;
527 val = (*this)[idx].GetRSigmaChargeErr();
528 break;
529 case 7:
530 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
531 return kFALSE;
532 if ((*this)[idx].GetRSigmaCharge() == -1.)
533 return kFALSE;
534 val = (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
535 break;
536 case 8:
537 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
538 return kFALSE;
539 if ((*this)[idx].GetRSigmaCharge() == -1.)
540 return kFALSE;
541 // relative error RsigmaCharge square
542 val = (*this)[idx].GetRSigmaChargeErr()* (*this)[idx].GetRSigmaChargeErr()
543 / ((*this)[idx].GetRSigmaCharge() * (*this)[idx].GetRSigmaCharge() );
544 // relative error Charge square
545 val += (*this)[idx].GetMeanChargeErr() * (*this)[idx].GetMeanChargeErr()
546 / ((*this)[idx].GetMeanCharge() * (*this)[idx].GetMeanCharge() );
547 // calculate relative error out of squares
548 val = TMath::Sqrt(val) ;
549 // multiply with value to get absolute error
550 val *= (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
551 break;
552 case 9:
553 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
554 return kFALSE;
555 val = (*this)[idx].GetPheFFactorMethod();
556 break;
557 case 10:
558 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
559 return kFALSE;
560 val = (*this)[idx].GetPheFFactorMethodErr();
561 break;
562 case 11:
563 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
564 return kFALSE;
565 val = (*this)[idx].GetMeanConversionFFactorMethod();
566 break;
567 case 12:
568 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
569 return kFALSE;
570 val = (*this)[idx].GetConversionFFactorMethodErr();
571 break;
572 case 13:
573 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
574 return kFALSE;
575 val = (*this)[idx].GetTotalFFactorFFactorMethod();
576 break;
577 case 14:
578 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsFFactorMethodValid())
579 return kFALSE;
580 val = (*this)[idx].GetTotalFFactorErrFFactorMethod();
581 break;
582 case 15:
583 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
584 return kFALSE;
585 val = fBlindPixel->GetMeanFluxInsidePlexiglass()*area;
586 break;
587 case 16:
588 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
589 return kFALSE;
590 val = fBlindPixel->GetMeanFluxErrInsidePlexiglass()*area;
591 break;
592 case 17:
593 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
594 return kFALSE;
595 val = (*this)[idx].GetMeanConversionBlindPixelMethod();
596 break;
597 case 18:
598 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
599 return kFALSE;
600 val = (*this)[idx].GetConversionBlindPixelMethodErr();
601 break;
602 case 19:
603 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
604 return kFALSE;
605 val = (*this)[idx].GetTotalFFactorBlindPixelMethod();
606 break;
607 case 20:
608 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsBlindPixelMethodValid())
609 return kFALSE;
610 val = (*this)[idx].GetTotalFFactorErrBlindPixelMethod();
611 break;
612 case 21:
613 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
614 return kFALSE;
615 val = fPINDiode->GetMeanFluxOutsidePlexiglass()*area;
616 break;
617 case 22:
618 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
619 return kFALSE;
620 val = fPINDiode->GetMeanFluxErrOutsidePlexiglass()*area;
621 break;
622 case 23:
623 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
624 return kFALSE;
625 val = (*this)[idx].GetMeanConversionPINDiodeMethod();
626 break;
627 case 24:
628 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
629 return kFALSE;
630 val = (*this)[idx].GetConversionPINDiodeMethodErr();
631 break;
632 case 25:
633 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
634 return kFALSE;
635 val = (*this)[idx].GetTotalFFactorPINDiodeMethod();
636 break;
637 case 26:
638 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid() || !(*this)[idx].IsPINDiodeMethodValid())
639 return kFALSE;
640 val = (*this)[idx].GetTotalFFactorErrPINDiodeMethod();
641 break;
642 case 27:
643 if ((*this)[idx].IsExcluded())
644 val = 1.;
645 else
646 return kFALSE;
647 break;
648 case 28:
649 if ((*this)[idx].IsExcluded())
650 return kFALSE;
651 if (!(*this)[idx].IsFitted())
652 val = 1;
653 else
654 return kFALSE;
655 break;
656 case 29:
657 if ((*this)[idx].IsExcluded())
658 return kFALSE;
659 if (!(*this)[idx].IsChargeValid())
660 val = 1;
661 else
662 return kFALSE;
663 break;
664 case 30:
665 if ((*this)[idx].IsExcluded())
666 return kFALSE;
667 if ((*this)[idx].IsHiGainOscillating())
668 val = 1;
669 else
670 return kFALSE;
671 break;
672 case 31:
673 if ((*this)[idx].IsExcluded())
674 return kFALSE;
675 if ((*this)[idx].IsLoGainOscillating())
676 val = 1;
677 else
678 return kFALSE;
679 break;
680 case 32:
681 if ((*this)[idx].IsExcluded())
682 return kFALSE;
683 val = (*this)[idx].GetHiGainNumPickup();
684 break;
685 case 33:
686 if ((*this)[idx].IsExcluded())
687 return kFALSE;
688 val = (*this)[idx].GetLoGainNumPickup();
689 break;
690 case 34:
691 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
692 return kFALSE;
693 if ((*this)[idx].IsHiGainSaturation())
694 val = 1;
695 else
696 return kFALSE;
697 break;
698 case 35:
699 if ((*this)[idx].IsExcluded())
700 return kFALSE;
701 if ((*this)[idx].IsFFactorMethodValid())
702 val = 1;
703 else
704 return kFALSE;
705 break;
706 case 36:
707 if ((*this)[idx].IsExcluded())
708 return kFALSE;
709 if ((*this)[idx].IsBlindPixelMethodValid())
710 val = 1;
711 else
712 return kFALSE;
713 break;
714 case 37:
715 if ((*this)[idx].IsExcluded())
716 return kFALSE;
717 if ((*this)[idx].IsPINDiodeMethodValid())
718 val = 1;
719 else
720 return kFALSE;
721 break;
722 case 38:
723 if ((*this)[idx].IsExcluded())
724 return kFALSE;
725 val = (*this)[idx].GetPed();
726 break;
727 case 39:
728 if ((*this)[idx].IsExcluded())
729 return kFALSE;
730 val = (*this)[idx].GetPedErr();
731 break;
732 case 40:
733 if ((*this)[idx].IsExcluded())
734 return kFALSE;
735 val = (*this)[idx].GetPedRms();
736 break;
737 case 41:
738 if ((*this)[idx].IsExcluded())
739 return kFALSE;
740 val = (*this)[idx].GetPedErr()/2.;
741 break;
742 case 42:
743 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
744 return kFALSE;
745 val = (*this)[idx].GetAbsTimeMean();
746 break;
747 case 43:
748 if ((*this)[idx].IsExcluded() || !(*this)[idx].IsChargeValid())
749 return kFALSE;
750 val = (*this)[idx].GetAbsTimeRms();
751 break;
752 default:
753 return kFALSE;
754 }
755 return val!=-1.;
756}
757
758// --------------------------------------------------------------------------
759//
760// What MHCamera needs in order to draw an individual pixel in the camera
761//
762void MCalibrationChargeCam::DrawPixelContent(Int_t idx) const
763{
764 if (idx == -1)
765 fAverageInnerPix->DrawClone();
766 if (idx == -2)
767 fAverageOuterPix->DrawClone();
768
769 (*this)[idx].DrawClone();
770}
771
772void MCalibrationChargeCam::ApplyBlindPixelCalibration()
773{
774
775 Float_t flux = fBlindPixel->GetMeanFluxInsidePlexiglass();
776 Float_t fluxerr = fBlindPixel->GetMeanFluxErrInsidePlexiglass();
777
778 TIter Next(fPixels);
779 MCalibrationChargePix *pix;
780 while ((pix=(MCalibrationChargePix*)Next()))
781 {
782
783 if(pix->IsChargeValid())
784 {
785
786 const Int_t idx = pix->GetPixId();
787
788 const Float_t charge = pix->GetMeanCharge();
789 const Float_t area = (*fGeomCam)[idx].GetA();
790 const Float_t chargeerr = pix->GetMeanChargeErr();
791
792 const Float_t nphot = flux * area;
793 const Float_t nphoterr = fluxerr * area;
794 const Float_t conversion = nphot/charge;
795 Float_t conversionerr;
796
797 conversionerr = nphoterr/charge
798 * nphoterr/charge ;
799 conversionerr += chargeerr/charge
800 * chargeerr/charge
801 * conversion*conversion;
802 conversionerr = TMath::Sqrt(conversionerr);
803
804 const Float_t conversionsigma = 0.;
805
806 pix->SetConversionBlindPixelMethod(conversion, conversionerr, conversionsigma);
807
808 if (conversionerr/conversion < 0.1)
809 pix->SetBlindPixelMethodValid();
810 }
811 }
812}
813
814
815void MCalibrationChargeCam::ApplyPINDiodeCalibration()
816{
817
818 Float_t flux = fPINDiode->GetMeanFluxOutsidePlexiglass();
819 Float_t fluxerr = fPINDiode->GetMeanFluxErrOutsidePlexiglass();
820
821 TIter Next(fPixels);
822 MCalibrationChargePix *pix;
823 while ((pix=(MCalibrationChargePix*)Next()))
824 {
825
826 if (pix->IsChargeValid())
827 {
828
829 const Int_t idx = pix->GetPixId();
830
831 const Float_t charge = pix->GetMeanCharge();
832 const Float_t area = (*fGeomCam)[idx].GetA();
833 const Float_t chargeerr = pix->GetMeanChargeErr();
834
835 const Float_t nphot = flux * area;
836 const Float_t nphoterr = fluxerr * area;
837 const Float_t conversion = nphot/charge;
838
839 Float_t conversionerr;
840
841 conversionerr = nphoterr/charge
842 * nphoterr/charge ;
843 conversionerr += chargeerr/charge
844 * chargeerr/charge
845 * conversion*conversion;
846 if (conversionerr > 0.)
847 conversionerr = TMath::Sqrt(conversionerr);
848
849 const Float_t conversionsigma = 0.;
850
851 pix->SetConversionPINDiodeMethod(conversion, conversionerr, conversionsigma);
852
853 if (conversionerr/conversion < 0.1)
854 pix->SetPINDiodeMethodValid();
855
856 }
857 }
858}
859
860
861
862Bool_t MCalibrationChargeCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
863{
864
865 mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
866 err = (*this)[ipx].GetConversionBlindPixelMethodErr();
867 sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
868
869 return kTRUE;
870}
871
872
873Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
874{
875
876 Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
877
878 if (conv < 0.)
879 return kFALSE;
880
881 mean = conv;
882 err = (*this)[ipx].GetConversionFFactorMethodErr();
883 sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
884
885 return kTRUE;
886}
887
888
889//-----------------------------------------------------------------------------------
890//
891// Calculates the conversion factor between the integral of FADCs slices
892// (as defined in the signal extractor MExtractSignal.cc)
893// and the number of photons reaching the plexiglass for one Inner Pixel
894//
895// FIXME: The PINDiode is still not working and so is the code
896//
897Bool_t MCalibrationChargeCam::GetConversionFactorPINDiode(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
898{
899
900 mean = (*this)[ipx].GetMeanConversionPINDiodeMethod();
901 err = (*this)[ipx].GetConversionPINDiodeMethodErr();
902 sigma = (*this)[ipx].GetSigmaConversionPINDiodeMethod();
903
904 return kFALSE;
905
906}
907
908//-----------------------------------------------------------------------------------
909//
910// Calculates the best combination of the three used methods possible
911// between the integral of FADCs slices
912// (as defined in the signal extractor MExtractSignal.cc)
913// and the number of photons reaching one Inner Pixel.
914// The procedure is not yet defined.
915//
916// FIXME: The PINDiode is still not working and so is the code
917//
918Bool_t MCalibrationChargeCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
919{
920 return kFALSE;
921
922}
923
924/*
925void MCalibrationChargeCam::DrawHiLoFits()
926{
927
928 if (!fOffsets)
929 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
930 if (!fSlopes)
931 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
932 if (!fOffvsSlope)
933 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
934
935 TIter Next(fPixels);
936 MCalibrationPix *pix;
937 MHCalibrationPixel *hist;
938 while ((pix=(MCalibrationPix*)Next()))
939 {
940 hist = pix->GetHist();
941 hist->FitHiGainvsLoGain();
942 fOffsets->Fill(hist->GetOffset(),1.);
943 fSlopes->Fill(hist->GetSlope(),1.);
944 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
945 }
946
947 TCanvas *c1 = new TCanvas();
948
949 c1->Divide(1,3);
950 c1->cd(1);
951 fOffsets->Draw();
952 gPad->Modified();
953 gPad->Update();
954
955 c1->cd(2);
956 fSlopes->Draw();
957 gPad->Modified();
958 gPad->Update();
959
960 c1->cd(3);
961 fOffvsSlope->Draw("col1");
962 gPad->Modified();
963 gPad->Update();
964}
965
966*/
Note: See TracBrowser for help on using the repository browser.