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