source: tags/Mars-V0.8.6/mcalib/MCalibrationChargeCalc.cc

Last change on this file was 5156, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 74.9 KB
Line 
1
2/* ======================================================================== *\
3!
4! *
5! * This file is part of MARS, the MAGIC Analysis and Reconstruction
6! * Software. It is distributed to you in the hope that it can be a useful
7! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
8! * It is distributed WITHOUT ANY WARRANTY.
9! *
10! * Permission to use, copy, modify and distribute this software and its
11! * documentation for any purpose is hereby granted without fee,
12! * provided that the above copyright notice appear in all copies and
13! * that both that copyright notice and this permission notice appear
14! * in supporting documentation. It is provided "as is" without express
15! * or implied warranty.
16! *
17!
18!
19! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25//////////////////////////////////////////////////////////////////////////////
26//
27// MCalibrationChargeCalc
28//
29// Task to calculate the calibration conversion factors and quantum efficiencies
30// from the fit results to the summed FADC slice distributions delivered by
31// MCalibrationChargeCam, MCalibrationChargePix, MCalibrationChargeBlindPix and
32// MCalibrationChargePINDiode, calculated and filled by MHCalibrationChargeCam,
33// MHCalibrationChargePix, MHCalibrationChargeBlindPix and MHCalibrationChargePINDiode.
34//
35// PreProcess(): Initialize pointers to MCalibrationChargeCam, MCalibrationChargeBlindPix
36// MCalibrationChargePINDiode and MCalibrationQECam
37//
38// Initialize pulser light wavelength
39//
40// ReInit(): MCalibrationCam::InitSize(NumPixels) is called from MGeomApply (which allocates
41// memory in a TClonesArray of type MCalibrationChargePix)
42// Initializes pointer to MBadPixelsCam
43//
44// Process(): Nothing to be done, histograms getting filled by MHCalibrationChargeCam
45//
46// PostProcess(): - FinalizePedestals()
47// - FinalizeCharges()
48// - FinalizeFFactorMethod()
49// - FinalizeBadPixels()
50// - FinalizeBlindCam()
51// - FinalizePINDiode()
52// - FinalizeFFactorQECam()
53// - FinalizeBlindPixelQECam()
54// - FinalizePINDiodeQECam()
55//
56// Input Containers:
57// MCalibrationChargeCam
58// MCalibrationChargeBlindPix
59// MCalibrationChargePINDiode
60// MCalibrationQECam
61// MPedestalCam
62// MBadPixelsCam
63// MGeomCam
64// MTime
65//
66// Output Containers:
67// MCalibrationChargeCam
68// MCalibrationChargeBlindPix
69// MCalibrationChargePINDiode
70// MCalibrationQECam
71// MBadPixelsCam
72//
73//
74// Preliminary description of the calibration in photons (email from 12/02/04)
75//
76// Why calibrating in photons:
77// ===========================
78//
79// At the Barcelona meeting in 2002, we decided to calibrate the camera in
80// photons. This for the following reasons:
81//
82// * The physical quantity arriving at the camera are photons. This is
83// the direct physical information from the air shower. The photons
84// have a flux and a spectrum.
85//
86// * The photon fluxes depend mostly on the shower energy (with
87// corrections deriving from the observation conditions), while the photon
88// spectra depend mostly on the observation conditions: zenith angle,
89// quality of the air, also the impact parameter of the shower.
90//
91// * The photomultiplier, in turn, has different response properties
92// (quantum efficiencies) for photons of different colour. (Moreover,
93// different pixels have slightly different quantum efficiencies).
94// The resulting number of photo-electrons is then amplified (linearly)
95// with respect to the photo-electron flux.
96//
97// * In the ideal case, one would like to disentagle the effects
98// of the observation conditions from the primary particle energy (which
99// one likes to measure). To do so, one needs:
100//
101// 1) A reliable calibration relating the FADC counts to the photo-electron
102// flux -> This is accomplished with the F-Factor method.
103//
104// 2) A reliable calibration of the wavelength-dependent quantum efficiency
105// -> This is accomplished with the combination of the three methods,
106// together with QE-measurements performed by David in order to do
107// the interpolation.
108//
109// 3) A reliable calibration of the observation conditions. This means:
110// - Tracing the atmospheric conditions -> LIDAR
111// - Tracing the observation zenith angle -> Drive System
112//
113// 4) Some knowlegde about the impact parameter:
114// - This is the only part which cannot be accomplished well with a
115// single telescope. We would thus need to convolute the spectrum
116// over the distribution of impact parameters.
117//
118//
119// How an ideal calibration would look like:
120// =========================================
121//
122// We know from the combined PIN-Diode and Blind-Pixel Method the response of
123// each pixel to well-measured light fluxes in three representative
124// wavelengths (green, blue, UV). We also know the response to these light
125// fluxes in photo-electrons. Thus, we can derive:
126//
127// - conversion factors to photo-electrons
128// - conversion factors to photons in three wavelengths.
129//
130// Together with David's measurements and some MC-simulation, we should be
131// able to derive tables for typical Cherenkov-photon spectra - convoluted
132// with the impact parameters and depending on the athmospheric conditions
133// and the zenith angle (the "outer parameters").
134//
135// From these tables we can create "calibration tables" containing some
136// effective quantum efficiency depending on these outer parameters and which
137// are different for each pixel.
138//
139// In an ideal MCalibrate, one would thus have to convert first the FADC
140// slices to Photo-electrons and then, depending on the outer parameters,
141// look up the effective quantum efficiency and get the mean number of
142// photons which is then used for the further analysis.
143//
144// How the (first) MAGIC calibration should look like:
145// ===================================================
146//
147// For the moment, we have only one reliable calibration method, although
148// with very large systematic errors. This is the F-Factor method. Knowing
149// that the light is uniform over the whole camera (which I would not at all
150// guarantee in the case of the CT1 pulser), one could in principle already
151// perform a relative calibration of the quantum efficiencies in the UV.
152// However, the spread in QE at UV is about 10-15% (according to the plot
153// that Abelardo sent around last time. The spread in photo-electrons is 15%
154// for the inner pixels, but much larger (40%) for the outer ones.
155//
156// I'm not sure if we can already say that we have measured the relative
157// difference in quantum efficiency for the inner pixels and produce a first
158// QE-table for each pixel. To so, I would rather check in other wavelengths
159// (which we can do in about one-two weeks when the optical transmission of
160// the calibration trigger is installed).
161//
162// Thus, for the moment being, I would join Thomas proposal to calibrate in
163// photo-electrons and apply one stupid average quantum efficiency for all
164// pixels. This keeping in mind that we will have much preciser information
165// in about one to two weeks.
166//
167//
168// What MCalibrate should calculate and what should be stored:
169// ===========================================================
170//
171// It is clear that in the end, MCerPhotEvt will store photons.
172// MCalibrationCam stores the conversionfactors to photo-electrons and also
173// some tables of how to apply the conversion to photons, given the outer
174// parameters. This is not yet implemented and not even discussed.
175//
176// To start, I would suggest that we define the "average quantum efficiency"
177// (maybe something like 25+-3%) and apply them equally to all
178// photo-electrons. Later, this average factor can be easily replaced by a
179// pixel-dependent factor and later by a (pixel-dependent) table.
180//
181//
182//
183//////////////////////////////////////////////////////////////////////////////
184#include "MCalibrationChargeCalc.h"
185
186#include <TSystem.h>
187#include <TH1.h>
188#include <TF1.h>
189
190#include "MLog.h"
191#include "MLogManip.h"
192
193#include "MParList.h"
194
195#include "MStatusDisplay.h"
196
197#include "MRawEvtHeader.h"
198
199#include "MGeomCam.h"
200#include "MGeomPix.h"
201#include "MHCamera.h"
202
203#include "MPedestalCam.h"
204#include "MPedestalPix.h"
205
206#include "MCalibrationIntensityChargeCam.h"
207#include "MCalibrationIntensityQECam.h"
208#include "MCalibrationIntensityBlindCam.h"
209
210#include "MHCalibrationChargeCam.h"
211#include "MHCalibrationChargeBlindCam.h"
212
213#include "MCalibrationChargeCam.h"
214#include "MCalibrationChargePix.h"
215#include "MCalibrationChargePINDiode.h"
216#include "MCalibrationBlindPix.h"
217#include "MCalibrationBlindCam.h"
218
219#include "MExtractedSignalCam.h"
220#include "MExtractedSignalPix.h"
221#include "MExtractedSignalBlindPixel.h"
222#include "MExtractedSignalPINDiode.h"
223
224#include "MBadPixelsIntensityCam.h"
225#include "MBadPixelsCam.h"
226
227#include "MCalibrationQECam.h"
228#include "MCalibrationQEPix.h"
229
230ClassImp(MCalibrationChargeCalc);
231
232using namespace std;
233
234const Float_t MCalibrationChargeCalc::fgChargeLimit = 2.5;
235const Float_t MCalibrationChargeCalc::fgChargeErrLimit = 0.;
236const Float_t MCalibrationChargeCalc::fgChargeRelErrLimit = 1.;
237const Float_t MCalibrationChargeCalc::fgLambdaErrLimit = 0.2;
238const Float_t MCalibrationChargeCalc::fgLambdaCheckLimit = 0.5;
239const Float_t MCalibrationChargeCalc::fgPheErrLimit = 4.5;
240const Float_t MCalibrationChargeCalc::fgFFactorErrLimit = 4.5;
241// --------------------------------------------------------------------------
242//
243// Default constructor.
244//
245// Sets the pointer to fQECam and fGeom to NULL
246//
247// Calls AddToBranchList for:
248// - MRawEvtData.fHiGainPixId
249// - MRawEvtData.fLoGainPixId
250// - MRawEvtData.fHiGainFadcSamples
251// - MRawEvtData.fLoGainFadcSamples
252//
253// Initializes:
254// - fChargeLimit to fgChargeLimit
255// - fChargeErrLimit to fgChargeErrLimit
256// - fChargeRelErrLimit to fgChargeRelErrLimit
257// - fFFactorErrLimit to fgFFactorErrLimit
258// - fLambdaCheckLimit to fgLambdaCheckLimit
259// - fLambdaErrLimit to fgLambdaErrLimit
260// - fPheErrLimit to fgPheErrLimit
261// - fPulserColor to MCalibrationCam::kCT1
262// - fOutputPath to "."
263// - fOutputFile to "ChargeCalibStat.txt"
264// - flag debug to kFALSE
265//
266// Sets all checks
267//
268// Calls:
269// - Clear()
270//
271MCalibrationChargeCalc::MCalibrationChargeCalc(const char *name, const char *title)
272 : fGeom(NULL), fSignal(NULL), fHeader(NULL)
273{
274
275 fName = name ? name : "MCalibrationChargeCalc";
276 fTitle = title ? title : "Task to calculate the calibration constants and MCalibrationCam ";
277
278 AddToBranchList("MRawEvtData.fHiGainPixId");
279 AddToBranchList("MRawEvtData.fLoGainPixId");
280 AddToBranchList("MRawEvtData.fHiGainFadcSamples");
281 AddToBranchList("MRawEvtData.fLoGainFadcSamples");
282
283 SetChargeLimit ();
284 SetChargeErrLimit ();
285 SetChargeRelErrLimit ();
286 SetFFactorErrLimit ();
287 SetLambdaCheckLimit ();
288 SetLambdaErrLimit ();
289 SetPheErrLimit ();
290 SetOutputPath ();
291 SetOutputFile ();
292 SetDebug ( kFALSE );
293
294 SetCheckDeadPixels ();
295 SetCheckDeviatingBehavior();
296 SetCheckExtractionWindow ();
297 SetCheckHistOverflow ();
298 SetCheckOscillations ();
299
300 Clear();
301
302}
303
304// --------------------------------------------------------------------------
305//
306// Sets:
307// - all variables to 0.,
308// - all flags to kFALSE
309// - all pointers to NULL
310// - the pulser colour to kNONE
311// - fBlindPixelFlags to 0
312// - fPINDiodeFlags to 0
313//
314void MCalibrationChargeCalc::Clear(const Option_t *o)
315{
316
317 fNumHiGainSamples = 0.;
318 fNumLoGainSamples = 0.;
319 fSqrtHiGainSamples = 0.;
320 fSqrtLoGainSamples = 0.;
321 fNumInnerFFactorMethodUsed = 0;
322
323 fIntensBad = NULL;
324 fBadPixels = NULL;
325 fIntensCam = NULL;
326 fCam = NULL;
327 fHCam = NULL;
328 fIntensQE = NULL;
329 fQECam = NULL;
330 fIntensBlind = NULL;
331 fBlindCam = NULL;
332 fHBlindCam = NULL;
333 fPINDiode = NULL;
334 fPedestals = NULL;
335
336 fPulserPattern = 0;
337 SetPulserColor ( MCalibrationCam::kNONE );
338
339 fBlindPixelFlags.Set(0);
340 fPINDiodeFlags .Set(0);
341 fResultFlags .Set(0);
342}
343
344
345// -----------------------------------------------------------------------------------
346//
347// The following container are searched for and execution aborted if not in MParList:
348// - MRawEvtHeader
349// - MPedestalCam
350// - MExtractedSignalCam
351//
352Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
353{
354
355 fHeader = (MRawEvtHeader*)pList->FindObject("MRawEvtHeader");
356 if (!fHeader)
357 {
358 *fLog << err << "MRawEvtHeader not found... abort." << endl;
359 return kFALSE;
360 }
361
362 //
363 // Containers that have to be there.
364 //
365 fPedestals = (MPedestalCam*)pList->FindObject("MPedestalCam");
366 if (!fPedestals)
367 {
368 *fLog << err << "MPedestalCam not found... aborting" << endl;
369 return kFALSE;
370 }
371
372 fSignal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
373 if (!fSignal)
374 {
375 *fLog << err << "MExtractedSignalCam not found... aborting" << endl;
376 return kFALSE;
377 }
378
379
380 return kTRUE;
381}
382
383
384// --------------------------------------------------------------------------
385//
386// Search for the following input containers and abort if not existing:
387// - MGeomCam
388// - MCalibrationIntensityChargeCam or MCalibrationChargeCam
389// - MCalibrationIntensityQECam or MCalibrationQECam
390// - MBadPixelsIntensityCam or MBadPixelsCam
391//
392// Search for the following input containers and give a warning if not existing:
393// - MCalibrationBlindPix
394// - MCalibrationChargePINDiode
395//
396// It retrieves the following variables from MCalibrationChargeCam:
397//
398// - fNumHiGainSamples
399// - fNumLoGainSamples
400//
401// It defines the PixId of every pixel in:
402//
403// - MCalibrationIntensityChargeCam
404// - MCalibrationChargeCam
405// - MCalibrationQECam
406//
407// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
408//
409// - MCalibrationChargePix
410// - MCalibrationQEPix
411//
412// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
413//
414// - MCalibrationChargeCam
415// - MCalibrationBlindPix (if existing)
416// - MCalibrationChargePINDiode (if existing)
417//
418Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
419{
420
421 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
422 if (!fGeom)
423 {
424 *fLog << err << "No MGeomCam found... aborting." << endl;
425 return kFALSE;
426 }
427
428 fIntensCam = (MCalibrationIntensityChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityChargeCam"));
429 if (fIntensCam)
430 *fLog << inf << "Found MCalibrationIntensityChargeCam ... " << endl;
431 else
432 {
433 fCam = (MCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
434 if (!fCam)
435 {
436 *fLog << err << "Cannot find MCalibrationChargeCam ... abort." << endl;
437 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
438 return kFALSE;
439 }
440 }
441
442 fHCam = (MHCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeCam"));
443 if (!fHCam)
444 {
445 *fLog << err << "Cannot find MHCalibrationChargeCam ... abort." << endl;
446 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
447 return kFALSE;
448 }
449
450 fIntensQE = (MCalibrationIntensityQECam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityQECam"));
451 if (fIntensQE)
452 *fLog << inf << "Found MCalibrationIntensityQECam ... " << endl;
453 else
454 {
455 fQECam = (MCalibrationQECam*)pList->FindObject(AddSerialNumber("MCalibrationQECam"));
456 if (!fQECam)
457 {
458 *fLog << err << "Cannot find MCalibrationQECam ... abort." << endl;
459 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationQECam before..." << endl;
460 return kFALSE;
461 }
462 }
463
464 fIntensBlind = (MCalibrationIntensityBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityBlindCam"));
465 if (fIntensBlind)
466 *fLog << inf << "Found MCalibrationIntensityBlindCam ... " << endl;
467 else
468 {
469 fBlindCam = (MCalibrationBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationBlindCam"));
470 if (!fBlindCam)
471 {
472 *fLog << endl;
473 *fLog << warn << GetDescriptor()
474 << ": No MCalibrationBlindCam found... no Blind Pixel method! " << endl;
475 }
476 }
477
478 fHBlindCam = (MHCalibrationChargeBlindCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeBlindCam"));
479 if (!fHBlindCam)
480 {
481 *fLog << endl;
482 *fLog << warn << GetDescriptor()
483 << ": No MHCalibrationChargeBlindCam found... no Blind Pixel method! " << endl;
484 }
485
486 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
487 if (fIntensBad)
488 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
489 else
490 {
491 fBadPixels = (MBadPixelsCam*)pList->FindObject(AddSerialNumber("MBadPixelsCam"));
492 if (!fBadPixels)
493 {
494 *fLog << err << "Cannot find MBadPixelsCam ... abort." << endl;
495 return kFALSE;
496 }
497 }
498
499 //
500 // Optional Containers
501 //
502 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
503 if (!fPINDiode)
504 {
505 *fLog << endl;
506 *fLog << warn << GetDescriptor()
507 << ": MCalibrationChargePINDiode not found... no PIN Diode method! " << endl;
508 }
509
510 fNumHiGainSamples = fSignal->GetNumUsedHiGainFADCSlices();
511 fNumLoGainSamples = fSignal->GetNumUsedLoGainFADCSlices();
512
513 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
514 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
515
516
517 MCalibrationQECam *qecam = fIntensQE
518 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
519 MCalibrationChargeCam *chargecam = fIntensCam
520 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
521 MBadPixelsCam *badcam = fIntensBad
522 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
523
524 UInt_t npixels = fGeom->GetNumPixels();
525
526 for (UInt_t i=0; i<npixels; i++)
527 {
528
529 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
530 MCalibrationQEPix &pqe = (MCalibrationQEPix&) (*qecam) [i];
531 MBadPixelsPix &bad = (*badcam) [i];
532
533 pix.SetPixId(i);
534 pqe.SetPixId(i);
535
536 if (bad.IsBad())
537 {
538 pix.SetExcluded();
539 pqe.SetExcluded();
540 continue;
541 }
542
543 if (IsDebug())
544 pix.SetDebug();
545 }
546
547 return kTRUE;
548}
549
550// ----------------------------------------------------------------------------------
551//
552// Set the correct colour to the charge containers
553//
554Int_t MCalibrationChargeCalc::Process()
555{
556
557 const UInt_t pattern = fHeader->GetCalibrationPattern();
558
559 if (pattern == fPulserPattern)
560 return kTRUE;
561
562 enum ColorCode_t
563 {
564 kSlot1Green = BIT(0),
565 kSlot2Green = BIT(1),
566 kSlot3Blue = BIT(2),
567 kSlot4UV = BIT(3),
568 kSlot5UV = BIT(4),
569 kSlot6Blue = BIT(5),
570 kSlot7Blue = BIT(6),
571 kSlot8Blue = BIT(7),
572 kSlot9AttBlue = BIT(8),
573 kSlot10Blue = BIT(9),
574 kSlot11Blue = BIT(10),
575 kSlot12UV = BIT(11),
576 kSlot13UV = BIT(12),
577 kSlot14Blue = BIT(13),
578 kSlot15Green = BIT(14),
579 kSlot16AttGreen = BIT(15),
580 kCT1Pulser = BIT(16),
581 kAnyGreen = kSlot1Green | kSlot2Green | kSlot15Green | kSlot16AttGreen,
582 kAnyUV = kSlot4UV | kSlot5UV | kSlot12UV | kSlot13UV,
583 kAnyBlue = kSlot3Blue | kSlot6Blue | kSlot7Blue | kSlot8Blue
584 | kSlot9AttBlue| kSlot10Blue | kSlot11Blue | kSlot14Blue
585 };
586
587 //
588 // The pattern has changed, we have to initialize everything new!!!
589 //
590 *fLog << inf << GetDescriptor() << "- New pulser pattern: " ;
591 for (Int_t i=16; i>= 0; i--)
592 *fLog << (pattern >> i & 1);
593 *fLog << endl;
594 fPulserPattern = pattern;
595
596 //
597 // Now retrieve the colour and check if not various colours have been used
598 //
599 fPulserColor = MCalibrationCam::kNONE;
600
601 if (fPulserPattern & kAnyGreen)
602 fPulserColor = MCalibrationCam::kGREEN;
603
604 if ((fPulserPattern & kAnyBlue ||
605 fPulserPattern & kAnyUV ||
606 fPulserPattern & kCT1Pulser) && fPulserColor != MCalibrationCam::kNONE)
607 {
608 *fLog << err << "Multiple colours used simultaneously in calibration file. Will skip this part!" << endl;
609 return kCONTINUE;
610 }
611
612 if (fPulserColor == MCalibrationCam::kNONE)
613 {
614 if (fPulserPattern & kAnyBlue)
615 fPulserColor = MCalibrationCam::kBLUE;
616 if (fPulserPattern & kAnyUV)
617 fPulserColor = MCalibrationCam::kUV;
618 if (fPulserPattern & kCT1Pulser)
619 fPulserColor = MCalibrationCam::kCT1;
620 }
621
622 *fLog << inf << "Found new colour ... " << flush;
623
624 switch (fPulserColor)
625 {
626 case MCalibrationCam::kGREEN: *fLog << "Green."; break;
627 case MCalibrationCam::kBLUE: *fLog << "Blue."; break;
628 case MCalibrationCam::kUV: *fLog << "UV."; break;
629 case MCalibrationCam::kCT1: *fLog << "CT1."; break;
630 default: break;
631 }
632 *fLog << endl;
633
634 MCalibrationBlindCam *blindcam = fIntensBlind
635 ? (MCalibrationBlindCam*) fIntensBlind->GetCam() : fBlindCam;
636 MCalibrationChargeCam *chargecam = fIntensCam
637 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
638
639 //
640 // Initialize the pulser colours
641 //
642 chargecam ->SetPulserColor( fPulserColor );
643 blindcam ->SetPulserColor( fPulserColor );
644 fHCam ->SetColor ( fPulserColor );
645 fHBlindCam->SetColor ( fPulserColor );
646
647 if (fPINDiode)
648 fPINDiode->SetColor( fPulserColor );
649
650 return kTRUE;
651}
652
653// -----------------------------------------------------------------------
654//
655// Analogue to the MTask::CallPostProcess, but without the manipulation
656// of the bit fIsPreProcessed. Needed in order to call PostProcess multiple
657// times in the intensity calibration.
658//
659Int_t MCalibrationChargeCalc::CallPostProcess()
660{
661
662 // if (!fIsPreprocessed)
663 // return kTRUE;
664
665 *fLog << all << fName << "... " << flush;
666 if (fDisplay)
667 fDisplay->SetStatusLine2(*this);
668
669 return PostProcess();
670}
671
672// -----------------------------------------------------------------------
673//
674// Return if number of executions is null.
675//
676// First loop over pixels, average areas and sectors, call:
677// - FinalizePedestals()
678// - FinalizeCharges()
679// for every entry. Count number of valid pixels in loop and return kFALSE
680// if there are none (the "Michele check").
681//
682// Call FinalizeBadPixels()
683//
684// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
685//
686// Call FinalizeBlindCam()
687// Call FinalizePINDiode()
688//
689// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
690// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
691// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
692//
693// Call FinalizeUnsuitablePixels()
694//
695// Call MParContainer::SetReadyToSave() for fIntensCam, fCam, fQECam, fBadPixels and
696// fBlindCam and fPINDiode if they exist
697//
698// Print out some statistics
699//
700Int_t MCalibrationChargeCalc::PostProcess()
701{
702
703 if (GetNumExecutions()==0)
704 return kFALSE;
705
706 *fLog << endl;
707
708 if (fPINDiode)
709 if (!fPINDiode->IsValid())
710 {
711 *fLog << warn << GetDescriptor()
712 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
713 fPINDiode = NULL;
714 }
715
716 *fLog << endl;
717
718 MCalibrationBlindCam *blindcam = fIntensBlind
719 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
720 MCalibrationQECam *qecam = fIntensQE
721 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
722 MCalibrationChargeCam *chargecam = fIntensCam
723 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
724 MBadPixelsCam *badcam = fIntensBad
725 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
726
727 //
728 // First loop over pixels, call FinalizePedestals and FinalizeCharges
729 //
730 Int_t nvalid = 0;
731
732 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
733 {
734
735 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[pixid];
736 //
737 // Check if the pixel has been excluded from the fits
738 //
739 if (pix.IsExcluded())
740 continue;
741
742 MPedestalPix &ped = (*fPedestals)[pixid];
743 MBadPixelsPix &bad = (*badcam) [pixid];
744
745 const Int_t aidx = (*fGeom)[pixid].GetAidx();
746
747 FinalizePedestals(ped,pix,aidx);
748
749 if (FinalizeCharges(pix,bad,"pixel "))
750 nvalid++;
751 }
752
753 *fLog << endl;
754 //
755 // The Michele check ...
756 //
757 if (nvalid == 0)
758 {
759 *fLog << err << GetDescriptor() << ": All pixels have non-valid calibration. "
760 << "Did you forget to fill the histograms "
761 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
762 *fLog << err << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
763 << "instead of a calibration run " << endl;
764 return kFALSE;
765 }
766
767 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
768 {
769
770 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
771 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
772
773 FinalizePedestals(ped,pix,aidx);
774 FinalizeCharges(pix,
775 fIntensCam ? fIntensCam->GetAverageBadArea(aidx) : fCam->GetAverageBadArea(aidx),
776 "area id");
777 }
778
779 *fLog << endl;
780
781 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
782 {
783
784 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
785
786 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
787 FinalizePedestals(ped,pix, 0);
788 }
789
790 *fLog << endl;
791
792 //
793 // Finalize Bad Pixels
794 //
795 FinalizeBadPixels();
796
797 //
798 // Finalize F-Factor method
799 //
800 if (FinalizeFFactorMethod())
801 chargecam->SetFFactorMethodValid(kTRUE);
802 else
803 {
804 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
805 chargecam->SetFFactorMethodValid(kFALSE);
806 return kFALSE;
807 }
808
809 *fLog << endl;
810 //
811 // Finalize Blind Pixel
812 //
813 qecam->SetBlindPixelMethodValid(FinalizeBlindCam());
814
815 //
816 // Finalize PIN Diode
817 //
818 qecam->SetBlindPixelMethodValid(FinalizePINDiode());
819
820 //
821 // Finalize QE Cam
822 //
823 FinalizeFFactorQECam();
824 FinalizeBlindPixelQECam();
825 FinalizePINDiodeQECam();
826 FinalizeCombinedQECam();
827
828 //
829 // Re-direct the output to an ascii-file from now on:
830 //
831 MLog asciilog;
832 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
833 SetLogStream(&asciilog);
834 //
835 // Finalize calibration statistics
836 //
837 FinalizeUnsuitablePixels();
838
839 chargecam->SetReadyToSave();
840 qecam ->SetReadyToSave();
841 badcam ->SetReadyToSave();
842
843 if (blindcam)
844 blindcam->SetReadyToSave();
845 if (fPINDiode)
846 fPINDiode->SetReadyToSave();
847
848 *fLog << inf << endl;
849 *fLog << GetDescriptor() << ": Fatal errors statistics:" << endl;
850
851 PrintUncalibrated(MBadPixelsPix::kChargeIsPedestal,
852 Form("%s%2.1f%s","Signal less than ",fChargeLimit," Pedestal RMS: "));
853 PrintUncalibrated(MBadPixelsPix::kChargeRelErrNotValid,
854 Form("%s%2.1f%s","Signal Error bigger than ",fChargeRelErrLimit," times Mean Signal: "));
855 PrintUncalibrated(MBadPixelsPix::kLoGainSaturation,
856 "Pixels with Low Gain Saturation: ");
857 PrintUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin,
858 Form("%s%2.1f%s","Mean Abs. Arr. Time in First ",1.," Bin(s): "));
859 PrintUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins,
860 Form("%s%2.1f%s","Mean Abs. Arr. Time in Last ",2.," Bin(s): "));
861 PrintUncalibrated(MBadPixelsPix::kDeviatingNumPhes,
862 "Pixels with deviating number of phes: ");
863 PrintUncalibrated(MBadPixelsPix::kDeviatingFFactor,
864 "Pixels with deviating F-Factor: ");
865
866 *fLog << inf << endl;
867 *fLog << GetDescriptor() << ": Unreliable errors statistics:" << endl;
868
869 PrintUncalibrated(MBadPixelsPix::kChargeSigmaNotValid,
870 "Signal Sigma smaller than Pedestal RMS: ");
871 PrintUncalibrated(MBadPixelsPix::kHiGainOscillating,
872 "Pixels with changing Hi Gain signal over time: ");
873 PrintUncalibrated(MBadPixelsPix::kLoGainOscillating,
874 "Pixels with changing Lo Gain signal over time: ");
875 PrintUncalibrated(MBadPixelsPix::kHiGainNotFitted,
876 "Pixels with unsuccesful Gauss fit to the Hi Gain: ");
877 PrintUncalibrated(MBadPixelsPix::kLoGainNotFitted,
878 "Pixels with unsuccesful Gauss fit to the Lo Gain: ");
879
880 SetLogStream(&gLog);
881
882 return kTRUE;
883}
884
885// ----------------------------------------------------------------------------------
886//
887// Retrieves pedestal and pedestal RMS from MPedestalPix
888// Retrieves total entries from MPedestalCam
889// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
890// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
891//
892// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
893// - MCalibrationChargePix::CalcLoGainPedestal()
894//
895void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
896{
897
898 //
899 // get the pedestals
900 //
901 const Float_t pedes = ped.GetPedestal();
902 const Float_t prms = ped.GetPedestalRms();
903 const Int_t num = fPedestals->GetTotalEntries();
904
905 //
906 // RMS error set by PedCalcFromLoGain, 0 in case MPedCalcPedRun was used.
907 //
908 const Float_t prmserr = num>0 ? prms/TMath::Sqrt(2.*num) : ped.GetPedestalRmsError();
909
910 //
911 // set them in the calibration camera
912 //
913 if (cal.IsHiGainSaturation())
914 {
915 cal.SetPedestal(pedes * fNumLoGainSamples,
916 prms * fSqrtLoGainSamples,
917 prmserr * fSqrtLoGainSamples);
918 cal.CalcLoGainPedestal((Float_t)fNumLoGainSamples, aidx);
919 }
920 else
921 {
922 cal.SetPedestal(pedes * fNumHiGainSamples,
923 prms * fSqrtHiGainSamples,
924 prmserr * fSqrtHiGainSamples);
925 }
926
927}
928
929// ----------------------------------------------------------------------------------------------------
930//
931// Check fit results validity. Bad Pixels flags are set if:
932//
933// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
934// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
935// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
936// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
937// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
938//
939// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
940//
941// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
942// and returns kFALSE if not succesful.
943//
944// Calls MCalibrationChargePix::CalcFFactor() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
945// and returns kFALSE if not succesful.
946//
947// Calls MCalibrationChargePix::CalcConvFFactor()and sets flag: MBadPixelsPix::kDeviatingNumPhes)
948// and returns kFALSE if not succesful.
949//
950Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad, const char* what)
951{
952
953 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
954 return kFALSE;
955
956 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
957 {
958 *fLog << warn
959 << Form("Fitted Charge: %5.2f < %2.1f",cal.GetMean(),fChargeLimit)
960 << Form(" * Pedestal RMS %5.2f in %s%3i",cal.GetPedRms(),what,cal.GetPixId()) << endl;
961 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
962 }
963
964 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
965 {
966 *fLog << warn
967 << Form("Fitted Charge: %4.2f < %2.1f",cal.GetMean(),fChargeRelErrLimit)
968 << Form(" * its error %4.2f in %s%3i",cal.GetMeanErr(),what,cal.GetPixId()) << endl;
969 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
970 }
971
972 if (cal.GetSigma() < cal.GetPedRms())
973 {
974 *fLog << warn
975 << Form("Sigma of Fitted Charge: %6.2f <",cal.GetSigma())
976 << Form(" Ped. RMS=%5.2f in %s%3i",cal.GetPedRms(),what,cal.GetPixId()) << endl;
977 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
978 return kFALSE;
979 }
980
981 if (!cal.CalcReducedSigma())
982 {
983 *fLog << warn
984 << Form("Could not calculate the reduced sigma in %s: ",what)
985 << Form(" %4i",cal.GetPixId())
986 << endl;
987 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
988 return kFALSE;
989 }
990
991 if (!cal.CalcFFactor())
992 {
993 *fLog << warn
994 << Form("Could not calculate the F-Factor in %s: ",what)
995 << Form(" %4i",cal.GetPixId())
996 << endl;
997 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
998 return kFALSE;
999 }
1000
1001 if (!cal.CalcConvFFactor())
1002 {
1003 *fLog << warn
1004 << Form("Could not calculate the Conv. FADC counts to Phes in %s: ",what)
1005 << Form(" %4i",cal.GetPixId())
1006 << endl;
1007 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1008 return kFALSE;
1009 }
1010
1011 return kTRUE;
1012}
1013
1014
1015
1016// -----------------------------------------------------------------------------------
1017//
1018// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
1019// - MBadPixelsPix::kChargeIsPedestal
1020// - MBadPixelsPix::kChargeRelErrNotValid
1021// - MBadPixelsPix::kMeanTimeInFirstBin
1022// - MBadPixelsPix::kMeanTimeInLast2Bins
1023// - MBadPixelsPix::kDeviatingNumPhes
1024// - MBadPixelsPix::kHiGainOverFlow
1025// - MBadPixelsPix::kLoGainOverFlow
1026//
1027// - Call MCalibrationPix::SetExcluded() for the bad pixels
1028//
1029// Sets pixel to MBadPixelsPix::kUnreliableRun, if one of the following flags is set:
1030// - MBadPixelsPix::kChargeSigmaNotValid
1031//
1032void MCalibrationChargeCalc::FinalizeBadPixels()
1033{
1034
1035 MBadPixelsCam *badcam = fIntensBad
1036 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1037 MCalibrationChargeCam *chargecam = fIntensCam
1038 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1039
1040 for (Int_t i=0; i<badcam->GetSize(); i++)
1041 {
1042
1043 MBadPixelsPix &bad = (*badcam) [i];
1044 MCalibrationPix &pix = (*chargecam)[i];
1045
1046 if (IsCheckDeadPixels())
1047 {
1048 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
1049 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1050
1051 if (bad.IsUncalibrated( MBadPixelsPix::kChargeErrNotValid ))
1052 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1053
1054 if (bad.IsUncalibrated( MBadPixelsPix::kChargeRelErrNotValid ))
1055 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1056 }
1057
1058 if (IsCheckExtractionWindow())
1059 {
1060 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
1061 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1062
1063 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
1064 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1065 }
1066
1067 if (IsCheckDeviatingBehavior())
1068 {
1069 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
1070 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1071 }
1072
1073 if (IsCheckHistOverflow())
1074 {
1075 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOverFlow ))
1076 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1077
1078 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOverFlow ))
1079 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1080 }
1081
1082 if (bad.IsUnsuitable( MBadPixelsPix::kUnsuitableRun ))
1083 pix.SetExcluded();
1084
1085 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
1086 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
1087 }
1088}
1089
1090// ------------------------------------------------------------------------
1091//
1092//
1093// First loop: Calculate a mean and mean RMS of photo-electrons per area index
1094// Include only pixels which are not MBadPixelsPix::kUnsuitableRun nor
1095// MBadPixelsPix::kChargeSigmaNotValid (see FinalizeBadPixels()) and set
1096// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
1097//
1098// Second loop: Get mean number of photo-electrons and its RMS including
1099// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
1100// and further exclude those deviating by more than fPheErrLimit mean
1101// sigmas from the mean (obtained in first loop). Set
1102// MBadPixelsPix::kDeviatingNumPhes if excluded.
1103//
1104// For the suitable pixels with flag MBadPixelsPix::kChargeSigmaNotValid
1105// set the number of photo-electrons as the mean number of photo-electrons
1106// calculated in that area index.
1107//
1108// Set weighted mean and variance of photo-electrons per area index in:
1109// average area pixels of MCalibrationChargeCam (obtained from:
1110// MCalibrationChargeCam::GetAverageArea() )
1111//
1112// Set weighted mean and variance of photo-electrons per sector in:
1113// average sector pixels of MCalibrationChargeCam (obtained from:
1114// MCalibrationChargeCam::GetAverageSector() )
1115//
1116//
1117// Third loop: Set mean number of photo-electrons and its RMS in the pixels
1118// only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1119//
1120Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
1121{
1122
1123 MBadPixelsCam *badcam = fIntensBad
1124 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1125 MCalibrationChargeCam *chargecam = fIntensCam
1126 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1127
1128 const Int_t npixels = fGeom->GetNumPixels();
1129 const Int_t nareas = fGeom->GetNumAreas();
1130 const Int_t nsectors = fGeom->GetNumSectors();
1131
1132 TArrayF lowlim (nareas);
1133 TArrayF upplim (nareas);
1134 TArrayD areavars (nareas);
1135 TArrayD areaweights (nareas);
1136 TArrayD sectorweights (nsectors);
1137 TArrayD areaphes (nareas);
1138 TArrayD sectorphes (nsectors);
1139 TArrayI numareavalid (nareas);
1140 TArrayI numsectorvalid(nsectors);
1141
1142 //
1143 // First loop: Get mean number of photo-electrons and the RMS
1144 // The loop is only to recognize later pixels with very deviating numbers
1145 //
1146 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
1147
1148 for (Int_t i=0; i<npixels; i++)
1149 {
1150
1151 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1152 MBadPixelsPix &bad = (*badcam)[i];
1153
1154 if (!pix.IsFFactorMethodValid())
1155 continue;
1156
1157 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1158 {
1159 pix.SetFFactorMethodValid(kFALSE);
1160 continue;
1161 }
1162
1163 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1164 continue;
1165
1166 const Float_t nphe = pix.GetPheFFactorMethod();
1167 const Int_t aidx = (*fGeom)[i].GetAidx();
1168
1169 camphes.Fill(i,nphe);
1170 camphes.SetUsed(i);
1171
1172 areaphes [aidx] += nphe;
1173 areavars [aidx] += nphe*nphe;
1174 numareavalid[aidx] ++;
1175 }
1176
1177 for (Int_t i=0; i<nareas; i++)
1178 {
1179 if (numareavalid[i] == 0)
1180 {
1181 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
1182 << "in area index: " << i << endl;
1183 continue;
1184 }
1185
1186 if (numareavalid[i] == 1)
1187 areavars[i] = 0.;
1188 else if (numareavalid[i] == 0)
1189 {
1190 areaphes[i] = -1.;
1191 areaweights[i] = -1.;
1192 }
1193 else
1194 {
1195 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1);
1196 areaphes[i] = areaphes[i] / numareavalid[i];
1197 }
1198
1199 if (areavars[i] < 0.)
1200 {
1201 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of photo-electrons found "
1202 << "in area index: " << i << endl;
1203 continue;
1204 }
1205
1206 lowlim [i] = areaphes[i] - fPheErrLimit*TMath::Sqrt(areavars[i]);
1207 upplim [i] = areaphes[i] + fPheErrLimit*TMath::Sqrt(areavars[i]);
1208
1209 TH1D *hist = camphes.ProjectionS(TArrayI(),TArrayI(1,&i),"_py",100);
1210 hist->Fit("gaus","Q");
1211 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1212 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1213 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1214
1215 if (IsDebug())
1216 hist->DrawClone();
1217
1218 if (ndf < 2)
1219 {
1220 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1221 << "in the camera with area index: " << i << endl;
1222 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1223 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1224 delete hist;
1225 continue;
1226 }
1227
1228 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1229
1230 if (prob < 0.001)
1231 {
1232 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1233 << "in the camera with area index: " << i << endl;
1234 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1235 << " is smaller than 0.001 " << endl;
1236 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1237 delete hist;
1238 continue;
1239 }
1240
1241 *fLog << inf << GetDescriptor() << ": Mean number of photo-electrons "
1242 << "with area idx " << i << ": "
1243 << Form("%7.2f+-%6.2f",mean,sigma) << endl;
1244
1245 lowlim [i] = mean - fPheErrLimit*sigma;
1246 upplim [i] = mean + fPheErrLimit*sigma;
1247
1248 delete hist;
1249 }
1250
1251 *fLog << endl;
1252
1253 numareavalid.Reset();
1254 areaphes .Reset();
1255 areavars .Reset();
1256 //
1257 // Second loop: Get mean number of photo-electrons and its RMS excluding
1258 // pixels deviating by more than fPheErrLimit sigma.
1259 // Set the conversion factor FADC counts to photo-electrons
1260 //
1261 for (Int_t i=0; i<npixels; i++)
1262 {
1263
1264 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1265
1266 if (!pix.IsFFactorMethodValid())
1267 continue;
1268
1269 MBadPixelsPix &bad = (*badcam)[i];
1270
1271 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1272 continue;
1273
1274 const Float_t nvar = pix.GetPheFFactorMethodVar();
1275
1276 if (nvar <= 0.)
1277 {
1278 pix.SetFFactorMethodValid(kFALSE);
1279 continue;
1280 }
1281
1282 const Int_t aidx = (*fGeom)[i].GetAidx();
1283 const Int_t sector = (*fGeom)[i].GetSector();
1284 const Float_t area = (*fGeom)[i].GetA();
1285 const Float_t nphe = pix.GetPheFFactorMethod();
1286
1287 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1288 {
1289 *fLog << warn << GetDescriptor() << ": Number of phes: "
1290 << Form("%7.2f out of %3.1f sigma limit: ",nphe,fPheErrLimit)
1291 << Form("[%7.2f,%7.2f] pixel%4i",lowlim[aidx],upplim[aidx],i) << endl;
1292 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1293 if (IsCheckDeviatingBehavior())
1294 {
1295 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1296 pix.SetFFactorMethodValid(kFALSE);
1297 }
1298 continue;
1299 }
1300
1301 areaweights [aidx] += nphe*nphe;
1302 areaphes [aidx] += nphe;
1303 numareavalid [aidx] ++;
1304
1305 if (aidx == 0)
1306 fNumInnerFFactorMethodUsed++;
1307
1308 sectorweights [sector] += nphe*nphe/area/area;
1309 sectorphes [sector] += nphe/area;
1310 numsectorvalid[sector] ++;
1311 }
1312
1313 *fLog << endl;
1314
1315 for (Int_t aidx=0; aidx<nareas; aidx++)
1316 {
1317
1318 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1319
1320 if (numareavalid[aidx] == 1)
1321 areaweights[aidx] = 0.;
1322 else if (numareavalid[aidx] == 0)
1323 {
1324 areaphes[aidx] = -1.;
1325 areaweights[aidx] = -1.;
1326 }
1327 else
1328 {
1329 areaweights[aidx] = (areaweights[aidx] - areaphes[aidx]*areaphes[aidx]/numareavalid[aidx])
1330 / (numareavalid[aidx]-1);
1331 areaphes[aidx] /= numareavalid[aidx];
1332 }
1333
1334 if (areaweights[aidx] < 0. || areaphes[aidx] <= 0.)
1335 {
1336 *fLog << warn << GetDescriptor()
1337 << ": Mean number phes from area index " << aidx << " could not be calculated: "
1338 << " Mean: " << areaphes[aidx]
1339 << " Variance: " << areaweights[aidx] << endl;
1340 apix.SetFFactorMethodValid(kFALSE);
1341 continue;
1342 }
1343
1344 *fLog << inf << GetDescriptor()
1345 << ": Average total number phes in area idx " << aidx << ": "
1346 << Form("%7.2f%s%6.2f",areaphes[aidx]," +- ",TMath::Sqrt(areaweights[aidx])) << endl;
1347
1348 apix.SetPheFFactorMethod ( areaphes[aidx] );
1349 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1350 apix.SetFFactorMethodValid ( kTRUE );
1351
1352 }
1353
1354 *fLog << endl;
1355
1356 for (Int_t sector=0; sector<nsectors; sector++)
1357 {
1358
1359 if (numsectorvalid[sector] == 1)
1360 sectorweights[sector] = 0.;
1361 else if (numsectorvalid[sector] == 0)
1362 {
1363 sectorphes[sector] = -1.;
1364 sectorweights[sector] = -1.;
1365 }
1366 else
1367 {
1368 sectorweights[sector] = (sectorweights[sector]
1369 - sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector]
1370 )
1371 / (numsectorvalid[sector]-1.);
1372 sectorphes[sector] /= numsectorvalid[sector];
1373 }
1374
1375 MCalibrationChargePix &spix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
1376
1377 if (sectorweights[sector] < 0. || sectorphes[sector] <= 0.)
1378 {
1379 *fLog << warn << GetDescriptor()
1380 <<": Mean number phes per area for sector " << sector << " could not be calculated: "
1381 << " Mean: " << sectorphes[sector]
1382 << " Variance: " << sectorweights[sector] << endl;
1383 spix.SetFFactorMethodValid(kFALSE);
1384 continue;
1385 }
1386
1387 *fLog << inf << GetDescriptor()
1388 << ": Average number phes per area in sector " << sector << ": "
1389 << Form("%5.3f+-%4.3f [phe/mm^2]",sectorphes[sector],TMath::Sqrt(sectorweights[sector]))
1390 << endl;
1391
1392 spix.SetPheFFactorMethod ( sectorphes[sector] );
1393 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1394 spix.SetFFactorMethodValid ( kTRUE );
1395
1396 }
1397
1398 //
1399 // Third loop: Set mean number of photo-electrons and its RMS in the pixels
1400 // only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1401 //
1402 for (Int_t i=0; i<npixels; i++)
1403 {
1404
1405 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1406 MBadPixelsPix &bad = (*badcam)[i];
1407
1408 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1409 continue;
1410
1411 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1412 {
1413 const Int_t aidx = (*fGeom)[i].GetAidx();
1414 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1415
1416 pix.SetPheFFactorMethod ( apix.GetPheFFactorMethod() );
1417 pix.SetPheFFactorMethodVar( apix.GetPheFFactorMethodVar() );
1418 if (!pix.CalcConvFFactor())
1419 {
1420 *fLog << warn << GetDescriptor()
1421 << ": Could not calculate the Conv. FADC counts to Phes in pixel: "
1422 << Form(" %4i",pix.GetPixId())
1423 << endl;
1424 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1425 if (IsCheckDeviatingBehavior())
1426 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1427 }
1428 }
1429 }
1430
1431 return kTRUE;
1432}
1433
1434
1435
1436// ------------------------------------------------------------------------
1437//
1438// Returns kFALSE if pointer to MCalibrationBlindCam is NULL
1439//
1440// The check returns kFALSE if:
1441//
1442// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1443// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1444//
1445// Calls:
1446// - MCalibrationBlindPix::CalcFluxInsidePlexiglass()
1447//
1448Bool_t MCalibrationChargeCalc::FinalizeBlindCam()
1449{
1450
1451 MCalibrationBlindCam *blindcam = fIntensBlind
1452 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
1453
1454 if (!blindcam)
1455 return kFALSE;
1456
1457 Int_t nvalid = 0;
1458
1459 for (Int_t i=0; i<blindcam->GetSize(); i++)
1460 {
1461
1462 MCalibrationBlindPix &blindpix = (MCalibrationBlindPix&)(*blindcam)[i];
1463
1464 if (!blindpix.IsValid())
1465 continue;
1466
1467 const Float_t lambda = blindpix.GetLambda();
1468 const Float_t lambdaerr = blindpix.GetLambdaErr();
1469 const Float_t lambdacheck = blindpix.GetLambdaCheck();
1470
1471 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1472 {
1473 *fLog << warn << GetDescriptor()
1474 << Form("%s%4.2f%s%4.2f%s%4.2f%s%2i",": Lambda: ",lambda," and Lambda-Check: ",
1475 lambdacheck," differ by more than ",fLambdaCheckLimit," in the Blind Pixel Nr.",i)
1476 << endl;
1477 blindpix.SetValid(kFALSE);
1478 continue;
1479 }
1480
1481 if (lambdaerr > fLambdaErrLimit)
1482 {
1483 *fLog << warn << GetDescriptor()
1484 << Form("%s%4.2f%s%4.2f%s%2i",": Error of Fitted Lambda: ",lambdaerr," is greater than ",
1485 fLambdaErrLimit," in Blind Pixel Nr.",i) << endl;
1486 blindpix.SetValid(kFALSE);
1487 continue;
1488 }
1489
1490 if (!blindpix.CalcFluxInsidePlexiglass())
1491 {
1492 *fLog << warn << "Could not calculate the flux of photons from Blind Pixel Nr." << i << endl;
1493 blindpix.SetValid(kFALSE);
1494 continue;
1495 }
1496
1497 nvalid++;
1498 }
1499
1500 if (!nvalid)
1501 return kFALSE;
1502
1503 return kTRUE;
1504}
1505
1506// ------------------------------------------------------------------------
1507//
1508// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1509//
1510// The check returns kFALSE if:
1511//
1512// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1513// 2) PINDiode has a fit error smaller than fChargeErrLimit
1514// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1515// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1516//
1517// Calls:
1518// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1519//
1520Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1521{
1522
1523 if (!fPINDiode)
1524 return kFALSE;
1525
1526 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1527 {
1528 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1529 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1530 return kFALSE;
1531 }
1532
1533 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1534 {
1535 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1536 << fChargeErrLimit << " in PINDiode " << endl;
1537 return kFALSE;
1538 }
1539
1540 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1541 {
1542 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1543 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1544 return kFALSE;
1545 }
1546
1547 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1548 {
1549 *fLog << warn << GetDescriptor()
1550 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1551 return kFALSE;
1552 }
1553
1554
1555 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1556 {
1557 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1558 << "will skip PIN Diode Calibration " << endl;
1559 return kFALSE;
1560 }
1561
1562 return kTRUE;
1563}
1564
1565// ------------------------------------------------------------------------
1566//
1567// Calculate the average number of photons outside the plexiglass with the
1568// formula:
1569//
1570// av.Num.photons(area index) = av.Num.Phes(area index)
1571// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1572// / MCalibrationQEPix::GetPMTCollectionEff()
1573// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1574// / MCalibrationQECam::GetPlexiglassQE()
1575//
1576// Calculate the variance on the average number of photons assuming that the error on the
1577// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1578// values keeps it ordinary error since it is systematic.
1579//
1580// Loop over pixels:
1581//
1582// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1583// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1584//
1585// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1586// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1587//
1588// - Calculate the quantum efficiency with the formula:
1589//
1590// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1591//
1592// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1593//
1594// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1595// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1596//
1597// - Call MCalibrationQEPix::UpdateFFactorMethod()
1598//
1599void MCalibrationChargeCalc::FinalizeFFactorQECam()
1600{
1601
1602 if (fNumInnerFFactorMethodUsed < 2)
1603 {
1604 *fLog << warn << GetDescriptor()
1605 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1606 return;
1607 }
1608
1609 MCalibrationQECam *qecam = fIntensQE
1610 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1611 MCalibrationChargeCam *chargecam = fIntensCam
1612 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1613 MBadPixelsCam *badcam = fIntensBad
1614 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1615
1616 MCalibrationChargePix &avpix = (MCalibrationChargePix&)chargecam->GetAverageArea(0);
1617 MCalibrationQEPix &qepix = (MCalibrationQEPix&) qecam->GetAverageArea(0);
1618
1619 const Float_t avphotons = avpix.GetPheFFactorMethod()
1620 / qepix.GetDefaultQE(fPulserColor)
1621 / qepix.GetPMTCollectionEff()
1622 / qepix.GetLightGuidesEff(fPulserColor)
1623 / qecam->GetPlexiglassQE();
1624
1625 const Float_t avphotrelvar = avpix.GetPheFFactorMethodRelVar()
1626 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1627 + qepix.GetPMTCollectionEffRelVar()
1628 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1629 + qecam->GetPlexiglassQERelVar();
1630
1631 const UInt_t nareas = fGeom->GetNumAreas();
1632
1633 //
1634 // Set the results in the MCalibrationChargeCam
1635 //
1636 chargecam->SetNumPhotonsFFactorMethod (avphotons);
1637
1638 if (avphotrelvar > 0.)
1639 chargecam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1640
1641 TArrayF lowlim (nareas);
1642 TArrayF upplim (nareas);
1643 TArrayD avffactorphotons (nareas);
1644 TArrayD avffactorphotvar (nareas);
1645 TArrayI numffactor (nareas);
1646
1647 const UInt_t npixels = fGeom->GetNumPixels();
1648
1649 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1650
1651 for (UInt_t i=0; i<npixels; i++)
1652 {
1653
1654 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1655 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1656 MBadPixelsPix &bad = (*badcam) [i];
1657
1658 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1659 continue;
1660
1661 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1662 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1663
1664 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1665
1666 qepix.SetQEFFactor ( qe , fPulserColor );
1667 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1668 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1669
1670 if (!qepix.UpdateFFactorMethod())
1671 *fLog << warn << GetDescriptor()
1672 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1673
1674 //
1675 // The following pixels are those with deviating sigma, but otherwise OK,
1676 // probably those with stars during the pedestal run, but not the cal. run.
1677 //
1678 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1679 {
1680 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1681 if (IsCheckDeviatingBehavior())
1682 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1683 continue;
1684 }
1685
1686 const Int_t aidx = (*fGeom)[i].GetAidx();
1687 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1688
1689 camffactor.Fill(i,ffactor);
1690 camffactor.SetUsed(i);
1691
1692 avffactorphotons[aidx] += ffactor;
1693 avffactorphotvar[aidx] += ffactor*ffactor;
1694 numffactor[aidx]++;
1695 }
1696
1697 for (UInt_t i=0; i<nareas; i++)
1698 {
1699
1700 if (numffactor[i] == 0)
1701 {
1702 *fLog << warn << GetDescriptor() << ": No pixels with valid total F-Factor found "
1703 << "in area index: " << i << endl;
1704 continue;
1705 }
1706
1707 avffactorphotvar[i] = (avffactorphotvar[i] - avffactorphotons[i]*avffactorphotons[i]/numffactor[i])
1708 / (numffactor[i]-1.);
1709 avffactorphotons[i] = avffactorphotons[i] / numffactor[i];
1710
1711 if (avffactorphotvar[i] < 0.)
1712 {
1713 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of total F-Factor found "
1714 << "in area index: " << i << endl;
1715 continue;
1716 }
1717
1718 lowlim [i] = 1.1; // Lowest known F-Factor of a PMT
1719 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1720
1721 TArrayI area(1);
1722 area[0] = i;
1723
1724 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1725 hist->Fit("gaus","Q");
1726 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1727 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1728 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1729
1730 if (IsDebug())
1731 camffactor.DrawClone();
1732
1733 if (ndf < 2)
1734 {
1735 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1736 << "in the camera with area index: " << i << endl;
1737 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1738 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1739 delete hist;
1740 continue;
1741 }
1742
1743 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1744
1745 if (prob < 0.001)
1746 {
1747 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1748 << "in the camera with area index: " << i << endl;
1749 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1750 << " is smaller than 0.001 " << endl;
1751 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1752 delete hist;
1753 continue;
1754 }
1755
1756 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1757 << "with area index #" << i << ": "
1758 << Form("%4.2f+-%4.2f",mean,sigma) << endl;
1759
1760 lowlim [i] = 1.1;
1761 upplim [i] = mean + fFFactorErrLimit*sigma;
1762
1763 delete hist;
1764 }
1765
1766 *fLog << endl;
1767
1768 for (UInt_t i=0; i<npixels; i++)
1769 {
1770
1771 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1772 MBadPixelsPix &bad = (*badcam) [i];
1773
1774 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1775 continue;
1776
1777 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1778 const Int_t aidx = (*fGeom)[i].GetAidx();
1779
1780 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1781 {
1782 *fLog << warn << GetDescriptor() << ": Overall F-Factor "
1783 << Form("%5.2f",ffactor) << " out of range ["
1784 << Form("%5.2f,%5.2f",lowlim[aidx],upplim[aidx]) << "] pixel " << i << endl;
1785
1786 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1787 if (IsCheckDeviatingBehavior())
1788 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1789 }
1790 }
1791
1792 for (UInt_t i=0; i<npixels; i++)
1793 {
1794
1795 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1796 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1797 MBadPixelsPix &bad = (*badcam) [i];
1798
1799 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1800 {
1801 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1802 pix.SetFFactorMethodValid(kFALSE);
1803 pix.SetExcluded();
1804 continue;
1805 }
1806 }
1807}
1808
1809
1810// ------------------------------------------------------------------------
1811//
1812// Loop over pixels:
1813//
1814// - Continue, if not MCalibrationBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1815// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1816//
1817// - Calculate the quantum efficiency with the formula:
1818//
1819// QE = Num.Phes / MCalibrationBlindPix::GetFluxInsidePlexiglass()
1820// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1821//
1822// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1823// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1824// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1825//
1826// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1827//
1828void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1829{
1830
1831
1832 MCalibrationBlindCam *blindcam = fIntensBlind
1833 ? (MCalibrationBlindCam*) fIntensBlind->GetCam(): fBlindCam;
1834 MBadPixelsCam *badcam = fIntensBad
1835 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1836 MCalibrationQECam *qecam = fIntensQE
1837 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1838 MCalibrationChargeCam *chargecam = fIntensCam
1839 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1840
1841 //
1842 // Set the results in the MCalibrationChargeCam
1843 //
1844 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
1845 {
1846
1847 const Float_t photons = blindcam->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
1848 / qecam->GetPlexiglassQE();
1849 chargecam->SetNumPhotonsBlindPixelMethod(photons);
1850
1851 const Float_t photrelvar = blindcam->GetFluxInsidePlexiglassRelVar()
1852 + qecam->GetPlexiglassQERelVar();
1853
1854 if (photrelvar > 0.)
1855 chargecam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
1856 }
1857
1858 //
1859 // With the knowledge of the overall photon flux, calculate the
1860 // quantum efficiencies after the Blind Pixel and PIN Diode method
1861 //
1862 const UInt_t npixels = fGeom->GetNumPixels();
1863 for (UInt_t i=0; i<npixels; i++)
1864 {
1865
1866 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1867
1868 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
1869 {
1870 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1871 continue;
1872 }
1873
1874 MBadPixelsPix &bad = (*badcam) [i];
1875
1876 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1877 {
1878 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1879 continue;
1880 }
1881
1882 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1883 MGeomPix &geo = (*fGeom) [i];
1884
1885 const Float_t qe = pix.GetPheFFactorMethod()
1886 / blindcam->GetFluxInsidePlexiglass()
1887 / geo.GetA()
1888 * qecam->GetPlexiglassQE();
1889
1890 const Float_t qerelvar = blindcam->GetFluxInsidePlexiglassRelVar()
1891 + qecam->GetPlexiglassQERelVar()
1892 + pix.GetPheFFactorMethodRelVar();
1893
1894 qepix.SetQEBlindPixel ( qe , fPulserColor );
1895 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
1896 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
1897
1898 if (!qepix.UpdateBlindPixelMethod())
1899 *fLog << warn << GetDescriptor()
1900 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
1901 }
1902}
1903
1904// ------------------------------------------------------------------------
1905//
1906// Loop over pixels:
1907//
1908// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
1909// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
1910//
1911// - Calculate the quantum efficiency with the formula:
1912//
1913// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
1914//
1915// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
1916// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
1917// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
1918//
1919// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
1920//
1921void MCalibrationChargeCalc::FinalizePINDiodeQECam()
1922{
1923
1924 const UInt_t npixels = fGeom->GetNumPixels();
1925
1926 MCalibrationQECam *qecam = fIntensQE
1927 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1928 MCalibrationChargeCam *chargecam = fIntensCam
1929 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1930 MBadPixelsCam *badcam = fIntensBad
1931 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1932 //
1933 // With the knowledge of the overall photon flux, calculate the
1934 // quantum efficiencies after the PIN Diode method
1935 //
1936 for (UInt_t i=0; i<npixels; i++)
1937 {
1938
1939 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1940
1941 if (!fPINDiode)
1942 {
1943 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1944 continue;
1945 }
1946
1947 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
1948 {
1949 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1950 continue;
1951 }
1952
1953 MBadPixelsPix &bad = (*badcam) [i];
1954
1955 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1956 {
1957 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1958 continue;
1959 }
1960
1961 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1962 MGeomPix &geo = (*fGeom) [i];
1963
1964 const Float_t qe = pix.GetPheFFactorMethod()
1965 / fPINDiode->GetFluxOutsidePlexiglass()
1966 / geo.GetA();
1967
1968 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
1969
1970 qepix.SetQEPINDiode ( qe , fPulserColor );
1971 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
1972 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
1973
1974 if (!qepix.UpdatePINDiodeMethod())
1975 *fLog << warn << GetDescriptor()
1976 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
1977 }
1978}
1979
1980// ------------------------------------------------------------------------
1981//
1982// Loop over pixels:
1983//
1984// - Call MCalibrationQEPix::UpdateCombinedMethod()
1985//
1986void MCalibrationChargeCalc::FinalizeCombinedQECam()
1987{
1988
1989 const UInt_t npixels = fGeom->GetNumPixels();
1990
1991 MCalibrationQECam *qecam = fIntensQE
1992 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1993 MBadPixelsCam *badcam = fIntensBad
1994 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1995
1996 for (UInt_t i=0; i<npixels; i++)
1997 {
1998
1999 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2000 MBadPixelsPix &bad = (*badcam) [i];
2001
2002 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2003 {
2004 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2005 continue;
2006 }
2007
2008 qepix.UpdateCombinedMethod();
2009 }
2010}
2011
2012// -----------------------------------------------------------------------------------------------
2013//
2014// - Print out statistics about BadPixels of type UnsuitableType_t
2015// - store numbers of bad pixels of each type in fCam or fIntensCam
2016//
2017void MCalibrationChargeCalc::FinalizeUnsuitablePixels()
2018{
2019
2020 *fLog << inf << endl;
2021 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
2022 *fLog << dec << setfill(' ');
2023
2024 const Int_t nareas = fGeom->GetNumAreas();
2025
2026 TArrayI counts(nareas);
2027
2028 MBadPixelsCam *badcam = fIntensBad
2029 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2030 MCalibrationChargeCam *chargecam = fIntensCam
2031 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
2032
2033 for (Int_t i=0; i<badcam->GetSize(); i++)
2034 {
2035 MBadPixelsPix &bad = (*badcam)[i];
2036 if (!bad.IsBad())
2037 {
2038 const Int_t aidx = (*fGeom)[i].GetAidx();
2039 counts[aidx]++;
2040 }
2041 }
2042
2043 if (fGeom->InheritsFrom("MGeomCamMagic"))
2044 *fLog << " " << setw(7) << "Successfully calibrated Pixels: "
2045 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2046
2047 counts.Reset();
2048
2049 for (Int_t i=0; i<badcam->GetSize(); i++)
2050 {
2051 MBadPixelsPix &bad = (*badcam)[i];
2052
2053 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
2054 {
2055 const Int_t aidx = (*fGeom)[i].GetAidx();
2056 counts[aidx]++;
2057 }
2058 }
2059
2060 for (Int_t aidx=0; aidx<nareas; aidx++)
2061 chargecam->SetNumUnsuitable(counts[aidx], aidx);
2062
2063 if (fGeom->InheritsFrom("MGeomCamMagic"))
2064 *fLog << " " << setw(7) << "Uncalibrated Pixels: "
2065 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2066
2067 counts.Reset();
2068
2069 for (Int_t i=0; i<badcam->GetSize(); i++)
2070 {
2071
2072 MBadPixelsPix &bad = (*badcam)[i];
2073
2074 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
2075 {
2076 const Int_t aidx = (*fGeom)[i].GetAidx();
2077 counts[aidx]++;
2078 }
2079 }
2080
2081 for (Int_t aidx=0; aidx<nareas; aidx++)
2082 chargecam->SetNumUnreliable(counts[aidx], aidx);
2083
2084 *fLog << " " << setw(7) << "Unreliable Pixels: "
2085 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2086
2087}
2088
2089// -----------------------------------------------------------------------------------------------
2090//
2091// Print out statistics about BadPixels of type UncalibratedType_t
2092//
2093void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
2094{
2095
2096 UInt_t countinner = 0;
2097 UInt_t countouter = 0;
2098
2099 MBadPixelsCam *badcam = fIntensBad
2100 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2101
2102 for (Int_t i=0; i<badcam->GetSize(); i++)
2103 {
2104 MBadPixelsPix &bad = (*badcam)[i];
2105
2106 if (bad.IsUncalibrated(typ))
2107 {
2108 if (fGeom->GetPixRatio(i) == 1.)
2109 countinner++;
2110 else
2111 countouter++;
2112 }
2113 }
2114
2115 *fLog << " " << setw(7) << text
2116 << Form("%s%3i%s%3i","Inner: ",countinner," Outer: ",countouter) << endl;
2117}
2118
2119// --------------------------------------------------------------------------
2120//
2121// Set the path for output file
2122//
2123void MCalibrationChargeCalc::SetOutputPath(TString path)
2124{
2125 fOutputPath = path;
2126 if (fOutputPath.EndsWith("/"))
2127 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
2128}
2129
2130// --------------------------------------------------------------------------
2131//
2132// Set the output file
2133//
2134void MCalibrationChargeCalc::SetOutputFile(TString file)
2135{
2136 fOutputFile = file;
2137}
2138
2139// --------------------------------------------------------------------------
2140//
2141// Get the output file
2142//
2143const char* MCalibrationChargeCalc::GetOutputFile()
2144{
2145 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
2146}
2147
2148// --------------------------------------------------------------------------
2149//
2150// Read the environment for the following data members:
2151// - fChargeLimit
2152// - fChargeErrLimit
2153// - fChargeRelErrLimit
2154// - fDebug
2155// - fFFactorErrLimit
2156// - fLambdaErrLimit
2157// - fLambdaCheckErrLimit
2158// - fPheErrLimit
2159//
2160Int_t MCalibrationChargeCalc::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
2161{
2162
2163 Bool_t rc = kFALSE;
2164 if (IsEnvDefined(env, prefix, "ChargeLimit", print))
2165 {
2166 SetChargeLimit(GetEnvValue(env, prefix, "ChargeLimit", fChargeLimit));
2167 rc = kTRUE;
2168 }
2169 if (IsEnvDefined(env, prefix, "ChargeErrLimit", print))
2170 {
2171 SetChargeErrLimit(GetEnvValue(env, prefix, "ChargeErrLimit", fChargeErrLimit));
2172 rc = kTRUE;
2173 }
2174 if (IsEnvDefined(env, prefix, "ChargeRelErrLimit", print))
2175 {
2176 SetChargeRelErrLimit(GetEnvValue(env, prefix, "ChargeRelErrLimit", fChargeRelErrLimit));
2177 rc = kTRUE;
2178 }
2179 if (IsEnvDefined(env, prefix, "Debug", print))
2180 {
2181 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
2182 rc = kTRUE;
2183 }
2184 if (IsEnvDefined(env, prefix, "FFactorErrLimit", print))
2185 {
2186 SetFFactorErrLimit(GetEnvValue(env, prefix, "FFactorErrLimit", fFFactorErrLimit));
2187 rc = kTRUE;
2188 }
2189 if (IsEnvDefined(env, prefix, "LambdaErrLimit", print))
2190 {
2191 SetLambdaErrLimit(GetEnvValue(env, prefix, "LambdaErrLimit", fLambdaErrLimit));
2192 rc = kTRUE;
2193 }
2194 if (IsEnvDefined(env, prefix, "LambdaCheckLimit", print))
2195 {
2196 SetLambdaCheckLimit(GetEnvValue(env, prefix, "LambdaCheckLimit", fLambdaCheckLimit));
2197 rc = kTRUE;
2198 }
2199 if (IsEnvDefined(env, prefix, "PheErrLimit", print))
2200 {
2201 SetPheErrLimit(GetEnvValue(env, prefix, "PheErrLimit", fPheErrLimit));
2202 rc = kTRUE;
2203 }
2204 if (IsEnvDefined(env, prefix, "CheckDeadPixels", print))
2205 {
2206 SetCheckDeadPixels(GetEnvValue(env, prefix, "CheckDeadPixels", IsCheckDeadPixels()));
2207 rc = kTRUE;
2208 }
2209 if (IsEnvDefined(env, prefix, "CheckDeviatingBehavior", print))
2210 {
2211 SetCheckDeviatingBehavior(GetEnvValue(env, prefix, "CheckDeviatingBehavior", IsCheckDeviatingBehavior()));
2212 rc = kTRUE;
2213 }
2214 if (IsEnvDefined(env, prefix, "CheckExtractionWindow", print))
2215 {
2216 SetCheckExtractionWindow(GetEnvValue(env, prefix, "CheckExtractionWindow", IsCheckExtractionWindow()));
2217 rc = kTRUE;
2218 }
2219 if (IsEnvDefined(env, prefix, "CheckHistOverflow", print))
2220 {
2221 SetCheckHistOverflow(GetEnvValue(env, prefix, "CheckHistOverflow", IsCheckHistOverflow()));
2222 rc = kTRUE;
2223 }
2224 if (IsEnvDefined(env, prefix, "CheckOscillations", print))
2225 {
2226 SetCheckOscillations(GetEnvValue(env, prefix, "CheckOscillations", IsCheckOscillations()));
2227 rc = kTRUE;
2228 }
2229
2230 return rc;
2231}
Note: See TracBrowser for help on using the repository browser.