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

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