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 MCalibrationPix 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 | // Fitted values:
|
---|
41 | // ==============
|
---|
42 | //
|
---|
43 | // 0: Fitted Charge
|
---|
44 | // 1: Error of fitted Charge
|
---|
45 | // 2: Sigma of fitted Charge
|
---|
46 | // 3: Error of Sigma of fitted Charge
|
---|
47 | //
|
---|
48 | // Useful variables derived from the fit results:
|
---|
49 | // =============================================
|
---|
50 | //
|
---|
51 | // 4: Returned probability of Gauss fit to Charge distribution
|
---|
52 | // 5: Reduced Sigma of fitted Charge --> sqrt(sigma_Q^2 - PedRMS^2)
|
---|
53 | // 6: Error Reduced Sigma of fitted Charge
|
---|
54 | // 7: Reduced Sigma per Charge
|
---|
55 | // 8: Error of Reduced Sigma per Charge
|
---|
56 | //
|
---|
57 | // Results of the different calibration methods:
|
---|
58 | // =============================================
|
---|
59 | //
|
---|
60 | // 9: Number of Photo-electrons obtained with the F-Factor method
|
---|
61 | // 10: Error on Number of Photo-electrons obtained with the F-Factor method
|
---|
62 | // 11: Mean conversion factor obtained with the F-Factor method
|
---|
63 | // 12: Error on the mean conversion factor obtained with the F-Factor method
|
---|
64 | // 13: Overall F-Factor of the readout obtained with the F-Factor method
|
---|
65 | // 14: Error on Overall F-Factor of the readout obtained with the F-Factor method
|
---|
66 | // 15: Number of Photons inside Plexiglass obtained with the Blind Pixel method
|
---|
67 | // 16: Error on Number of Photons inside Plexiglass obtained with the Blind Pixel method
|
---|
68 | // 17: Mean conversion factor obtained with the Blind Pixel method
|
---|
69 | // 18: Error on the mean conversion factor obtained with the Blind Pixel method
|
---|
70 | // 19: Overall F-Factor of the readout obtained with the Blind Pixel method
|
---|
71 | // 20: Error on Overall F-Factor of the readout obtained with the Blind Pixel method
|
---|
72 | // 21: Number of Photons outside Plexiglass obtained with the PIN Diode method
|
---|
73 | // 22: Error on Number of Photons outside Plexiglass obtained with the PIN Diode method
|
---|
74 | // 23: Mean conversion factor obtained with the PIN Diode method
|
---|
75 | // 24: Error on the mean conversion factor obtained with the PIN Diode method
|
---|
76 | // 25: Overall F-Factor of the readout obtained with the PIN Diode method
|
---|
77 | // 26: Error on Overall F-Factor of the readout obtained with the PIN Diode method
|
---|
78 | //
|
---|
79 | // Localized defects:
|
---|
80 | // ==================
|
---|
81 | //
|
---|
82 | // 27: Excluded Pixels
|
---|
83 | // 28: Pixels where the fit did not succeed --> results obtained only from the histograms
|
---|
84 | // 29: Pixels with succeeded fit, but apparently wrong results
|
---|
85 | // 30: Pixels with un-expected behavior in the fourier spectrum (e.g. oscillations)
|
---|
86 | //
|
---|
87 | // Other classifications of pixels:
|
---|
88 | // ================================
|
---|
89 | //
|
---|
90 | // 31: Pixels with saturated Hi-Gain
|
---|
91 | //
|
---|
92 | // Classification of validity of the calibrations:
|
---|
93 | // ===============================================
|
---|
94 | //
|
---|
95 | // 32: Pixels with valid calibration by the F-Factor-Method
|
---|
96 | // 33: Pixels with valid calibration by the Blind Pixel-Method
|
---|
97 | // 34: Pixels with valid calibration by the PIN Diode-Method
|
---|
98 | //
|
---|
99 | // Used Pedestals:
|
---|
100 | // ===============
|
---|
101 | //
|
---|
102 | // 35: Mean Pedestal over the entire range of signal extraction
|
---|
103 | // 36: Error on the Mean Pedestal over the entire range of signal extraction
|
---|
104 | // 37: Pedestal RMS over the entire range of signal extraction
|
---|
105 | // 38: Error on the Pedestal RMS over the entire range of signal extraction
|
---|
106 | //
|
---|
107 | // Calculated absolute arrival times (very low precision!):
|
---|
108 | // ========================================================
|
---|
109 | //
|
---|
110 | // 39: Absolute Arrival time of the signal
|
---|
111 | // 40: Error on the Absolute Arrival time of the signal
|
---|
112 | // 41: RMS of the Absolute Arrival time of the signal
|
---|
113 | // 42: Error on the RMS of the Absolute Arrival time of the signal
|
---|
114 | //
|
---|
115 | /////////////////////////////////////////////////////////////////////////////
|
---|
116 | #include "MCalibrationChargeCam.h"
|
---|
117 |
|
---|
118 | #include <TH2.h>
|
---|
119 | #include <TCanvas.h>
|
---|
120 | #include <TClonesArray.h>
|
---|
121 |
|
---|
122 | #include "MLog.h"
|
---|
123 | #include "MLogManip.h"
|
---|
124 |
|
---|
125 | #include "MGeomCam.h"
|
---|
126 | #include "MGeomPix.h"
|
---|
127 |
|
---|
128 | #include "MCalibrationChargePix.h"
|
---|
129 | #include "MCalibrationBlindPix.h"
|
---|
130 | #include "MCalibrationChargePINDiode.h"
|
---|
131 |
|
---|
132 |
|
---|
133 | ClassImp(MCalibrationChargeCam);
|
---|
134 |
|
---|
135 | using namespace std;
|
---|
136 |
|
---|
137 | const Int_t MCalibrationChargeCam::gkBlindPixelId = 559;
|
---|
138 | const Int_t MCalibrationChargeCam::gkPINDiodeId = 9999;
|
---|
139 | const Float_t MCalibrationChargeCam::gkBlindPixelArea = 100;
|
---|
140 | // Average QE of Blind Pixel (three colours)
|
---|
141 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEGreen = 0.154;
|
---|
142 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEBlue = 0.226;
|
---|
143 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEUV = 0.247;
|
---|
144 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQECT1 = 0.247;
|
---|
145 | // Average QE Error of Blind Pixel (three colours)
|
---|
146 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEGreenErr = 0.015;
|
---|
147 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEBlueErr = 0.02;
|
---|
148 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQEUVErr = 0.02;
|
---|
149 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelQECT1Err = 0.02;
|
---|
150 | // Attenuation factor Blind Pixel (three colours)
|
---|
151 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelAttGreen = 1.97;
|
---|
152 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelAttBlue = 1.96;
|
---|
153 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelAttUV = 1.95;
|
---|
154 | const Float_t MCalibrationChargeCam::gkCalibrationBlindPixelAttCT1 = 1.95;
|
---|
155 |
|
---|
156 | // --------------------------------------------------------------------------
|
---|
157 | //
|
---|
158 | // Default constructor.
|
---|
159 | //
|
---|
160 | // Creates a TClonesArray of MCalibrationPix containers, initialized to 1 entry
|
---|
161 | // Later, a call to MCalibrationChargeCam::InitSize(Int_t size) has to be performed
|
---|
162 | //
|
---|
163 | // Creates an MCalibrationBlindPix container
|
---|
164 | //
|
---|
165 | MCalibrationChargeCam::MCalibrationChargeCam(const char *name, const char *title)
|
---|
166 | : fBlindPixel(NULL),
|
---|
167 | fPINDiode(NULL),
|
---|
168 | fOffsets(NULL),
|
---|
169 | fSlopes(NULL),
|
---|
170 | fOffvsSlope(NULL)
|
---|
171 | {
|
---|
172 | fName = name ? name : "MCalibrationChargeCam";
|
---|
173 | fTitle = title ? title : "Storage container for the Calibration Information in the camera";
|
---|
174 |
|
---|
175 | fPixels = new TClonesArray("MCalibrationChargePix",1);
|
---|
176 | fBlindPixel = new MCalibrationBlindPix();
|
---|
177 |
|
---|
178 | Clear();
|
---|
179 | }
|
---|
180 |
|
---|
181 | // --------------------------------------------------------------------------
|
---|
182 | //
|
---|
183 | // Delete the TClonesArray of MCalibrationPix containers
|
---|
184 | // Delete the MCalibrationPINDiode and the MCalibrationBlindPix
|
---|
185 | //
|
---|
186 | // Delete the histograms if they exist
|
---|
187 | //
|
---|
188 | MCalibrationChargeCam::~MCalibrationChargeCam()
|
---|
189 | {
|
---|
190 |
|
---|
191 | //
|
---|
192 | // delete fPixels should delete all Objects stored inside
|
---|
193 | //
|
---|
194 | delete fPixels;
|
---|
195 | delete fBlindPixel;
|
---|
196 |
|
---|
197 | if (fOffsets)
|
---|
198 | delete fOffsets;
|
---|
199 | if (fSlopes)
|
---|
200 | delete fSlopes;
|
---|
201 | if (fOffvsSlope)
|
---|
202 | delete fOffvsSlope;
|
---|
203 |
|
---|
204 | }
|
---|
205 |
|
---|
206 | // -------------------------------------------------------------------
|
---|
207 | //
|
---|
208 | // This function simply allocates memory via the ROOT command:
|
---|
209 | // (TObject**) TStorage::ReAlloc(fCont, newSize * sizeof(TObject*),
|
---|
210 | // fSize * sizeof(TObject*));
|
---|
211 | // newSize corresponds to size in our case
|
---|
212 | // fSize is the old size (in most cases: 1)
|
---|
213 | //
|
---|
214 | void MCalibrationChargeCam::InitSize(const UInt_t i)
|
---|
215 | {
|
---|
216 |
|
---|
217 | //
|
---|
218 | // check if we have already initialized to size
|
---|
219 | //
|
---|
220 | if (CheckBounds(i))
|
---|
221 | return;
|
---|
222 |
|
---|
223 | fPixels->ExpandCreate(i);
|
---|
224 |
|
---|
225 | }
|
---|
226 |
|
---|
227 | // --------------------------------------------------------------------------
|
---|
228 | //
|
---|
229 | // This function returns the current size of the TClonesArray
|
---|
230 | // independently if the MCalibrationPix is filled with values or not.
|
---|
231 | //
|
---|
232 | // It is the size of the array fPixels.
|
---|
233 | //
|
---|
234 | Int_t MCalibrationChargeCam::GetSize() const
|
---|
235 | {
|
---|
236 | return fPixels->GetEntriesFast();
|
---|
237 | }
|
---|
238 |
|
---|
239 | // --------------------------------------------------------------------------
|
---|
240 | //
|
---|
241 | // Check if position i is inside the current bounds of the TClonesArray
|
---|
242 | //
|
---|
243 | Bool_t MCalibrationChargeCam::CheckBounds(Int_t i) const
|
---|
244 | {
|
---|
245 | return i < GetSize();
|
---|
246 | }
|
---|
247 |
|
---|
248 |
|
---|
249 | // --------------------------------------------------------------------------
|
---|
250 | //
|
---|
251 | // Get i-th pixel (pixel number)
|
---|
252 | //
|
---|
253 | MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i)
|
---|
254 | {
|
---|
255 | return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
|
---|
256 | }
|
---|
257 |
|
---|
258 | // --------------------------------------------------------------------------
|
---|
259 | //
|
---|
260 | // Get i-th pixel (pixel number)
|
---|
261 | //
|
---|
262 | const MCalibrationChargePix &MCalibrationChargeCam::operator[](UInt_t i) const
|
---|
263 | {
|
---|
264 | return *static_cast<MCalibrationChargePix*>(fPixels->UncheckedAt(i));
|
---|
265 | }
|
---|
266 |
|
---|
267 |
|
---|
268 | // --------------------------------------
|
---|
269 | //
|
---|
270 | void MCalibrationChargeCam::Clear(Option_t *o)
|
---|
271 | {
|
---|
272 |
|
---|
273 | fPixels->ForEach(TObject, Clear)();
|
---|
274 | fBlindPixel->Clear();
|
---|
275 |
|
---|
276 | fMeanFluxInsidePlexiglass = -1.;
|
---|
277 | fMeanFluxErrInsidePlexiglass = -1.;
|
---|
278 |
|
---|
279 | fNumExcludedPixels = 0;
|
---|
280 |
|
---|
281 | CLRBIT(fFlags,kBlindPixelMethodValid);
|
---|
282 | CLRBIT(fFlags,kPINDiodeMethodValid);
|
---|
283 | CLRBIT(fFlags,kFluxInsidePlexiglassAvailable);
|
---|
284 |
|
---|
285 | return;
|
---|
286 | }
|
---|
287 |
|
---|
288 | void MCalibrationChargeCam::SetBlindPixelMethodValid(const Bool_t b)
|
---|
289 | {
|
---|
290 | b ? SETBIT(fFlags, kBlindPixelMethodValid) : CLRBIT(fFlags, kBlindPixelMethodValid);
|
---|
291 | }
|
---|
292 |
|
---|
293 | void MCalibrationChargeCam::SetPINDiodeMethodValid(const Bool_t b)
|
---|
294 | {
|
---|
295 | b ? SETBIT(fFlags, kPINDiodeMethodValid) : CLRBIT(fFlags, kPINDiodeMethodValid);
|
---|
296 | }
|
---|
297 |
|
---|
298 | Bool_t MCalibrationChargeCam::IsBlindPixelMethodValid() const
|
---|
299 | {
|
---|
300 | return TESTBIT(fFlags,kBlindPixelMethodValid);
|
---|
301 | }
|
---|
302 |
|
---|
303 | Bool_t MCalibrationChargeCam::IsPINDiodeMethodValid() const
|
---|
304 | {
|
---|
305 | return TESTBIT(fFlags,kPINDiodeMethodValid);
|
---|
306 | }
|
---|
307 |
|
---|
308 |
|
---|
309 | Bool_t MCalibrationChargeCam::IsFluxInsidePlexiglassAvailable() const
|
---|
310 | {
|
---|
311 | return TESTBIT(fFlags,kFluxInsidePlexiglassAvailable);
|
---|
312 | }
|
---|
313 |
|
---|
314 |
|
---|
315 | // --------------------------------------------------------------------------
|
---|
316 | //
|
---|
317 | // Print first the well fitted pixels
|
---|
318 | // and then the ones which are not FitValid
|
---|
319 | //
|
---|
320 | void MCalibrationChargeCam::Print(Option_t *o) const
|
---|
321 | {
|
---|
322 |
|
---|
323 | *fLog << all << GetDescriptor() << ":" << endl;
|
---|
324 | int id = 0;
|
---|
325 |
|
---|
326 | *fLog << all << "Succesfully calibrated pixels:" << endl;
|
---|
327 | *fLog << all << endl;
|
---|
328 |
|
---|
329 | TIter Next(fPixels);
|
---|
330 | MCalibrationChargePix *pix;
|
---|
331 | while ((pix=(MCalibrationChargePix*)Next()))
|
---|
332 | {
|
---|
333 |
|
---|
334 | if (pix->IsChargeValid() && !pix->IsExcluded() && !pix->IsOscillating())
|
---|
335 | {
|
---|
336 |
|
---|
337 | *fLog << all << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- "
|
---|
338 | << pix->GetPedRms() << " Reduced Charge: " << pix->GetMeanCharge() << " +- "
|
---|
339 | << pix->GetSigmaCharge() << " Reduced Sigma: " << pix->GetRSigmaCharge()
|
---|
340 | << " Nr Phe's: " << pix->GetPheFFactorMethod() << endl;
|
---|
341 | id++;
|
---|
342 | }
|
---|
343 | }
|
---|
344 |
|
---|
345 | *fLog << all << id << " succesful pixels :-))" << endl;
|
---|
346 | id = 0;
|
---|
347 |
|
---|
348 | *fLog << all << endl;
|
---|
349 | *fLog << all << "Pixels with errors:" << endl;
|
---|
350 | *fLog << all << endl;
|
---|
351 |
|
---|
352 | TIter Next2(fPixels);
|
---|
353 | while ((pix=(MCalibrationChargePix*)Next2()))
|
---|
354 | {
|
---|
355 |
|
---|
356 | if (!pix->IsExcluded() && !pix->IsChargeValid())
|
---|
357 | {
|
---|
358 |
|
---|
359 | *fLog << all << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- "
|
---|
360 | << pix->GetPedRms() << " Reduced Charge: " << pix->GetMeanCharge() << " +- "
|
---|
361 | << pix->GetSigmaCharge() << " Reduced Sigma: " << pix->GetRSigmaCharge() << 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 ( pix->IsOscillating() && !pix->IsExcluded())
|
---|
378 | {
|
---|
379 |
|
---|
380 | *fLog << all << pix->GetPixId() << " Pedestals: " << pix->GetPed() << " +- "
|
---|
381 | << pix->GetPedRms() << " Reduced Charge: " << pix->GetMeanCharge() << " +- "
|
---|
382 | << pix->GetSigmaCharge() << " Reduced Sigma: " << pix->GetRSigmaCharge() << endl;
|
---|
383 | id++;
|
---|
384 | }
|
---|
385 | }
|
---|
386 | *fLog << all << id << " Oscillating pixels :-((" << endl;
|
---|
387 |
|
---|
388 |
|
---|
389 | *fLog << all << endl;
|
---|
390 | *fLog << all << "Excluded pixels:" << endl;
|
---|
391 | *fLog << all << endl;
|
---|
392 |
|
---|
393 | id = 0;
|
---|
394 |
|
---|
395 | TIter Next4(fPixels);
|
---|
396 | while ((pix=(MCalibrationChargePix*)Next4()))
|
---|
397 | {
|
---|
398 | if (pix->IsExcluded())
|
---|
399 | {
|
---|
400 | *fLog << all << pix->GetPixId() << endl;
|
---|
401 | id++;
|
---|
402 | }
|
---|
403 | }
|
---|
404 | *fLog << all << id << " Excluded pixels " << endl;
|
---|
405 | }
|
---|
406 |
|
---|
407 |
|
---|
408 | // --------------------------------------------------------------------------
|
---|
409 | //
|
---|
410 | // Sets the user ranges of all histograms such that
|
---|
411 | // empty bins at the edges are not used. Additionally, it rebins the
|
---|
412 | // histograms such that in total, 50 bins are used.
|
---|
413 | //
|
---|
414 | void MCalibrationChargeCam::CutEdges()
|
---|
415 | {
|
---|
416 |
|
---|
417 | fBlindPixel->GetHist()->CutAllEdges();
|
---|
418 |
|
---|
419 | return;
|
---|
420 | }
|
---|
421 |
|
---|
422 | // --------------------------------------------------------------------------
|
---|
423 | //
|
---|
424 | // The types are as follows:
|
---|
425 | //
|
---|
426 | // Fitted values:
|
---|
427 | // ==============
|
---|
428 | //
|
---|
429 | // 0: Fitted Charge
|
---|
430 | // 1: Error of fitted Charge
|
---|
431 | // 2: Sigma of fitted Charge
|
---|
432 | // 3: Error of Sigma of fitted Charge
|
---|
433 | //
|
---|
434 | // Useful variables derived from the fit results:
|
---|
435 | // =============================================
|
---|
436 | //
|
---|
437 | // 4: Returned probability of Gauss fit to Charge distribution
|
---|
438 | // 5: Reduced Sigma of fitted Charge --> sqrt(sigma_Q^2 - PedRMS^2)
|
---|
439 | // 6: Error Reduced Sigma of fitted Charge
|
---|
440 | // 7: Reduced Sigma per Charge
|
---|
441 | // 8: Error of Reduced Sigma per Charge
|
---|
442 | //
|
---|
443 | // Results of the different calibration methods:
|
---|
444 | // =============================================
|
---|
445 | //
|
---|
446 | // 9: Number of Photo-electrons obtained with the F-Factor method
|
---|
447 | // 10: Error on Number of Photo-electrons obtained with the F-Factor method
|
---|
448 | // 11: Mean conversion factor obtained with the F-Factor method
|
---|
449 | // 12: Error on the mean conversion factor obtained with the F-Factor method
|
---|
450 | // 13: Overall F-Factor of the readout obtained with the F-Factor method
|
---|
451 | // 14: Error on Overall F-Factor of the readout obtained with the F-Factor method
|
---|
452 | // 15: Number of Photons inside Plexiglass obtained with the Blind Pixel method
|
---|
453 | // 16: Error on Number of Photons inside Plexiglass obtained with the Blind Pixel method
|
---|
454 | // 17: Mean conversion factor obtained with the Blind Pixel method
|
---|
455 | // 18: Error on the mean conversion factor obtained with the Blind Pixel method
|
---|
456 | // 19: Overall F-Factor of the readout obtained with the Blind Pixel method
|
---|
457 | // 20: Error on Overall F-Factor of the readout obtained with the Blind Pixel method
|
---|
458 | // 21: Number of Photons outside Plexiglass obtained with the PIN Diode method
|
---|
459 | // 22: Error on Number of Photons outside Plexiglass obtained with the PIN Diode method
|
---|
460 | // 23: Mean conversion factor obtained with the PIN Diode method
|
---|
461 | // 24: Error on the mean conversion factor obtained with the PIN Diode method
|
---|
462 | // 25: Overall F-Factor of the readout obtained with the PIN Diode method
|
---|
463 | // 26: Error on Overall F-Factor of the readout obtained with the PIN Diode method
|
---|
464 | //
|
---|
465 | // Localized defects:
|
---|
466 | // ==================
|
---|
467 | //
|
---|
468 | // 27: Excluded Pixels
|
---|
469 | // 28: Pixels where the fit did not succeed --> results obtained only from the histograms
|
---|
470 | // 29: Pixels with succeeded fit, but apparently wrong results
|
---|
471 | // 30: Pixels with un-expected behavior in the Hi Gain fourier spectrum (e.g. oscillations)
|
---|
472 | // 31: Pixels with un-expected behavior in the Lo Gain fourier spectrum (e.g. oscillations)a
|
---|
473 | // 32: Number of probable pickup events in the Hi Gain
|
---|
474 | // 33: Number of probable pickup events in the Lo Gain
|
---|
475 | //
|
---|
476 | // Other classifications of pixels:
|
---|
477 | // ================================
|
---|
478 | //
|
---|
479 | // 34: Pixels with saturated Hi-Gain
|
---|
480 | //
|
---|
481 | // Classification of validity of the calibrations:
|
---|
482 | // ===============================================
|
---|
483 | //
|
---|
484 | // 35: Pixels with valid calibration by the F-Factor-Method
|
---|
485 | // 36: Pixels with valid calibration by the Blind Pixel-Method
|
---|
486 | // 37: Pixels with valid calibration by the PIN Diode-Method
|
---|
487 | //
|
---|
488 | // Used Pedestals:
|
---|
489 | // ===============
|
---|
490 | //
|
---|
491 | // 38: Mean Pedestal over the entire range of signal extraction
|
---|
492 | // 39: Error on the Mean Pedestal over the entire range of signal extraction
|
---|
493 | // 40: Pedestal RMS over the entire range of signal extraction
|
---|
494 | // 41: Error on the Pedestal RMS over the entire range of signal extraction
|
---|
495 | //
|
---|
496 | // Calculated absolute arrival times (very low precision!):
|
---|
497 | // ========================================================
|
---|
498 | //
|
---|
499 | // 42: Absolute Arrival time of the signal
|
---|
500 | // 43: RMS of the Absolute Arrival time of the signal
|
---|
501 | //
|
---|
502 | Bool_t MCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
|
---|
503 | {
|
---|
504 |
|
---|
505 | if (idx > GetSize())
|
---|
506 | return kFALSE;
|
---|
507 |
|
---|
508 | Float_t area = cam[idx].GetA();
|
---|
509 |
|
---|
510 | if (area == 0)
|
---|
511 | return kFALSE;
|
---|
512 |
|
---|
513 | switch (type)
|
---|
514 | {
|
---|
515 | case 0:
|
---|
516 | if ((*this)[idx].IsExcluded())
|
---|
517 | return kFALSE;
|
---|
518 | val = (*this)[idx].GetMeanCharge();
|
---|
519 | break;
|
---|
520 | case 1:
|
---|
521 | if ((*this)[idx].IsExcluded())
|
---|
522 | return kFALSE;
|
---|
523 | val = (*this)[idx].GetMeanChargeErr();
|
---|
524 | break;
|
---|
525 | case 2:
|
---|
526 | if ((*this)[idx].IsExcluded())
|
---|
527 | return kFALSE;
|
---|
528 | val = (*this)[idx].GetSigmaCharge();
|
---|
529 | break;
|
---|
530 | case 3:
|
---|
531 | if ((*this)[idx].IsExcluded())
|
---|
532 | return kFALSE;
|
---|
533 | val = (*this)[idx].GetSigmaChargeErr();
|
---|
534 | break;
|
---|
535 | case 4:
|
---|
536 | if ((*this)[idx].IsExcluded())
|
---|
537 | return kFALSE;
|
---|
538 | val = (*this)[idx].GetChargeProb();
|
---|
539 | break;
|
---|
540 | case 5:
|
---|
541 | if ((*this)[idx].IsExcluded())
|
---|
542 | return kFALSE;
|
---|
543 | val = (*this)[idx].GetRSigmaCharge();
|
---|
544 | break;
|
---|
545 | case 6:
|
---|
546 | if ((*this)[idx].IsExcluded())
|
---|
547 | return kFALSE;
|
---|
548 | val = (*this)[idx].GetRSigmaChargeErr();
|
---|
549 | break;
|
---|
550 | case 7:
|
---|
551 | if ((*this)[idx].IsExcluded())
|
---|
552 | return kFALSE;
|
---|
553 | val = (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
|
---|
554 | break;
|
---|
555 | case 8:
|
---|
556 | if ((*this)[idx].IsExcluded())
|
---|
557 | return kFALSE;
|
---|
558 | // relative error RsigmaCharge square
|
---|
559 | val = (*this)[idx].GetRSigmaChargeErr()* (*this)[idx].GetRSigmaChargeErr()
|
---|
560 | / ((*this)[idx].GetRSigmaCharge() * (*this)[idx].GetRSigmaCharge() );
|
---|
561 | // relative error Charge square
|
---|
562 | val += (*this)[idx].GetMeanChargeErr() * (*this)[idx].GetMeanChargeErr()
|
---|
563 | / ((*this)[idx].GetMeanCharge() * (*this)[idx].GetMeanCharge() );
|
---|
564 | // calculate relative error out of squares
|
---|
565 | val = TMath::Sqrt(val) ;
|
---|
566 | // multiply with value to get absolute error
|
---|
567 | val *= (*this)[idx].GetRSigmaCharge() / (*this)[idx].GetMeanCharge();
|
---|
568 | break;
|
---|
569 | case 9:
|
---|
570 | if ((*this)[idx].IsExcluded())
|
---|
571 | return kFALSE;
|
---|
572 | val = (*this)[idx].GetPheFFactorMethod();
|
---|
573 | break;
|
---|
574 | case 10:
|
---|
575 | if ((*this)[idx].IsExcluded())
|
---|
576 | return kFALSE;
|
---|
577 | val = (*this)[idx].GetPheFFactorMethodErr();
|
---|
578 | break;
|
---|
579 | case 11:
|
---|
580 | if ((*this)[idx].IsExcluded())
|
---|
581 | return kFALSE;
|
---|
582 | val = (*this)[idx].GetMeanConversionFFactorMethod();
|
---|
583 | break;
|
---|
584 | case 12:
|
---|
585 | if ((*this)[idx].IsExcluded())
|
---|
586 | return kFALSE;
|
---|
587 | val = (*this)[idx].GetConversionFFactorMethodErr();
|
---|
588 | break;
|
---|
589 | case 13:
|
---|
590 | if ((*this)[idx].IsExcluded())
|
---|
591 | return kFALSE;
|
---|
592 | val = (*this)[idx].GetTotalFFactorFFactorMethod();
|
---|
593 | break;
|
---|
594 | case 14:
|
---|
595 | if ((*this)[idx].IsExcluded())
|
---|
596 | return kFALSE;
|
---|
597 | val = (*this)[idx].GetTotalFFactorErrFFactorMethod();
|
---|
598 | break;
|
---|
599 | case 15:
|
---|
600 | if ((*this)[idx].IsExcluded())
|
---|
601 | return kFALSE;
|
---|
602 | val = GetMeanFluxInsidePlexiglass()*area;
|
---|
603 | break;
|
---|
604 | case 16:
|
---|
605 | if ((*this)[idx].IsExcluded())
|
---|
606 | return kFALSE;
|
---|
607 | val = GetMeanFluxInsidePlexiglass()*area;
|
---|
608 | break;
|
---|
609 | case 17:
|
---|
610 | if ((*this)[idx].IsExcluded())
|
---|
611 | return kFALSE;
|
---|
612 | val = (*this)[idx].GetMeanConversionBlindPixelMethod();
|
---|
613 | break;
|
---|
614 | case 18:
|
---|
615 | if ((*this)[idx].IsExcluded())
|
---|
616 | return kFALSE;
|
---|
617 | val = (*this)[idx].GetConversionBlindPixelMethodErr();
|
---|
618 | break;
|
---|
619 | case 19:
|
---|
620 | if ((*this)[idx].IsExcluded())
|
---|
621 | return kFALSE;
|
---|
622 | val = (*this)[idx].GetTotalFFactorBlindPixelMethod();
|
---|
623 | break;
|
---|
624 | case 20:
|
---|
625 | if ((*this)[idx].IsExcluded())
|
---|
626 | return kFALSE;
|
---|
627 | val = (*this)[idx].GetTotalFFactorErrBlindPixelMethod();
|
---|
628 | break;
|
---|
629 | case 21:
|
---|
630 | if ((*this)[idx].IsExcluded())
|
---|
631 | return kFALSE;
|
---|
632 | val = fPINDiode->GetMeanFluxOutsidePlexiglass()*area;
|
---|
633 | break;
|
---|
634 | case 22:
|
---|
635 | if ((*this)[idx].IsExcluded())
|
---|
636 | return kFALSE;
|
---|
637 | val = fPINDiode->GetMeanFluxOutsidePlexiglass()*area;
|
---|
638 | break;
|
---|
639 | case 23:
|
---|
640 | if ((*this)[idx].IsExcluded())
|
---|
641 | return kFALSE;
|
---|
642 | val = (*this)[idx].GetMeanConversionPINDiodeMethod();
|
---|
643 | break;
|
---|
644 | case 24:
|
---|
645 | if ((*this)[idx].IsExcluded())
|
---|
646 | return kFALSE;
|
---|
647 | val = (*this)[idx].GetConversionPINDiodeMethodErr();
|
---|
648 | break;
|
---|
649 | case 25:
|
---|
650 | if ((*this)[idx].IsExcluded())
|
---|
651 | return kFALSE;
|
---|
652 | val = (*this)[idx].GetTotalFFactorBlindPixelMethod();
|
---|
653 | break;
|
---|
654 | case 26:
|
---|
655 | if ((*this)[idx].IsExcluded())
|
---|
656 | return kFALSE;
|
---|
657 | val = (*this)[idx].GetTotalFFactorErrBlindPixelMethod();
|
---|
658 | break;
|
---|
659 | case 27:
|
---|
660 | if ((*this)[idx].IsExcluded())
|
---|
661 | val = 1.;
|
---|
662 | else
|
---|
663 | return kFALSE;
|
---|
664 | break;
|
---|
665 | case 28:
|
---|
666 | if ((*this)[idx].IsExcluded())
|
---|
667 | return kFALSE;
|
---|
668 | if (!(*this)[idx].IsFitted())
|
---|
669 | val = 1;
|
---|
670 | else
|
---|
671 | return kFALSE;
|
---|
672 | break;
|
---|
673 | case 29:
|
---|
674 | if ((*this)[idx].IsExcluded())
|
---|
675 | return kFALSE;
|
---|
676 | if (!(*this)[idx].IsFitted())
|
---|
677 | return kFALSE;
|
---|
678 | if (!(*this)[idx].IsChargeValid())
|
---|
679 | val = 1;
|
---|
680 | else
|
---|
681 | return kFALSE;
|
---|
682 | break;
|
---|
683 | case 30:
|
---|
684 | if ((*this)[idx].IsExcluded())
|
---|
685 | return kFALSE;
|
---|
686 | if ((*this)[idx].IsHiGainOscillating())
|
---|
687 | val = 1;
|
---|
688 | else
|
---|
689 | return kFALSE;
|
---|
690 | break;
|
---|
691 | case 31:
|
---|
692 | if ((*this)[idx].IsExcluded())
|
---|
693 | return kFALSE;
|
---|
694 | if ((*this)[idx].IsLoGainOscillating())
|
---|
695 | val = 1;
|
---|
696 | else
|
---|
697 | return kFALSE;
|
---|
698 | break;
|
---|
699 | case 32:
|
---|
700 | if ((*this)[idx].IsExcluded())
|
---|
701 | return kFALSE;
|
---|
702 | val = (*this)[idx].GetHiGainNumPickup();
|
---|
703 | break;
|
---|
704 | case 33:
|
---|
705 | if ((*this)[idx].IsExcluded())
|
---|
706 | return kFALSE;
|
---|
707 | val = (*this)[idx].GetLoGainNumPickup();
|
---|
708 | break;
|
---|
709 | case 34:
|
---|
710 | if ((*this)[idx].IsExcluded())
|
---|
711 | return kFALSE;
|
---|
712 | if ((*this)[idx].IsHiGainSaturation())
|
---|
713 | val = 1;
|
---|
714 | else
|
---|
715 | return kFALSE;
|
---|
716 | break;
|
---|
717 | case 35:
|
---|
718 | if ((*this)[idx].IsExcluded())
|
---|
719 | return kFALSE;
|
---|
720 | if ((*this)[idx].IsFFactorMethodValid())
|
---|
721 | val = 1;
|
---|
722 | else
|
---|
723 | return kFALSE;
|
---|
724 | break;
|
---|
725 | case 36:
|
---|
726 | if ((*this)[idx].IsExcluded())
|
---|
727 | return kFALSE;
|
---|
728 | if ((*this)[idx].IsBlindPixelMethodValid())
|
---|
729 | val = 1;
|
---|
730 | else
|
---|
731 | return kFALSE;
|
---|
732 | break;
|
---|
733 | case 37:
|
---|
734 | if ((*this)[idx].IsExcluded())
|
---|
735 | return kFALSE;
|
---|
736 | if ((*this)[idx].IsPINDiodeMethodValid())
|
---|
737 | val = 1;
|
---|
738 | else
|
---|
739 | return kFALSE;
|
---|
740 | break;
|
---|
741 | case 38:
|
---|
742 | if ((*this)[idx].IsExcluded())
|
---|
743 | return kFALSE;
|
---|
744 | val = (*this)[idx].GetPed();
|
---|
745 | break;
|
---|
746 | case 39:
|
---|
747 | if ((*this)[idx].IsExcluded())
|
---|
748 | return kFALSE;
|
---|
749 | val = (*this)[idx].GetPedErr();
|
---|
750 | break;
|
---|
751 | case 40:
|
---|
752 | if ((*this)[idx].IsExcluded())
|
---|
753 | return kFALSE;
|
---|
754 | val = (*this)[idx].GetPedRms();
|
---|
755 | break;
|
---|
756 | case 41:
|
---|
757 | if ((*this)[idx].IsExcluded())
|
---|
758 | return kFALSE;
|
---|
759 | val = (*this)[idx].GetPedErr()/2.;
|
---|
760 | break;
|
---|
761 | case 42:
|
---|
762 | if ((*this)[idx].IsExcluded())
|
---|
763 | return kFALSE;
|
---|
764 | val = (*this)[idx].GetAbsTimeMean();
|
---|
765 | break;
|
---|
766 | case 43:
|
---|
767 | if ((*this)[idx].IsExcluded())
|
---|
768 | return kFALSE;
|
---|
769 | val = (*this)[idx].GetAbsTimeRms();
|
---|
770 | break;
|
---|
771 | default:
|
---|
772 | return kFALSE;
|
---|
773 | }
|
---|
774 | return val!=-1.;
|
---|
775 | }
|
---|
776 |
|
---|
777 | // --------------------------------------------------------------------------
|
---|
778 | //
|
---|
779 | // What MHCamera needs in order to draw an individual pixel in the camera
|
---|
780 | //
|
---|
781 | void MCalibrationChargeCam::DrawPixelContent(Int_t idx) const
|
---|
782 | {
|
---|
783 | (*this)[idx].DrawClone();
|
---|
784 | }
|
---|
785 |
|
---|
786 |
|
---|
787 | // --------------------------------------------------------------------------
|
---|
788 | //
|
---|
789 | //
|
---|
790 | //
|
---|
791 | Bool_t MCalibrationChargeCam::CalcFluxInsidePlexiglass()
|
---|
792 | {
|
---|
793 |
|
---|
794 | if (!fBlindPixel->IsFitOK())
|
---|
795 | return kFALSE;
|
---|
796 |
|
---|
797 | const Float_t mean = fBlindPixel->GetLambda();
|
---|
798 | const Float_t merr = fBlindPixel->GetErrLambda();
|
---|
799 |
|
---|
800 | //
|
---|
801 | // Start calculation of number of photons
|
---|
802 | //
|
---|
803 | // The blind pixel has exactly 100 mm^2 area (with negligible error),
|
---|
804 | //
|
---|
805 | fMeanFluxInsidePlexiglass = mean*gkBlindPixelArea;
|
---|
806 |
|
---|
807 | // Start calculation of number of photons relative Variance (!!)
|
---|
808 | fMeanFluxErrInsidePlexiglass = merr*merr/mean/mean;
|
---|
809 |
|
---|
810 | switch (fColor)
|
---|
811 | {
|
---|
812 | case kEGreen:
|
---|
813 | fMeanFluxInsidePlexiglass /= gkCalibrationBlindPixelQEGreen;
|
---|
814 | fMeanFluxErrInsidePlexiglass += gkCalibrationBlindPixelQEGreenErr*gkCalibrationBlindPixelQEGreenErr
|
---|
815 | / gkCalibrationBlindPixelQEGreen / gkCalibrationBlindPixelQEGreen;
|
---|
816 |
|
---|
817 | fMeanFluxInsidePlexiglass *= TMath::Power(10,gkCalibrationBlindPixelAttGreen); // correct for absorption
|
---|
818 | // attenuation has negligible error
|
---|
819 | break;
|
---|
820 | case kEBlue:
|
---|
821 | fMeanFluxInsidePlexiglass /= gkCalibrationBlindPixelQEBlue;
|
---|
822 | fMeanFluxErrInsidePlexiglass += gkCalibrationBlindPixelQEBlueErr*gkCalibrationBlindPixelQEBlueErr
|
---|
823 | / gkCalibrationBlindPixelQEBlue / gkCalibrationBlindPixelQEBlue;
|
---|
824 |
|
---|
825 | fMeanFluxInsidePlexiglass *= TMath::Power(10,gkCalibrationBlindPixelAttBlue); // correct for absorption
|
---|
826 | // attenuation has negligible error
|
---|
827 | break;
|
---|
828 | case kEUV:
|
---|
829 | fMeanFluxInsidePlexiglass /= gkCalibrationBlindPixelQEUV;
|
---|
830 | fMeanFluxErrInsidePlexiglass += gkCalibrationBlindPixelQEUVErr*gkCalibrationBlindPixelQEUVErr
|
---|
831 | / gkCalibrationBlindPixelQEUV / gkCalibrationBlindPixelQEUV;
|
---|
832 |
|
---|
833 | fMeanFluxInsidePlexiglass *= TMath::Power(10,gkCalibrationBlindPixelAttUV); // correct for absorption
|
---|
834 | // attenuation has negligible error
|
---|
835 | break;
|
---|
836 | case kECT1:
|
---|
837 | default:
|
---|
838 | fMeanFluxInsidePlexiglass /= gkCalibrationBlindPixelQECT1;
|
---|
839 | fMeanFluxErrInsidePlexiglass += gkCalibrationBlindPixelQECT1Err*gkCalibrationBlindPixelQECT1Err
|
---|
840 | / gkCalibrationBlindPixelQECT1 / gkCalibrationBlindPixelQECT1;
|
---|
841 |
|
---|
842 | fMeanFluxInsidePlexiglass *= TMath::Power(10,gkCalibrationBlindPixelAttCT1); // correct for absorption
|
---|
843 | // attenuation has negligible error
|
---|
844 | break;
|
---|
845 | }
|
---|
846 |
|
---|
847 | *fLog << inf << endl;
|
---|
848 | *fLog << inf << " Photon flux [ph/mm^2] inside Plexiglass: "
|
---|
849 | << fMeanFluxInsidePlexiglass << endl;
|
---|
850 |
|
---|
851 | if (fMeanFluxInsidePlexiglass > 0.)
|
---|
852 | SETBIT(fFlags,kFluxInsidePlexiglassAvailable);
|
---|
853 | else
|
---|
854 | {
|
---|
855 | CLRBIT(fFlags,kFluxInsidePlexiglassAvailable);
|
---|
856 | return kFALSE;
|
---|
857 | }
|
---|
858 |
|
---|
859 | if (fMeanFluxErrInsidePlexiglass < 0.)
|
---|
860 | {
|
---|
861 | *fLog << warn << " Relative Variance on Photon flux inside Plexiglass: "
|
---|
862 | << fMeanFluxErrInsidePlexiglass << endl;
|
---|
863 | CLRBIT(fFlags,kFluxInsidePlexiglassAvailable);
|
---|
864 | return kFALSE;
|
---|
865 | }
|
---|
866 |
|
---|
867 | // Finish calculation of errors -> convert from relative variance to absolute error
|
---|
868 | fMeanFluxErrInsidePlexiglass = TMath::Sqrt(fMeanFluxErrInsidePlexiglass);
|
---|
869 | fMeanFluxErrInsidePlexiglass *= fMeanFluxInsidePlexiglass;
|
---|
870 |
|
---|
871 | *fLog << inf << " Error on photon flux [ph/mm^2] inside Plexiglass: "
|
---|
872 | << fMeanFluxErrInsidePlexiglass << endl;
|
---|
873 | *fLog << inf << endl;
|
---|
874 |
|
---|
875 | TIter Next(fPixels);
|
---|
876 | MCalibrationChargePix *pix;
|
---|
877 | while ((pix=(MCalibrationChargePix*)Next()))
|
---|
878 | {
|
---|
879 |
|
---|
880 | if(pix->IsChargeValid())
|
---|
881 | {
|
---|
882 |
|
---|
883 | const Int_t idx = pix->GetPixId();
|
---|
884 |
|
---|
885 | const Float_t charge = pix->GetMeanCharge();
|
---|
886 | const Float_t area = (*fGeomCam)[idx].GetA();
|
---|
887 | const Float_t chargeerr = pix->GetMeanChargeErr();
|
---|
888 |
|
---|
889 | const Float_t nphot = fMeanFluxInsidePlexiglass*area;
|
---|
890 | const Float_t nphoterr = fMeanFluxErrInsidePlexiglass*area;
|
---|
891 | const Float_t conversion = nphot/charge;
|
---|
892 | Float_t conversionerr;
|
---|
893 |
|
---|
894 | conversionerr = nphoterr/charge
|
---|
895 | * nphoterr/charge ;
|
---|
896 | conversionerr += chargeerr/charge
|
---|
897 | * chargeerr/charge
|
---|
898 | * conversion*conversion;
|
---|
899 | conversionerr = TMath::Sqrt(conversionerr);
|
---|
900 |
|
---|
901 | const Float_t conversionsigma = 0.;
|
---|
902 |
|
---|
903 | pix->SetConversionBlindPixelMethod(conversion, conversionerr, conversionsigma);
|
---|
904 |
|
---|
905 | if (conversionerr/conversion < 0.1)
|
---|
906 | pix->SetBlindPixelMethodValid();
|
---|
907 | }
|
---|
908 | }
|
---|
909 | return kTRUE;
|
---|
910 | }
|
---|
911 |
|
---|
912 |
|
---|
913 | void MCalibrationChargeCam::ApplyPINDiodeCalibration()
|
---|
914 | {
|
---|
915 |
|
---|
916 | Float_t flux = fPINDiode->GetMeanFluxOutsidePlexiglass();
|
---|
917 | Float_t fluxerr = fPINDiode->GetMeanFluxErrOutsidePlexiglass();
|
---|
918 |
|
---|
919 | TIter Next(fPixels);
|
---|
920 | MCalibrationChargePix *pix;
|
---|
921 | while ((pix=(MCalibrationChargePix*)Next()))
|
---|
922 | {
|
---|
923 |
|
---|
924 | if (pix->IsChargeValid())
|
---|
925 | {
|
---|
926 |
|
---|
927 | const Int_t idx = pix->GetPixId();
|
---|
928 |
|
---|
929 | const Float_t charge = pix->GetMeanCharge();
|
---|
930 | const Float_t area = (*fGeomCam)[idx].GetA();
|
---|
931 | const Float_t chargeerr = pix->GetMeanChargeErr();
|
---|
932 |
|
---|
933 | const Float_t nphot = flux * area;
|
---|
934 | const Float_t nphoterr = fluxerr * area;
|
---|
935 | const Float_t conversion = nphot/charge;
|
---|
936 |
|
---|
937 | Float_t conversionerr;
|
---|
938 |
|
---|
939 | conversionerr = nphoterr/charge
|
---|
940 | * nphoterr/charge ;
|
---|
941 | conversionerr += chargeerr/charge
|
---|
942 | * chargeerr/charge
|
---|
943 | * conversion*conversion;
|
---|
944 | if (conversionerr > 0.)
|
---|
945 | conversionerr = TMath::Sqrt(conversionerr);
|
---|
946 |
|
---|
947 | const Float_t conversionsigma = 0.;
|
---|
948 |
|
---|
949 | pix->SetConversionPINDiodeMethod(conversion, conversionerr, conversionsigma);
|
---|
950 |
|
---|
951 | if (conversionerr/conversion < 0.1)
|
---|
952 | pix->SetPINDiodeMethodValid();
|
---|
953 |
|
---|
954 | }
|
---|
955 | }
|
---|
956 | }
|
---|
957 |
|
---|
958 |
|
---|
959 |
|
---|
960 | Bool_t MCalibrationChargeCam::GetConversionFactorBlindPixel(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
|
---|
961 | {
|
---|
962 |
|
---|
963 |
|
---|
964 | if (!IsFluxInsidePlexiglassAvailable())
|
---|
965 | if (!CalcFluxInsidePlexiglass())
|
---|
966 | return kFALSE;
|
---|
967 |
|
---|
968 | mean = (*this)[ipx].GetMeanConversionBlindPixelMethod();
|
---|
969 | err = (*this)[ipx].GetConversionBlindPixelMethodErr();
|
---|
970 | sigma = (*this)[ipx].GetSigmaConversionBlindPixelMethod();
|
---|
971 |
|
---|
972 | return kTRUE;
|
---|
973 | }
|
---|
974 |
|
---|
975 |
|
---|
976 | Bool_t MCalibrationChargeCam::GetConversionFactorFFactor(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
|
---|
977 | {
|
---|
978 |
|
---|
979 | Float_t conv = (*this)[ipx].GetMeanConversionFFactorMethod();
|
---|
980 |
|
---|
981 | if (conv < 0.)
|
---|
982 | return kFALSE;
|
---|
983 |
|
---|
984 | mean = conv;
|
---|
985 | err = (*this)[ipx].GetConversionFFactorMethodErr();
|
---|
986 | sigma = (*this)[ipx].GetSigmaConversionFFactorMethod();
|
---|
987 |
|
---|
988 | return kTRUE;
|
---|
989 | }
|
---|
990 |
|
---|
991 |
|
---|
992 | //-----------------------------------------------------------------------------------
|
---|
993 | //
|
---|
994 | // Calculates the conversion factor between the integral of FADCs slices
|
---|
995 | // (as defined in the signal extractor MExtractSignal.cc)
|
---|
996 | // and the number of photons reaching the plexiglass for one Inner Pixel
|
---|
997 | //
|
---|
998 | // FIXME: The PINDiode is still not working and so is the code
|
---|
999 | //
|
---|
1000 | Bool_t MCalibrationChargeCam::GetConversionFactorPINDiode(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
|
---|
1001 | {
|
---|
1002 |
|
---|
1003 | mean = (*this)[ipx].GetMeanConversionPINDiodeMethod();
|
---|
1004 | err = (*this)[ipx].GetConversionPINDiodeMethodErr();
|
---|
1005 | sigma = (*this)[ipx].GetSigmaConversionPINDiodeMethod();
|
---|
1006 |
|
---|
1007 | return kFALSE;
|
---|
1008 |
|
---|
1009 | }
|
---|
1010 |
|
---|
1011 | //-----------------------------------------------------------------------------------
|
---|
1012 | //
|
---|
1013 | // Calculates the best combination of the three used methods possible
|
---|
1014 | // between the integral of FADCs slices
|
---|
1015 | // (as defined in the signal extractor MExtractSignal.cc)
|
---|
1016 | // and the number of photons reaching one Inner Pixel.
|
---|
1017 | // The procedure is not yet defined.
|
---|
1018 | //
|
---|
1019 | // FIXME: The PINDiode is still not working and so is the code
|
---|
1020 | //
|
---|
1021 | Bool_t MCalibrationChargeCam::GetConversionFactorCombined(Int_t ipx, Float_t &mean, Float_t &err, Float_t &sigma)
|
---|
1022 | {
|
---|
1023 | return kFALSE;
|
---|
1024 |
|
---|
1025 | }
|
---|
1026 |
|
---|
1027 | /*
|
---|
1028 | void MCalibrationChargeCam::DrawHiLoFits()
|
---|
1029 | {
|
---|
1030 |
|
---|
1031 | if (!fOffsets)
|
---|
1032 | fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
|
---|
1033 | if (!fSlopes)
|
---|
1034 | fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
|
---|
1035 | if (!fOffvsSlope)
|
---|
1036 | fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
|
---|
1037 |
|
---|
1038 | TIter Next(fPixels);
|
---|
1039 | MCalibrationPix *pix;
|
---|
1040 | MHCalibrationPixel *hist;
|
---|
1041 | while ((pix=(MCalibrationPix*)Next()))
|
---|
1042 | {
|
---|
1043 | hist = pix->GetHist();
|
---|
1044 | hist->FitHiGainvsLoGain();
|
---|
1045 | fOffsets->Fill(hist->GetOffset(),1.);
|
---|
1046 | fSlopes->Fill(hist->GetSlope(),1.);
|
---|
1047 | fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
|
---|
1048 | }
|
---|
1049 |
|
---|
1050 | TCanvas *c1 = new TCanvas();
|
---|
1051 |
|
---|
1052 | c1->Divide(1,3);
|
---|
1053 | c1->cd(1);
|
---|
1054 | fOffsets->Draw();
|
---|
1055 | gPad->Modified();
|
---|
1056 | gPad->Update();
|
---|
1057 |
|
---|
1058 | c1->cd(2);
|
---|
1059 | fSlopes->Draw();
|
---|
1060 | gPad->Modified();
|
---|
1061 | gPad->Update();
|
---|
1062 |
|
---|
1063 | c1->cd(3);
|
---|
1064 | fOffvsSlope->Draw("col1");
|
---|
1065 | gPad->Modified();
|
---|
1066 | gPad->Update();
|
---|
1067 | }
|
---|
1068 |
|
---|
1069 | */
|
---|