source: branches/Mars_McMismatchStudy/mcalib/MCalibrationChargeCalc.cc@ 18222

Last change on this file since 18222 was 10166, checked in by tbretz, 14 years ago
Removed the old obsolete cvs header line.
File size: 79.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! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19! Author(s): Thomas Bretz <mailto:tbretz@astro.uni-wuerzburg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2007
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// ClassVersion 2:
183// - Float_t fPheErrLimit;
184// + Float_t fPheErrLowerLimit; // Lower limit acceptance nr. phe's w.r.t. area idx mean (in sigmas)
185// + Float_t fPheErrUpperLimit; // Upper limit acceptance nr. phe's w.r.t. area idx mean (in sigmas)
186//
187//
188// ClassVersion 3:
189// + Bool_t fUseExtractorRes; // Include extractor resolution in F-Factor method
190//
191// ClassVersion 4:
192// + Float_t fUnsuitablesLimit; // Limit for relative number of unsuitable pixels
193// + Float_t fUnreliablesLimit; // Limit for relative number of unreliable pixels
194// + Float_t fExternalNumPhes; // External mean number of photo-electrons set from outside
195// + Float_t fExternalNumPhesRelVar; // External rel. var. number of photo-electrons set from outside
196//
197// ClassVersion 5:
198// - TString fOutputPath
199// - TString fOutputFile
200//
201// ClassVersion 6:
202// - Bool_t fContinousCalibration
203//
204//////////////////////////////////////////////////////////////////////////////
205#include "MCalibrationChargeCalc.h"
206
207#include <TSystem.h>
208#include <TH1.h>
209#include <TF1.h>
210
211#include "MLog.h"
212#include "MLogManip.h"
213
214#include "MMath.h"
215
216#include "MParList.h"
217
218#include "MCalibrationPattern.h"
219
220#include "MGeomCam.h"
221#include "MGeom.h"
222#include "MHCamera.h"
223
224#include "MPedestalCam.h"
225#include "MPedestalPix.h"
226
227#include "MHCalibrationChargeCam.h"
228#include "MHCalibrationChargeBlindCam.h"
229
230#include "MCalibrationChargeCam.h"
231#include "MCalibrationChargePix.h"
232#include "MCalibrationChargePINDiode.h"
233#include "MCalibrationBlindPix.h"
234#include "MCalibrationBlindCam.h"
235
236#include "MExtractedSignalCam.h"
237#include "MExtractedSignalPix.h"
238#include "MExtractedSignalBlindPixel.h"
239#include "MExtractedSignalPINDiode.h"
240
241#include "MBadPixelsCam.h"
242
243#include "MCalibrationQECam.h"
244#include "MCalibrationQEPix.h"
245
246ClassImp(MCalibrationChargeCalc);
247
248using namespace std;
249
250const Float_t MCalibrationChargeCalc::fgChargeLimit = 4.5;
251const Float_t MCalibrationChargeCalc::fgChargeErrLimit = 0.;
252const Float_t MCalibrationChargeCalc::fgChargeRelErrLimit = 1.;
253const Float_t MCalibrationChargeCalc::fgLambdaErrLimit = 0.2;
254const Float_t MCalibrationChargeCalc::fgLambdaCheckLimit = 0.5;
255const Float_t MCalibrationChargeCalc::fgPheErrLowerLimit = 6.0;
256const Float_t MCalibrationChargeCalc::fgPheErrUpperLimit = 5.5;
257const Float_t MCalibrationChargeCalc::fgFFactorErrLimit = 4.5;
258const Float_t MCalibrationChargeCalc::fgArrTimeRmsLimit = 5.0;
259const Float_t MCalibrationChargeCalc::fgUnsuitablesLimit = 0.1;
260const Float_t MCalibrationChargeCalc::fgUnreliablesLimit = 0.3;
261
262const char *MCalibrationChargeCalc::fgNamePedestalCam = "MPedestalCam";
263
264
265// --------------------------------------------------------------------------
266//
267// Default constructor.
268//
269// Sets the pointer to fQECam and fGeom to NULL
270//
271// Calls AddToBranchList for:
272// - MRawEvtData.fHiGainPixId
273// - MRawEvtData.fLoGainPixId
274// - MRawEvtData.fHiGainFadcSamples
275// - MRawEvtData.fLoGainFadcSamples
276//
277// Initializes:
278// - fArrTimeRmsLimit to fgArrTimeRmsLimit
279// - fChargeLimit to fgChargeLimit
280// - fChargeErrLimit to fgChargeErrLimit
281// - fChargeRelErrLimit to fgChargeRelErrLimit
282// - fFFactorErrLimit to fgFFactorErrLimit
283// - fLambdaCheckLimit to fgLambdaCheckLimit
284// - fLambdaErrLimit to fgLambdaErrLimit
285// - fNamePedestalCam to fgNamePedestalCam
286// - fPheErrLowerLimit to fgPheErrLowerLimit
287// - fPheErrUpperLimit to fgPheErrUpperLimit
288// - fPulserColor to MCalibrationCam::kCT1
289// - fOutputPath to "."
290// - fOutputFile to "ChargeCalibStat.txt"
291// - flag debug to kFALSE
292// - IsUseExtractorRes to kFALSE
293// - IsUseExternalNumPhes to kFALSE
294// - fExternalNumPhes to 0.
295// - fExternalNumPhesRelVar to 0.
296//
297// Sets all checks
298//
299// Calls:
300// - Clear()
301//
302MCalibrationChargeCalc::MCalibrationChargeCalc(const char *name, const char *title)
303 : fPulserColor(MCalibrationCam::kNONE), fContinousCalibration(kFALSE),
304 fGeom(NULL), fSignal(NULL), fCalibPattern(NULL), fExtractor(NULL)
305{
306
307 fName = name ? name : "MCalibrationChargeCalc";
308 fTitle = title ? title : "Task to calculate the calibration constants and MCalibrationCam ";
309
310 AddToBranchList("MRawEvtData.fHiGainPixId");
311 AddToBranchList("MRawEvtData.fLoGainPixId");
312 AddToBranchList("MRawEvtData.fHiGainFadcSamples");
313 AddToBranchList("MRawEvtData.fLoGainFadcSamples");
314
315 SetArrTimeRmsLimit ();
316 SetChargeLimit ();
317 SetChargeErrLimit ();
318 SetChargeRelErrLimit ();
319 SetDebug ( kFALSE );
320 SetFFactorErrLimit ();
321 SetLambdaCheckLimit ();
322 SetLambdaErrLimit ();
323 SetNamePedestalCam ();
324 SetPheErrLowerLimit ();
325 SetPheErrUpperLimit ();
326 SetUnsuitablesLimit ();
327 SetUnreliablesLimit ();
328 SetUseExtractorRes ();
329 SetUseUnreliables ();
330 SetUseExternalNumPhes(kFALSE);
331
332 SetExternalNumPhes ();
333 SetExternalNumPhesRelVar();
334
335 SetCheckArrivalTimes ();
336 SetCheckDeadPixels ();
337 SetCheckDeviatingBehavior();
338 SetCheckExtractionWindow ();
339 SetCheckHistOverflow ();
340 SetCheckOscillations ();
341
342 Clear();
343
344}
345
346// --------------------------------------------------------------------------
347//
348// Sets:
349// - all variables to 0.,
350// - all flags to kFALSE
351// - all pointers to NULL
352// - the pulser colour to kNONE
353// - fBlindPixelFlags to 0
354// - fPINDiodeFlags to 0
355//
356void MCalibrationChargeCalc::Clear(const Option_t *o)
357{
358
359 fNumHiGainSamples = 0.;
360 fNumLoGainSamples = 0.;
361 fSqrtHiGainSamples = 0.;
362 fSqrtLoGainSamples = 0.;
363 fNumInnerFFactorMethodUsed = 0;
364
365 fNumProcessed = 0;
366
367 fBadPixels = NULL;
368 fCam = NULL;
369 fHCam = NULL;
370 fQECam = NULL;
371 fBlindCam = NULL;
372 fHBlindCam = NULL;
373 fPINDiode = NULL;
374 fPedestals = NULL;
375
376 SetPulserColor ( MCalibrationCam::kNONE );
377
378 fStrength = 0.;
379 fBlindPixelFlags.Set(0);
380 fPINDiodeFlags .Set(0);
381 fResultFlags .Set(0);
382}
383
384
385// -----------------------------------------------------------------------------------
386//
387// The following container are searched for and execution aborted if not in MParList:
388// - MPedestalCam
389// - MCalibrationPattern
390// - MExtractedSignalCam
391//
392Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
393{
394
395 /*
396 if (IsInterlaced())
397 {
398 fTrigPattern = (MTriggerPattern*)pList->FindObject("MTriggerPattern");
399 if (!fTrigPattern)
400 {
401 *fLog << err << "MTriggerPattern not found... abort." << endl;
402 return kFALSE;
403 }
404 }
405 */
406
407 fCalibPattern = (MCalibrationPattern*)pList->FindObject("MCalibrationPattern");
408 if (!fCalibPattern)
409 {
410 *fLog << err << "MCalibrationPattern not found... abort." << endl;
411 return kFALSE;
412 }
413
414 //
415 // Containers that have to be there.
416 //
417 fSignal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
418 if (!fSignal)
419 {
420 *fLog << err << "MExtractedSignalCam not found... aborting" << endl;
421 return kFALSE;
422 }
423
424 if (fPedestals)
425 return kTRUE;
426
427 fPedestals = (MPedestalCam*)pList->FindObject(AddSerialNumber(fNamePedestalCam), "MPedestalCam");
428 if (!fPedestals)
429 {
430 *fLog << err << fNamePedestalCam << " [MPedestalCam] not found... aborting" << endl;
431 return kFALSE;
432 }
433
434 // Obsolete... it is already set in the constructor isn't it?
435 fPulserColor = MCalibrationCam::kNONE;
436
437 return kTRUE;
438}
439
440
441// --------------------------------------------------------------------------
442//
443// Search for the following input containers and abort if not existing:
444// - MGeomCam
445// - MCalibrationChargeCam
446// - MCalibrationQECam
447// - MBadPixelsCam
448//
449// Search for the following input containers and give a warning if not existing:
450// - MCalibrationBlindPix
451// - MCalibrationChargePINDiode
452//
453// It retrieves the following variables from MCalibrationChargeCam:
454//
455// - fNumHiGainSamples
456// - fNumLoGainSamples
457//
458// It defines the PixId of every pixel in:
459//
460// - MCalibrationChargeCam
461// - MCalibrationQECam
462//
463// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
464//
465// - MCalibrationChargePix
466// - MCalibrationQEPix
467//
468// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
469//
470// - MCalibrationChargeCam
471// - MCalibrationBlindPix (if existing)
472// - MCalibrationChargePINDiode (if existing)
473//
474Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
475{
476
477 //
478 // Geometry
479 //
480 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
481 if (!fGeom)
482 {
483 *fLog << err << "No MGeomCam found... aborting." << endl;
484 return kFALSE;
485 }
486
487 //
488 // Charge histogram
489 //
490 fHCam = (MHCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeCam"));
491 if (!fHCam)
492 {
493 *fLog << err << "Cannot find MHCalibrationChargeCam ... abort." << endl;
494 *fLog << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
495 return kFALSE;
496 }
497
498 //
499 // MCalibrationChargeCam
500 //
501 fCam = (MCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
502 if (!fCam)
503 {
504 *fLog << err << "Cannot find MCalibrationChargeCam ... abort." << endl;
505 *fLog << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
506 return kFALSE;
507 }
508
509 //
510 // Blind pixel calibration
511 //
512 fBlindCam = (MCalibrationBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationBlindCam"));
513 if (!fBlindCam)
514 *fLog << inf << "No MCalibrationBlindCam found... no Blind Pixel method!" << endl;
515
516 //
517 // Required containers
518 //
519 fQECam = (MCalibrationQECam*)pList->FindObject(AddSerialNumber("MCalibrationQECam"));
520 if (!fQECam)
521 {
522 *fLog << err << "Cannot find MCalibrationQECam ... abort." << endl;
523 return kFALSE;
524 }
525
526 fBadPixels = (MBadPixelsCam*)pList->FindObject(AddSerialNumber("MBadPixelsCam"));
527 if (!fBadPixels)
528 {
529 *fLog << err << "Cannot find MBadPixelsCam ... abort." << endl;
530 return kFALSE;
531 }
532
533 //
534 // Optional Containers
535 //
536 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
537 if (!fPINDiode)
538 *fLog << inf << "No MCalibrationChargePINDiode found... no PIN Diode method!" << endl;
539
540 fHBlindCam = (MHCalibrationChargeBlindCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeBlindCam"));
541 if (!fHBlindCam)
542 *fLog << inf << "No MHCalibrationChargeBlindCam found... no Blind Pixel method!" << endl;
543
544
545 //
546 // Further initialization
547 //
548 if (IsDebug())
549 {
550 const UInt_t npixels = fGeom->GetNumPixels();
551 for (UInt_t i=0; i<npixels; i++)
552 (*fCam)[i].SetDebug();
553 }
554
555 fResultFlags.Set(fGeom->GetNumAreas());
556
557 return kTRUE;
558}
559
560// ----------------------------------------------------------------------------------
561//
562// Set the correct colour to the charge containers
563//
564Int_t MCalibrationChargeCalc::Process()
565{
566
567 const MCalibrationCam::PulserColor_t col = fCalibPattern->GetPulserColor();
568
569 //
570 // If pulser color didn't change there is nothing to do
571 //
572 const Float_t strength = fCalibPattern->GetPulserStrength();
573 const Float_t strdiff = TMath::Abs(strength-fStrength);
574 if (col == fPulserColor && strdiff < 0.05 )
575 {
576 fNumProcessed++;
577 return kTRUE;
578 }
579
580 //
581 // If there is no color is nothing to do
582 //
583 if (col == MCalibrationCam::kNONE)
584 return kTRUE;
585
586 //
587 // Now retrieve the colour and check if not various colours have been used
588 //
589 if (!fContinousCalibration && fPulserColor != MCalibrationCam::kNONE)
590 {
591 *fLog << warn << "Multiple colours used simultaneously!" << flush;
592 fHCam->Finalize();
593 if (fHBlindCam)
594 fHBlindCam->Finalize();
595
596 Finalize();
597
598 fHCam->ResetHists();
599 if (fHBlindCam)
600 fHBlindCam->ResetHists();
601
602 *fLog << inf << "Starting next calibration... " << flush;
603
604 fHCam->SetColor(col);
605 if (fHBlindCam)
606 fHBlindCam->SetColor(col);
607
608 fCam->SetPulserColor(col);
609 if (fBlindCam)
610 fBlindCam->SetPulserColor(col);
611 }
612
613 //
614 // Keep pulser color and strength in mind and output some information
615 //
616 fPulserColor = col;
617 fStrength = strength;
618
619 *fLog << inf << "Found new colour ... " << flush;
620
621 switch (col)
622 {
623 case MCalibrationCam::kGREEN: *fLog << "Green"; break;
624 case MCalibrationCam::kBLUE: *fLog << "Blue"; break;
625 case MCalibrationCam::kUV: *fLog << "UV"; break;
626 case MCalibrationCam::kCT1: *fLog << "CT1"; break;
627 default: break;
628 }
629
630 *fLog << inf << " with strength: " << strength << endl;
631
632 //
633 // Distribute new color to all containers
634 //
635 fCam->SetPulserColor(col);
636 if (fBlindCam)
637 fBlindCam->SetPulserColor(col);
638
639 fHCam->SetColor(col);
640 if (fHBlindCam)
641 fHBlindCam->SetColor(col);
642 if (fPINDiode)
643 fPINDiode->SetColor(col);
644
645 //
646 // Reset number of processed events
647 //
648 fNumProcessed = 0;
649
650 return kTRUE;
651}
652
653// -----------------------------------------------------------------------
654//
655// Return if number of executions is null.
656//
657Int_t MCalibrationChargeCalc::PostProcess()
658{
659
660 if (GetNumExecutions() < 1)
661 return kTRUE;
662
663 if (fPulserColor == MCalibrationCam::kNONE)
664 {
665 *fLog << err << "ERROR - No event was flagged with a pulser color... no calibration possible." << endl;
666 return kFALSE;
667 }
668
669 if (fNumProcessed == 0)
670 return kTRUE;
671
672 *fLog << endl;
673
674 return Finalize();
675}
676
677// ----------------------------------------------------------------------------------------------------
678//
679// Check for outliers. They are marked with
680// MBadPixelsPix::kFluctuatingArrivalTimes
681//
682// see also MCalibrationRelTimeCalc::FinalizeRelTimes
683//
684void MCalibrationChargeCalc::FinalizeAbsTimes()
685{
686 const Int_t npixels = fGeom->GetNumPixels();
687 const Int_t nareas = fGeom->GetNumAreas();
688
689 // Create an array capable of holding all pixels
690 TArrayF arr(npixels);
691
692 for (Int_t aidx=0; aidx<nareas; aidx++)
693 {
694 Int_t n = 0;
695 for (Int_t i=0; i<npixels; i++)
696 {
697 // Check for this aidx only
698 if ((*fGeom)[i].GetAidx()!=aidx)
699 continue;
700
701 // check if pixel may not contain a valid value
702 if ((*fBadPixels)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
703 continue;
704
705 // check if it was excluded for some reason
706 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
707 if (pix.IsExcluded())
708 continue;
709
710 // if TimePrecision is valid fill it into array
711 if (pix.GetAbsTimeRms()>0)
712 arr[n++] = pix.GetAbsTimeRms();
713 }
714
715 // Check the ratio of valid entries to the ratio of pixels
716 const Float_t ratio = 100*n/fGeom->GetNumPixWithAidx(aidx);
717 if (3*ratio<2)
718 *fLog << warn << "Area " << setw(4) << aidx << ": Only " << ratio << "% pixels with valid time resolution found." << endl;
719
720 // Calculate median and median deviation
721 Double_t med;
722 const Double_t dev = MMath::MedianDev(n, arr.GetArray(), med);
723
724 // Calculate upper and lower limit
725 const Float_t lolim = TMath::Max(med-fArrTimeRmsLimit*dev, 0.);
726 const Float_t hilim = TMath::Max(med+fArrTimeRmsLimit*dev, 0.);
727
728 // Now find the outliers
729 for (Int_t i=0; i<npixels; i++)
730 {
731 // Search only within this aidx
732 if ((*fGeom)[i].GetAidx()!=aidx)
733 continue;
734
735 // skip pixels already known to be unsuitable
736 if ((*fBadPixels)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
737 continue;
738
739 // check if a pixel has been excluded. This
740 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
741
742 // Check if time precision is valid (might be invalid
743 // for example in cae of empty histograms)
744 const Float_t res = pix.GetAbsTimeRms();
745 if (res<0) //FIXME!!! How does this happen?
746 {
747 *fLog << warn << "Pixel " << setw(4) << i << ": Abs-time rms could not be calculated." << endl;
748 (*fBadPixels)[i].SetUncalibrated(MBadPixelsPix::kDeviatingAbsTimeResolution);
749 continue;
750 }
751
752 if (res==0)
753 {
754 *fLog << warn << "Pixel " << setw(4) << i << ": Abs-time rms==0." << endl;
755 (*fBadPixels)[i].SetUncalibrated(MBadPixelsPix::kDeviatingAbsTimeResolution);
756 continue;
757 }
758
759 // Now compare to a lower and upper limit
760 if (res<=lolim || res>=hilim)
761 {
762 *fLog << warn << "Pixel " << setw(4) << i << ": Deviation from abs-time rms: "
763 << Form("%5.2f", res) << " out of range "
764 << Form("[%4.2f,%4.2f]", lolim, hilim) << endl;
765
766 (*fBadPixels)[i].SetUncalibrated(MBadPixelsPix::kDeviatingAbsTimeResolution);
767 }
768 }
769 }
770}
771
772// -----------------------------------------------------------------------
773//
774// Return kTRUE if fPulserColor is kNONE
775//
776// First loop over pixels, average areas and sectors, call:
777// - FinalizePedestals()
778// - FinalizeCharges()
779// for every entry. Count number of valid pixels in loop and return kFALSE
780// if there are none (the "Michele check").
781//
782// Call FinalizeBadPixels()
783//
784// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
785//
786// Call FinalizeBlindCam()
787// Call FinalizePINDiode()
788//
789// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
790// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
791// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
792//
793// Call FinalizeUnsuitablePixels()
794//
795// Call MParContainer::SetReadyToSave() for fCam, fQECam, fBadPixels and
796// fBlindCam and fPINDiode if they exist
797//
798// Print out some statistics
799//
800Int_t MCalibrationChargeCalc::Finalize()
801{
802 // The number of used slices are just a mean value
803 // the real number might change from event to event.
804 // (up to 50%!)
805
806 fNumHiGainSamples = fSignal->GetNumUsedHiGainFADCSlices();
807 fNumLoGainSamples = fSignal->GetNumUsedLoGainFADCSlices();
808
809 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
810 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
811
812 if (fPINDiode)
813 if (!fPINDiode->IsValid())
814 {
815 *fLog << warn << GetDescriptor()
816 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
817 fPINDiode = NULL;
818 }
819
820 //
821 // First loop over pixels, call FinalizePedestals and FinalizeCharges
822 //
823 Int_t nvalid = 0;
824 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
825 {
826 //
827 // Check if the pixel has been excluded from the fits
828 //
829 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[pixid];
830 if (pix.IsExcluded())
831 continue;
832
833 FinalizePedestals((*fPedestals)[pixid], pix, (*fGeom)[pixid].GetAidx());
834
835 if (FinalizeCharges(pix, (*fBadPixels)[pixid], "Pixel "))
836 nvalid++;
837 }
838
839 FinalizeAbsTimes();
840
841 *fLog << endl;
842
843 //
844 // The Michele check ...
845 //
846 if (nvalid == 0)
847 {
848 if (!fContinousCalibration)
849 {
850 *fLog << warn << GetDescriptor() << ": All pixels have non-valid calibration. "
851 << "Did you forget to fill the histograms "
852 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
853 *fLog << warn << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
854 << "instead of a calibration run " << endl;
855 return kFALSE;
856 }
857 }
858
859 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
860 {
861
862 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
863 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
864 const MArrayI &arr = fHCam->GetAverageAreaNum();
865
866 FinalizePedestals(ped,pix,aidx);
867
868 //
869 // Correct for sqrt(number of valid pixels) in the pedestal RMS
870 // (already done for calibration sigma in MHCalibrationCam::CalcAverageSigma()
871 //
872 const Double_t sqrtnum = TMath::Sqrt((Double_t)arr[aidx]);
873
874 pix.SetPedRMS(pix.GetPedRms()*sqrtnum, pix.GetPedRmsErr()*sqrtnum);
875 pix.SetSigma (pix.GetSigma()/pix.GetFFactorFADC2Phe());
876
877 FinalizeCharges(pix, fCam->GetAverageBadArea(aidx),"area id");
878 }
879
880 *fLog << endl;
881
882 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
883 {
884
885 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
886
887 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
888 FinalizePedestals(ped,pix, 0);
889 }
890
891 *fLog << endl;
892
893 //
894 // Finalize Bad Pixels
895 //
896 FinalizeBadPixels();
897
898 //
899 // Finalize F-Factor method
900 //
901 if (FinalizeFFactorMethod())
902 fCam->SetFFactorMethodValid(kTRUE);
903 else
904 {
905 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
906 fCam->SetFFactorMethodValid(kFALSE);
907 if (!fContinousCalibration)
908 return kFALSE;
909 }
910
911 *fLog << endl;
912
913 //
914 // Finalize Blind Pixel
915 //
916 fQECam->SetBlindPixelMethodValid(FinalizeBlindCam());
917
918 //
919 // Finalize PIN Diode
920 //
921 fQECam->SetBlindPixelMethodValid(FinalizePINDiode());
922
923 //
924 // Finalize QE Cam
925 //
926 FinalizeFFactorQECam();
927 FinalizeBlindPixelQECam();
928 FinalizePINDiodeQECam();
929 FinalizeCombinedQECam();
930
931 //
932 // Re-direct the output to an ascii-file from now on:
933 //
934 *fLog << inf << endl;
935 *fLog << GetDescriptor() << ": Fatal errors statistics:" << endl;
936
937 PrintUncalibrated(MBadPixelsPix::kChargeIsPedestal,
938 Form("%s%2.1f%s","Signal less than ",fChargeLimit," Pedestal RMS: "));
939 PrintUncalibrated(MBadPixelsPix::kChargeRelErrNotValid,
940 Form("%s%2.1f%s","Signal Error bigger than ",fChargeRelErrLimit," times Mean Signal: "));
941 PrintUncalibrated(MBadPixelsPix::kLoGainSaturation,
942 "Low Gain Saturation: ");
943// PrintUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin,
944// Form("%s%2.1f%s","Mean Abs. Arr. Time in First ",1.," Bin(s): "));
945// PrintUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins,
946// Form("%s%2.1f%s","Mean Abs. Arr. Time in Last ",2.," Bin(s): "));
947 PrintUncalibrated(MBadPixelsPix::kHiGainOverFlow,
948 "Pixels with High Gain Overflow: ");
949 PrintUncalibrated(MBadPixelsPix::kLoGainOverFlow,
950 "Pixels with Low Gain Overflow : ");
951// PrintUncalibrated(MBadPixelsPix::kFluctuatingArrivalTimes,
952// "Fluctuating Pulse Arrival Times: ");
953 PrintUncalibrated(MBadPixelsPix::kDeviatingAbsTimeResolution,
954 Form("%s%2.1f%s ","Abs.time rms more than ", fArrTimeRmsLimit, " dev from median: "));
955 PrintUncalibrated(MBadPixelsPix::kDeadPedestalRms,
956 "Presumably dead from Pedestal Rms: ");
957 PrintUncalibrated(MBadPixelsPix::kDeviatingNumPhes,
958 "Deviating number of phes: ");
959 PrintUncalibrated(MBadPixelsPix::kLoGainBlackout,
960 "Too many blackout events in low gain: ");
961 PrintUncalibrated(MBadPixelsPix::kPreviouslyExcluded,
962 "Previously excluded: ");
963
964 *fLog << inf << endl;
965 *fLog << GetDescriptor() << ": Unreliable errors statistics:" << endl;
966
967 PrintUncalibrated(MBadPixelsPix::kChargeSigmaNotValid,
968 "Signal Sigma smaller than Pedestal RMS: ");
969 PrintUncalibrated(MBadPixelsPix::kHiGainOscillating,
970 "Changing Hi Gain signal over time: ");
971 PrintUncalibrated(MBadPixelsPix::kLoGainOscillating,
972 "Changing Lo Gain signal over time: ");
973 PrintUncalibrated(MBadPixelsPix::kHiGainNotFitted,
974 "Unsuccesful Gauss fit to the Hi Gain: ");
975 PrintUncalibrated(MBadPixelsPix::kLoGainNotFitted,
976 "Unsuccesful Gauss fit to the Lo Gain: ");
977 PrintUncalibrated(MBadPixelsPix::kDeviatingFFactor,
978 "Deviating F-Factor: ");
979
980 fCam->SetReadyToSave();
981 fQECam->SetReadyToSave();
982 fBadPixels->SetReadyToSave();
983
984 if (fBlindCam)
985 fBlindCam->SetReadyToSave();
986 if (fPINDiode)
987 fPINDiode->SetReadyToSave();
988
989 //
990 // Finalize calibration statistics
991 //
992 return FinalizeUnsuitablePixels();
993}
994
995// ----------------------------------------------------------------------------------
996//
997// Retrieves pedestal and pedestal RMS from MPedestalPix
998// Retrieves total entries from MPedestalCam
999// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
1000// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
1001//
1002// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
1003// - MCalibrationChargePix::CalcLoGainPedestal()
1004//
1005void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
1006{
1007
1008 //
1009 // get the pedestals
1010 //
1011 const Float_t pedes = ped.GetPedestal();
1012 const Float_t prms = ped.GetPedestalRms();
1013 const Int_t num = fPedestals->GetNumSlices()*ped.GetNumEvents();
1014
1015 //
1016 // RMS error set by PedCalcFromLoGain, 0 in case MPedCalcPedRun was used.
1017 //
1018 const Float_t prmserr = num>0 ? prms/TMath::Sqrt(2.*num) : ped.GetPedestalRmsError();
1019
1020 //
1021 // set them in the calibration camera
1022 //
1023 if (cal.IsHiGainSaturation())
1024 {
1025 cal.SetPedestal(pedes * fNumLoGainSamples,
1026 prms * fSqrtLoGainSamples,
1027 prmserr * fSqrtLoGainSamples);
1028 cal.CalcLoGainPedestal(fNumLoGainSamples);
1029 }
1030 else
1031 {
1032
1033 cal.SetPedestal(pedes * fNumHiGainSamples,
1034 prms * fSqrtHiGainSamples,
1035 prmserr * fSqrtHiGainSamples);
1036 }
1037
1038}
1039
1040// ----------------------------------------------------------------------------------------------------
1041//
1042// Check fit results validity. Bad Pixels flags are set if:
1043//
1044// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
1045// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
1046// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
1047// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
1048// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
1049//
1050// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
1051//
1052// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
1053// and returns kFALSE if not succesful.
1054//
1055// Calls MCalibrationChargePix::CalcFFactor() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
1056// and returns kFALSE if not succesful.
1057//
1058// Calls MCalibrationChargePix::CalcConvFFactor()and sets flag: MBadPixelsPix::kDeviatingNumPhes)
1059// and returns kFALSE if not succesful.
1060//
1061Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad, const char* what)
1062{
1063
1064 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1065 return kFALSE;
1066
1067 const TString desc = Form("%7s%4d: ", what, cal.GetPixId());
1068
1069 if (cal.GetMean()<0)
1070 {
1071 *fLog << warn << desc << "Charge not fitted." << endl;
1072 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
1073 return kFALSE;
1074 }
1075
1076 if (cal.GetSigma()<0)
1077 {
1078 *fLog << warn << desc << "Charge Sigma invalid." << endl;
1079 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
1080 return kFALSE;
1081 }
1082
1083 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
1084 {
1085 *fLog << warn << desc
1086 << Form("Fitted Charge: %5.2f < %2.1f",cal.GetMean(),fChargeLimit)
1087 << Form(" * Pedestal RMS %5.2f",cal.GetPedRms()) << endl;
1088 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
1089 }
1090
1091 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
1092 {
1093 *fLog << warn << desc
1094 << Form("Fitted Charge: %4.2f < %2.1f",cal.GetMean(),fChargeRelErrLimit)
1095 << Form(" * its error %4.2f",cal.GetMeanErr()) << endl;
1096 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
1097 }
1098
1099 if (cal.GetSigma() < cal.GetPedRms())
1100 {
1101 *fLog << warn << desc << "Sigma of Fitted Charge: "
1102 << Form("%6.2f <",cal.GetSigma()) << " Ped. RMS="
1103 << Form("%5.2f", cal.GetPedRms()) << endl;
1104 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1105 return kFALSE;
1106 }
1107
1108 if (!cal.CalcReducedSigma())
1109 {
1110 *fLog << warn << desc << "Could not calculate the reduced sigma" << endl;
1111 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1112 return kFALSE;
1113 }
1114
1115 if (!cal.CalcFFactor())
1116 {
1117 *fLog << warn << desc << "Could not calculate the F-Factor"<< endl;
1118 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1119 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1120 return kFALSE;
1121 }
1122
1123 if (cal.GetPheFFactorMethod() < 0.)
1124 {
1125 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1126 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1127 cal.SetFFactorMethodValid(kFALSE);
1128 return kFALSE;
1129 }
1130
1131 if (!cal.CalcConvFFactor())
1132 {
1133 *fLog << warn << desc << "Could not calculate the Conv. FADC counts to Phes"<< endl;
1134 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1135 return kFALSE;
1136 }
1137
1138 if (!IsUseExtractorRes())
1139 return kTRUE;
1140
1141 if (!fExtractor)
1142 {
1143 *fLog << err << "Extractor resolution has been chosen, but no extractor is set. Cannot calibrate." << endl;
1144 return kFALSE;
1145 }
1146
1147 const Float_t resinphes = cal.IsHiGainSaturation()
1148 ? cal.GetPheFFactorMethod()*fExtractor->GetResolutionPerPheLoGain()
1149 : cal.GetPheFFactorMethod()*fExtractor->GetResolutionPerPheHiGain();
1150
1151 Float_t resinfadc = cal.IsHiGainSaturation()
1152 ? resinphes/cal.GetMeanConvFADC2Phe()/cal.GetConversionHiLo()
1153 : resinphes/cal.GetMeanConvFADC2Phe();
1154
1155 if (resinfadc > 3.0*cal.GetPedRms() )
1156 {
1157 *fLog << warn << desc << "Extractor Resolution " << Form("%5.2f", resinfadc) << " bigger than 3 Pedestal RMS "
1158 << Form("%4.2f", cal.GetPedRms()) << endl;
1159 resinfadc = 3.0*cal.GetPedRms();
1160 }
1161
1162 if (!cal.CalcReducedSigma(resinfadc))
1163 {
1164 *fLog << warn << desc << "Could not calculate the reduced sigma" << endl;
1165 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1166 return kFALSE;
1167 }
1168
1169 if (!cal.CalcFFactor())
1170 {
1171 *fLog << warn << desc << "Could not calculate the F-Factor" << endl;
1172 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1173 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1174 return kFALSE;
1175 }
1176
1177 if (cal.GetPheFFactorMethod() < 0.)
1178 {
1179 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1180 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1181 cal.SetFFactorMethodValid(kFALSE);
1182 return kFALSE;
1183 }
1184
1185 if (!cal.CalcConvFFactor())
1186 {
1187 *fLog << warn << desc << "Could not calculate the conv. FADC cts to phes" << endl;
1188 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1189 return kFALSE;
1190 }
1191
1192 return kTRUE;
1193}
1194
1195// -----------------------------------------------------------------------------------
1196//
1197// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
1198// - MBadPixelsPix::kChargeIsPedestal
1199// - MBadPixelsPix::kChargeRelErrNotValid
1200// - MBadPixelsPix::kMeanTimeInFirstBin
1201// - MBadPixelsPix::kMeanTimeInLast2Bins
1202// - MBadPixelsPix::kDeviatingNumPhes
1203// - MBadPixelsPix::kHiGainOverFlow
1204// - MBadPixelsPix::kLoGainOverFlow
1205//
1206// Sets pixel to MBadPixelsPix::kUnreliableRun, if one of the following flags is set:
1207// - MBadPixelsPix::kChargeSigmaNotValid
1208//
1209void MCalibrationChargeCalc::FinalizeBadPixels()
1210{
1211
1212 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1213 {
1214
1215 MBadPixelsPix &bad = (*fBadPixels)[i];
1216
1217 if (IsCheckDeadPixels())
1218 {
1219 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
1220 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1221 }
1222/*
1223 if (IsCheckExtractionWindow())
1224 {
1225 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
1226 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1227
1228 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
1229 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1230 }
1231 */
1232 if (IsCheckDeviatingBehavior())
1233 {
1234 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
1235 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1236 }
1237
1238 if (IsCheckHistOverflow())
1239 {
1240 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOverFlow ))
1241 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1242
1243 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOverFlow ))
1244 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1245 }
1246
1247 if (IsCheckArrivalTimes())
1248 {
1249 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingAbsTimeResolution ))
1250 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1251 //if (bad.IsUncalibrated( MBadPixelsPix::kFluctuatingArrivalTimes ))
1252 // bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1253 }
1254
1255 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
1256 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
1257 }
1258}
1259
1260// ------------------------------------------------------------------------
1261//
1262//
1263// First loop: Calculate a mean and mean RMS of photo-electrons per area index
1264// Include only pixels which are not MBadPixelsPix::kUnsuitableRun nor
1265// MBadPixelsPix::kChargeSigmaNotValid (see FinalizeBadPixels()) and set
1266// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
1267//
1268// Second loop: Get mean number of photo-electrons and its RMS including
1269// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
1270// and further exclude those deviating by more than fPheErrLimit mean
1271// sigmas from the mean (obtained in first loop). Set
1272// MBadPixelsPix::kDeviatingNumPhes if excluded.
1273//
1274// For the suitable pixels with flag MBadPixelsPix::kChargeSigmaNotValid
1275// set the number of photo-electrons as the mean number of photo-electrons
1276// calculated in that area index.
1277//
1278// Set weighted mean and variance of photo-electrons per area index in:
1279// average area pixels of MCalibrationChargeCam (obtained from:
1280// MCalibrationChargeCam::GetAverageArea() )
1281//
1282// Set weighted mean and variance of photo-electrons per sector in:
1283// average sector pixels of MCalibrationChargeCam (obtained from:
1284// MCalibrationChargeCam::GetAverageSector() )
1285//
1286//
1287// Third loop: Set mean number of photo-electrons and its RMS in the pixels
1288// only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1289//
1290Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
1291{
1292 const Int_t npixels = fGeom->GetNumPixels();
1293 const Int_t nareas = fGeom->GetNumAreas();
1294 const Int_t nsectors = fGeom->GetNumSectors();
1295
1296 TArrayF lowlim (nareas);
1297 TArrayF upplim (nareas);
1298 TArrayD areavars (nareas);
1299 TArrayD areaweights (nareas);
1300 TArrayD sectorweights (nsectors);
1301 TArrayD areaphes (nareas);
1302 TArrayD sectorphes (nsectors);
1303 TArrayI numareavalid (nareas);
1304 TArrayI numsectorvalid(nsectors);
1305
1306 //
1307 // First loop: Get mean number of photo-electrons and the RMS
1308 // The loop is only to recognize later pixels with very deviating numbers
1309 //
1310 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
1311
1312 for (Int_t i=0; i<npixels; i++)
1313 {
1314 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1315
1316 MBadPixelsPix &bad = (*fBadPixels)[i];
1317
1318 if (!pix.IsFFactorMethodValid())
1319 continue;
1320
1321 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1322 {
1323 pix.SetFFactorMethodValid(kFALSE);
1324 continue;
1325 }
1326
1327 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1328 continue;
1329
1330 if (!IsUseUnreliables())
1331 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1332 continue;
1333
1334 const Float_t nphe = pix.GetPheFFactorMethod();
1335 const Int_t aidx = (*fGeom)[i].GetAidx();
1336 camphes.Fill(i,nphe);
1337 camphes.SetUsed(i);
1338 areaphes [aidx] += nphe;
1339 areavars [aidx] += nphe*nphe;
1340 numareavalid[aidx] ++;
1341 }
1342
1343 for (Int_t i=0; i<nareas; i++)
1344 {
1345 if (numareavalid[i] == 0)
1346 {
1347 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
1348 << "in area index: " << i << endl;
1349 continue;
1350 }
1351
1352 if (numareavalid[i] == 1)
1353 areavars[i] = 0.;
1354 else
1355 {
1356 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1);
1357 areaphes[i] /= numareavalid[i];
1358 }
1359
1360 if (areavars[i] < 0.)
1361 {
1362 *fLog << warn << "Area " << setw(4) << i << ": No pixels with valid variance of photo-electrons found." << endl;
1363 continue;
1364 }
1365
1366 // FIXME: WHAT IS THIS FOR? It is overwritten!
1367 lowlim[i] = areaphes[i] - fPheErrLowerLimit*TMath::Sqrt(areavars[i]);
1368 upplim[i] = areaphes[i] + fPheErrUpperLimit*TMath::Sqrt(areavars[i]);
1369
1370
1371 TH1D *hist = camphes.ProjectionS(TArrayI(),TArrayI(1,&i),"_py",100);
1372 hist->Fit("gaus","Q0");
1373 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1374 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1375 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1376
1377 if (IsDebug())
1378 hist->DrawClone();
1379
1380 if (ndf < 5)
1381 {
1382 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons " << endl;
1383 *fLog << " in the camera with area index: " << i << endl;
1384 *fLog << " Number of dof.: " << ndf << " is smaller than 5 " << endl;
1385 *fLog << " Will use the simple mean and rms " << endl;
1386 delete hist;
1387 SetPheFitOK(i,kFALSE);
1388 continue;
1389 }
1390
1391 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1392
1393 if (prob < 0.001)
1394 {
1395 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons " << endl;
1396 *fLog << " in the camera with area index: " << i << endl;
1397 *fLog << " Fit probability " << prob << " is smaller than 0.001 " << endl;
1398 *fLog << " Will use the simple mean and rms " << endl;
1399 delete hist;
1400 SetPheFitOK(i,kFALSE);
1401 continue;
1402 }
1403
1404 if (mean < 0.)
1405 {
1406 *fLog << inf << "Area " << setw(4) << i << ": Fitted mean number of phe smaller 0." << endl;
1407 *fLog << warn << "Area " << setw(4) << i << ": Will use the simple mean and rms " << endl;
1408 SetPheFitOK(i,kFALSE);
1409 delete hist;
1410 continue;
1411 }
1412
1413 *fLog << inf << "Area " << setw(4) << i << ": Mean number of phes: "
1414 << Form("%7.2f+-%6.2f",mean,sigma) << endl;
1415
1416 lowlim[i] = mean - fPheErrLowerLimit*sigma;
1417 upplim[i] = mean + fPheErrUpperLimit*sigma;
1418
1419 if (lowlim[i]<1)
1420 lowlim[i] = 1;
1421
1422 delete hist;
1423
1424 SetPheFitOK(i,kTRUE);
1425 }
1426
1427 *fLog << endl;
1428
1429 numareavalid.Reset();
1430 areaphes .Reset();
1431 areavars .Reset();
1432 //
1433 // Second loop: Get mean number of photo-electrons and its RMS excluding
1434 // pixels deviating by more than fPheErrLimit sigma.
1435 // Set the conversion factor FADC counts to photo-electrons
1436 //
1437 for (Int_t i=0; i<npixels; i++)
1438 {
1439
1440 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1441
1442 if (!pix.IsFFactorMethodValid())
1443 continue;
1444
1445 MBadPixelsPix &bad = (*fBadPixels)[i];
1446
1447 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1448 continue;
1449
1450 if (!IsUseUnreliables())
1451 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1452 continue;
1453
1454 const Float_t nvar = pix.GetPheFFactorMethodVar();
1455 if (nvar <= 0.)
1456 {
1457 pix.SetFFactorMethodValid(kFALSE);
1458 continue;
1459 }
1460
1461 const Int_t aidx = (*fGeom)[i].GetAidx();
1462 const Int_t sector = (*fGeom)[i].GetSector();
1463 const Float_t area = (*fGeom)[i].GetA();
1464 const Float_t nphe = pix.GetPheFFactorMethod();
1465
1466 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1467 {
1468 *fLog << warn << "Pixel " << setw(4) << i << ": Num of phe "
1469 << Form("%7.2f out of +%3.1f-%3.1f sigma limit ",nphe,fPheErrUpperLimit,fPheErrLowerLimit)
1470 << Form("[%7.2f,%7.2f]",lowlim[aidx],upplim[aidx]) << endl;
1471 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1472
1473 if (IsCheckDeviatingBehavior())
1474 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1475 continue;
1476 }
1477
1478 areaweights [aidx] += nphe*nphe;
1479 areaphes [aidx] += nphe;
1480 numareavalid [aidx] ++;
1481
1482 if (aidx == 0)
1483 fNumInnerFFactorMethodUsed++;
1484
1485 sectorweights [sector] += nphe*nphe/area/area;
1486 sectorphes [sector] += nphe/area;
1487 numsectorvalid[sector] ++;
1488 }
1489
1490 *fLog << endl;
1491
1492 for (Int_t aidx=0; aidx<nareas; aidx++)
1493 {
1494
1495 MCalibrationChargePix &apix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
1496
1497 if (numareavalid[aidx] == 1)
1498 areaweights[aidx] = 0.;
1499 else if (numareavalid[aidx] == 0)
1500 {
1501 areaphes[aidx] = -1.;
1502 areaweights[aidx] = -1.;
1503 }
1504 else
1505 {
1506 areaweights[aidx] = (areaweights[aidx] - areaphes[aidx]*areaphes[aidx]/numareavalid[aidx])
1507 / (numareavalid[aidx]-1);
1508 areaphes[aidx] /= numareavalid[aidx];
1509 }
1510
1511 if (areaweights[aidx] < 0. || areaphes[aidx] <= 0.)
1512 {
1513 *fLog << warn << "Area " << setw(4) << aidx << ": Mean number phes could not be calculated: "
1514 << " Mean: " << areaphes[aidx]
1515 << " Var: " << areaweights[aidx] << endl;
1516 apix.SetFFactorMethodValid(kFALSE);
1517 continue;
1518 }
1519
1520 *fLog << inf << "Area " << setw(4) << aidx << ": Average total phes: "
1521 << Form("%7.2f +- %6.2f",areaphes[aidx],TMath::Sqrt(areaweights[aidx])) << endl;
1522
1523 apix.SetPheFFactorMethod ( areaphes[aidx] );
1524 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1525 apix.SetFFactorMethodValid ( kTRUE );
1526
1527 }
1528
1529 *fLog << endl;
1530
1531 for (Int_t sector=0; sector<nsectors; sector++)
1532 {
1533
1534 if (numsectorvalid[sector] == 1)
1535 sectorweights[sector] = 0.;
1536 else if (numsectorvalid[sector] == 0)
1537 {
1538 sectorphes[sector] = -1.;
1539 sectorweights[sector] = -1.;
1540 }
1541 else
1542 {
1543 sectorweights[sector] = (sectorweights[sector]
1544 - sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector]
1545 )
1546 / (numsectorvalid[sector]-1.);
1547 sectorphes[sector] /= numsectorvalid[sector];
1548 }
1549
1550 MCalibrationChargePix &spix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
1551
1552 if (sectorweights[sector] < 0. || sectorphes[sector] <= 0.)
1553 {
1554 *fLog << warn << "Sector " << setw(4) << sector
1555 <<": Mean number phes/area could not be calculated:"
1556 << " Mean: " << sectorphes[sector] << " Var: " << sectorweights[sector] << endl;
1557 spix.SetFFactorMethodValid(kFALSE);
1558 continue;
1559 }
1560
1561 *fLog << inf << "Sector " << setw(4) << sector
1562 << ": Avg number phes/mm^2: "
1563 << Form("%5.3f+-%4.3f",sectorphes[sector],TMath::Sqrt(sectorweights[sector]))
1564 << endl;
1565
1566 spix.SetPheFFactorMethod ( sectorphes[sector] );
1567 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1568 spix.SetFFactorMethodValid ( kTRUE );
1569
1570 }
1571
1572 //
1573 // Third loop: Set mean number of photo-electrons and its RMS in the pixels
1574 // only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1575 //
1576 for (Int_t i=0; i<npixels; i++)
1577 {
1578
1579 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1580
1581 MBadPixelsPix &bad = (*fBadPixels)[i];
1582
1583 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1584 continue;
1585
1586 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1587 {
1588 const Int_t aidx = (*fGeom)[i].GetAidx();
1589 MCalibrationChargePix &apix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
1590
1591 pix.SetPheFFactorMethod ( apix.GetPheFFactorMethod() );
1592 pix.SetPheFFactorMethodVar( apix.GetPheFFactorMethodVar() );
1593
1594 if (!pix.CalcConvFFactor())
1595 {
1596 *fLog << warn << GetDescriptor()
1597 << "Pixel " << setw(4) << pix.GetPixId()
1598 <<": Could not calculate the Conv. FADC counts to Phes"
1599 << endl;
1600 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1601 if (IsCheckDeviatingBehavior())
1602 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1603 }
1604
1605 }
1606 }
1607
1608 return kTRUE;
1609}
1610
1611
1612
1613// ------------------------------------------------------------------------
1614//
1615// Returns kFALSE if pointer to MCalibrationBlindCam is NULL
1616//
1617// The check returns kFALSE if:
1618//
1619// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1620// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1621//
1622// Calls:
1623// - MCalibrationBlindPix::CalcFluxInsidePlexiglass()
1624//
1625Bool_t MCalibrationChargeCalc::FinalizeBlindCam()
1626{
1627 if (!fBlindCam)
1628 return kFALSE;
1629
1630 Int_t nvalid = 0;
1631
1632 for (Int_t i=0; i<fBlindCam->GetSize(); i++)
1633 {
1634
1635 MCalibrationBlindPix &blindpix = (MCalibrationBlindPix&)(*fBlindCam)[i];
1636
1637 if (!blindpix.IsValid())
1638 continue;
1639
1640 const Float_t lambda = blindpix.GetLambda();
1641 const Float_t lambdaerr = blindpix.GetLambdaErr();
1642 const Float_t lambdacheck = blindpix.GetLambdaCheck();
1643
1644 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1645 {
1646 *fLog << warn << "BlindPix " << i << ": Lambda="
1647 << Form("%4.2f", lambda) << " and Lambda-Check="
1648 << Form("%4.2f", lambdacheck) << " differ by more than "
1649 << Form("%4.2f", fLambdaCheckLimit) << endl;
1650 blindpix.SetValid(kFALSE);
1651 continue;
1652 }
1653
1654 if (lambdaerr > fLambdaErrLimit)
1655 {
1656 *fLog << warn << "BlindPix " << i << ": Error of Fitted Lambda="
1657 << Form("%4.2f", lambdaerr) << " is greater than "
1658 << Form("%4.2f", fLambdaErrLimit) << endl;
1659 blindpix.SetValid(kFALSE);
1660 continue;
1661 }
1662
1663 if (!blindpix.CalcFluxInsidePlexiglass())
1664 {
1665 *fLog << warn << "Could not calculate the flux of photons from Blind Pixel Nr." << i << endl;
1666 blindpix.SetValid(kFALSE);
1667 continue;
1668 }
1669
1670 nvalid++;
1671 }
1672
1673 if (!nvalid)
1674 return kFALSE;
1675
1676 return kTRUE;
1677}
1678
1679// ------------------------------------------------------------------------
1680//
1681// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1682//
1683// The check returns kFALSE if:
1684//
1685// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1686// 2) PINDiode has a fit error smaller than fChargeErrLimit
1687// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1688// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1689//
1690// Calls:
1691// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1692//
1693Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1694{
1695
1696 if (!fPINDiode)
1697 return kFALSE;
1698
1699 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1700 {
1701 *fLog << warn << "PINDiode : Fitted Charge is smaller than "
1702 << fChargeLimit << " Pedestal RMS." << endl;
1703 return kFALSE;
1704 }
1705
1706 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1707 {
1708 *fLog << warn << "PINDiode : Error of Fitted Charge is smaller than "
1709 << fChargeErrLimit << endl;
1710 return kFALSE;
1711 }
1712
1713 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1714 {
1715 *fLog << warn << "PINDiode : Fitted Charge is smaller than "
1716 << fChargeRelErrLimit << "* its error" << endl;
1717 return kFALSE;
1718 }
1719
1720 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1721 {
1722 *fLog << warn << "PINDiode : Sigma of Fitted Charge smaller than Pedestal RMS" << endl;
1723 return kFALSE;
1724 }
1725
1726
1727 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1728 {
1729 *fLog << warn << "PINDiode : Could not calculate the flux of photons, "
1730 << "will skip PIN Diode Calibration " << endl;
1731 return kFALSE;
1732 }
1733
1734 return kTRUE;
1735}
1736
1737// ------------------------------------------------------------------------
1738//
1739// Calculate the average number of photons outside the plexiglass with the
1740// formula:
1741//
1742// av.Num.photons(area index) = av.Num.Phes(area index)
1743// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1744// / MCalibrationQEPix::GetPMTCollectionEff()
1745// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1746// / MCalibrationQECam::GetPlexiglassQE()
1747//
1748// Calculate the variance on the average number of photons assuming that the error on the
1749// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1750// values keeps it ordinary error since it is systematic.
1751//
1752// Loop over pixels:
1753//
1754// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1755// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1756//
1757// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1758// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1759//
1760// - Calculate the quantum efficiency with the formula:
1761//
1762// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1763//
1764// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1765//
1766// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1767// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1768//
1769// - Call MCalibrationQEPix::UpdateFFactorMethod()
1770//
1771void MCalibrationChargeCalc::FinalizeFFactorQECam()
1772{
1773
1774 if (fNumInnerFFactorMethodUsed < 2)
1775 {
1776 *fLog << warn << GetDescriptor()
1777 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1778 return;
1779 }
1780
1781 MCalibrationChargePix &avpix = (MCalibrationChargePix&)fCam->GetAverageArea(0);
1782 MCalibrationQEPix &qepix = (MCalibrationQEPix&) fQECam->GetAverageArea(0);
1783
1784 if (IsDebug())
1785 *fLog << dbginf << "External Phes: " << fExternalNumPhes
1786 << " Internal Phes: " << avpix.GetPheFFactorMethod() << endl;
1787
1788 const Float_t avphotons = ( IsUseExternalNumPhes()
1789 ? fExternalNumPhes
1790 : avpix.GetPheFFactorMethod() )
1791 / qepix.GetDefaultQE(fPulserColor)
1792 / qepix.GetPMTCollectionEff()
1793 / qepix.GetLightGuidesEff(fPulserColor)
1794 / fQECam->GetPlexiglassQE();
1795
1796 const Float_t avphotrelvar = ( IsUseExternalNumPhes()
1797 ? fExternalNumPhesRelVar
1798 : avpix.GetPheFFactorMethodRelVar() )
1799 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1800 + qepix.GetPMTCollectionEffRelVar()
1801 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1802 + fQECam->GetPlexiglassQERelVar();
1803
1804 const UInt_t nareas = fGeom->GetNumAreas();
1805
1806 //
1807 // Set the results in the MCalibrationChargeCam
1808 //
1809 fCam->SetNumPhotonsFFactorMethod (avphotons);
1810
1811 if (avphotrelvar > 0.)
1812 fCam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1813
1814 TArrayF lowlim (nareas);
1815 TArrayF upplim (nareas);
1816 TArrayD avffactorphotons (nareas);
1817 TArrayD avffactorphotvar (nareas);
1818 TArrayI numffactor (nareas);
1819
1820 const UInt_t npixels = fGeom->GetNumPixels();
1821
1822 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1823
1824 for (UInt_t i=0; i<npixels; i++)
1825 {
1826
1827 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1828 MCalibrationQEPix &qpix = (MCalibrationQEPix&) (*fQECam)[i];
1829
1830 MBadPixelsPix &bad = (*fBadPixels)[i];
1831
1832 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1833 continue;
1834
1835 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1836 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1837
1838 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1839
1840 qpix.SetQEFFactor ( qe , fPulserColor );
1841 qpix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1842 qpix.SetFFactorMethodValid( kTRUE , fPulserColor );
1843
1844 if (!qpix.UpdateFFactorMethod( fQECam->GetPlexiglassQE() ))
1845 *fLog << warn << GetDescriptor()
1846 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1847
1848 //
1849 // The following pixels are those with deviating sigma, but otherwise OK,
1850 // probably those with stars during the pedestal run, but not the cal. run.
1851 //
1852 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1853 {
1854 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1855 if (IsCheckDeviatingBehavior())
1856 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1857 continue;
1858 }
1859
1860 const Int_t aidx = (*fGeom)[i].GetAidx();
1861 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1862
1863 camffactor.Fill(i,ffactor);
1864 camffactor.SetUsed(i);
1865
1866 avffactorphotons[aidx] += ffactor;
1867 avffactorphotvar[aidx] += ffactor*ffactor;
1868 numffactor[aidx]++;
1869 }
1870
1871 for (UInt_t i=0; i<nareas; i++)
1872 {
1873
1874 if (numffactor[i] == 0)
1875 {
1876 *fLog << warn << GetDescriptor() << ": No pixels with valid total F-Factor found "
1877 << "in area index: " << i << endl;
1878 continue;
1879 }
1880
1881 avffactorphotvar[i] = (avffactorphotvar[i] - avffactorphotons[i]*avffactorphotons[i]/numffactor[i])
1882 / (numffactor[i]-1.);
1883 avffactorphotons[i] = avffactorphotons[i] / numffactor[i];
1884
1885 if (avffactorphotvar[i] < 0.)
1886 {
1887 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of total F-Factor found "
1888 << "in area index: " << i << endl;
1889 continue;
1890 }
1891
1892 lowlim [i] = 1.; // Lowest known F-Factor of a PMT
1893 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1894
1895 TArrayI area(1);
1896 area[0] = i;
1897
1898 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1899 hist->Fit("gaus","Q0");
1900 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1901 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1902 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1903
1904 if (IsDebug())
1905 camffactor.DrawClone();
1906
1907 if (ndf < 2)
1908 {
1909 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor " << endl;
1910 *fLog << " in the camera with area index: " << i << endl;
1911 *fLog << " Number of dof.: " << ndf << " is smaller than 2 " << endl;
1912 *fLog << " Will use the simple mean and rms." << endl;
1913 delete hist;
1914 SetFFactorFitOK(i,kFALSE);
1915 continue;
1916 }
1917
1918 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1919
1920 if (prob < 0.001)
1921 {
1922 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor " << endl;
1923 *fLog << " in the camera with area index: " << i << endl;
1924 *fLog << " Fit probability " << prob << " is smaller than 0.001 " << endl;
1925 *fLog << " Will use the simple mean and rms." << endl;
1926 delete hist;
1927 SetFFactorFitOK(i,kFALSE);
1928 continue;
1929 }
1930
1931 *fLog << inf << "Area " << setw(4) << i <<": Mean F-Factor :"
1932 << Form("%4.2f+-%4.2f",mean,sigma) << endl;
1933
1934 lowlim [i] = 1.;
1935 upplim [i] = mean + fFFactorErrLimit*sigma;
1936
1937 delete hist;
1938
1939 SetFFactorFitOK(i,kTRUE);
1940 }
1941
1942 *fLog << endl;
1943
1944 for (UInt_t i=0; i<npixels; i++)
1945 {
1946
1947 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1948 MBadPixelsPix &bad = (*fBadPixels)[i];
1949
1950 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1951 continue;
1952
1953 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1954 const Int_t aidx = (*fGeom)[i].GetAidx();
1955
1956 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1957 {
1958 *fLog << warn << "Pixel " << setw(4) << i<< ": Overall F-Factor "
1959 << Form("%5.2f",ffactor) << " out of range ["
1960 << Form("%5.2f,%5.2f",lowlim[aidx],upplim[aidx]) << "]" << endl;
1961
1962 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1963 if (IsCheckDeviatingBehavior())
1964 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1965 }
1966 }
1967
1968 for (UInt_t i=0; i<npixels; i++)
1969 {
1970
1971 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1972 MCalibrationQEPix &qpix = (MCalibrationQEPix&)(*fQECam)[i];
1973
1974 if ((*fBadPixels)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1975 {
1976 qpix.SetFFactorMethodValid(kFALSE,fPulserColor);
1977 pix.SetFFactorMethodValid(kFALSE);
1978 continue;
1979 }
1980 }
1981}
1982
1983
1984// ------------------------------------------------------------------------
1985//
1986// Loop over pixels:
1987//
1988// - Continue, if not MCalibrationBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1989// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1990//
1991// - Calculate the quantum efficiency with the formula:
1992//
1993// QE = Num.Phes / MCalibrationBlindPix::GetFluxInsidePlexiglass()
1994// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1995//
1996// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1997// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1998// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1999//
2000// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
2001//
2002void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
2003{
2004
2005 if (!fBlindCam)
2006 return;
2007
2008 //
2009 // Set the results in the MCalibrationChargeCam
2010 //
2011 if (!fBlindCam->IsFluxInsidePlexiglassAvailable())
2012 {
2013
2014 const Float_t photons = fBlindCam->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
2015 / fQECam->GetPlexiglassQE();
2016 fCam->SetNumPhotonsBlindPixelMethod(photons);
2017
2018 const Float_t photrelvar = fBlindCam->GetFluxInsidePlexiglassRelVar()
2019 + fQECam->GetPlexiglassQERelVar();
2020
2021 if (photrelvar > 0.)
2022 fCam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
2023 }
2024
2025 //
2026 // With the knowledge of the overall photon flux, calculate the
2027 // quantum efficiencies after the Blind Pixel and PIN Diode method
2028 //
2029 const UInt_t npixels = fGeom->GetNumPixels();
2030 for (UInt_t i=0; i<npixels; i++)
2031 {
2032
2033 MCalibrationQEPix &qepix = (MCalibrationQEPix&)(*fQECam)[i];
2034
2035 if (!fBlindCam->IsFluxInsidePlexiglassAvailable())
2036 {
2037 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
2038 continue;
2039 }
2040
2041 MBadPixelsPix &bad = (*fBadPixels)[i];
2042
2043 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2044 {
2045 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
2046 continue;
2047 }
2048
2049 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
2050 MGeom &geo = (*fGeom)[i];
2051
2052 const Float_t qe = pix.GetPheFFactorMethod()
2053 / fBlindCam->GetFluxInsidePlexiglass()
2054 / geo.GetA()
2055 * fQECam->GetPlexiglassQE();
2056
2057 const Float_t qerelvar = fBlindCam->GetFluxInsidePlexiglassRelVar()
2058 + fQECam->GetPlexiglassQERelVar()
2059 + pix.GetPheFFactorMethodRelVar();
2060
2061 qepix.SetQEBlindPixel ( qe , fPulserColor );
2062 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
2063 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
2064
2065 if (!qepix.UpdateBlindPixelMethod( fQECam->GetPlexiglassQE()))
2066 *fLog << warn << GetDescriptor()
2067 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
2068 }
2069}
2070
2071// ------------------------------------------------------------------------
2072//
2073// Loop over pixels:
2074//
2075// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
2076// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
2077//
2078// - Calculate the quantum efficiency with the formula:
2079//
2080// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
2081//
2082// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
2083// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
2084// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
2085//
2086// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
2087//
2088void MCalibrationChargeCalc::FinalizePINDiodeQECam()
2089{
2090
2091 const UInt_t npixels = fGeom->GetNumPixels();
2092
2093 if (!fPINDiode)
2094 return;
2095
2096 //
2097 // With the knowledge of the overall photon flux, calculate the
2098 // quantum efficiencies after the PIN Diode method
2099 //
2100 for (UInt_t i=0; i<npixels; i++)
2101 {
2102
2103 MCalibrationQEPix &qepix = (MCalibrationQEPix&)(*fQECam)[i];
2104
2105 if (!fPINDiode)
2106 {
2107 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2108 continue;
2109 }
2110
2111 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
2112 {
2113 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2114 continue;
2115 }
2116
2117 MBadPixelsPix &bad = (*fBadPixels)[i];
2118
2119 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2120 {
2121 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2122 continue;
2123 }
2124
2125 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
2126 MGeom &geo = (*fGeom)[i];
2127
2128 const Float_t qe = pix.GetPheFFactorMethod()
2129 / fPINDiode->GetFluxOutsidePlexiglass()
2130 / geo.GetA();
2131
2132 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
2133
2134 qepix.SetQEPINDiode ( qe , fPulserColor );
2135 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
2136 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
2137
2138 if (!qepix.UpdatePINDiodeMethod())
2139 *fLog << warn << GetDescriptor()
2140 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
2141 }
2142}
2143
2144// ------------------------------------------------------------------------
2145//
2146// Loop over pixels:
2147//
2148// - Call MCalibrationQEPix::UpdateCombinedMethod()
2149//
2150void MCalibrationChargeCalc::FinalizeCombinedQECam()
2151{
2152
2153 const UInt_t npixels = fGeom->GetNumPixels();
2154
2155 for (UInt_t i=0; i<npixels; i++)
2156 {
2157
2158 MCalibrationQEPix &qepix = (MCalibrationQEPix&)(*fQECam)[i];
2159
2160 if (!(*fBadPixels)[i].IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2161 {
2162 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2163 continue;
2164 }
2165
2166 qepix.UpdateCombinedMethod();
2167 }
2168}
2169
2170// -----------------------------------------------------------------------------------------------
2171//
2172// - Print out statistics about BadPixels of type UnsuitableType_t
2173// - store numbers of bad pixels of each type in fCam
2174//
2175Bool_t MCalibrationChargeCalc::FinalizeUnsuitablePixels()
2176{
2177
2178 *fLog << inf << endl;
2179 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
2180 *fLog << dec << setfill(' ');
2181
2182 const Int_t nareas = fGeom->GetNumAreas();
2183
2184 TArrayI suit(nareas);
2185 TArrayI unsuit(nareas);
2186 TArrayI unrel(nareas);
2187
2188 Int_t unsuitcnt=0;
2189 Int_t unrelcnt =0;
2190
2191 // Count number of succesfully calibrated pixels
2192 for (Int_t aidx=0; aidx<nareas; aidx++)
2193 {
2194 suit[aidx] = fBadPixels->GetNumSuitable(MBadPixelsPix::kUnsuitableRun, fGeom, aidx);
2195 unsuit[aidx] = fBadPixels->GetNumUnsuitable(MBadPixelsPix::kUnsuitableRun, fGeom, aidx);
2196 unrel[aidx] = fBadPixels->GetNumUnsuitable(MBadPixelsPix::kUnreliableRun, fGeom, aidx);
2197
2198 unsuitcnt += unsuit[aidx];
2199 unrelcnt += unrel[aidx];
2200
2201 fCam->SetNumUnsuitable(unsuit[aidx], aidx);
2202 fCam->SetNumUnreliable(unrel[aidx], aidx);
2203 }
2204
2205 TArrayI counts(nareas);
2206 for (Int_t i=0; i<fCam->GetSize(); i++)
2207 {
2208 const MCalibrationPix &pix = (*fCam)[i];
2209 if (pix.IsHiGainSaturation())
2210 {
2211 const Int_t aidx = (*fGeom)[i].GetAidx();
2212 counts[aidx]++;
2213 }
2214 }
2215
2216 if (fGeom->InheritsFrom("MGeomCamMagic"))
2217 {
2218 *fLog << " Successfully calibrated Pixels: Inner: "
2219 << Form("%3i",suit[0]) << " Outer: " << Form("%3i",suit[1]) << endl;
2220 *fLog << " Uncalibrated Pixels: Inner: "
2221 << Form("%3i",unsuit[0]) << " Outer: " << Form("%3i",unsuit[1]) << endl;
2222 *fLog << " Unreliable Pixels: Inner: "
2223 << Form("%3i",unrel[0]) << " Outer: " << Form("%3i",unrel[1]) << endl;
2224 *fLog << " High-gain saturated Pixels: Inner: "
2225 << Form("%3i",counts[0]) << " Outer: " << Form("%3i",counts[1]) << endl;
2226 *fLog << endl;
2227 }
2228
2229 if (unsuitcnt > fUnsuitablesLimit*fGeom->GetNumPixels())
2230 {
2231 *fLog << err << "Number of unsuitable pixels: " << 100.*unsuitcnt/fGeom->GetNumPixels()
2232 << "% exceeds limit of " << fUnsuitablesLimit*100 << "%" << endl;
2233 return kFALSE;
2234 }
2235
2236 if (unrelcnt > fUnreliablesLimit*fGeom->GetNumPixels())
2237 {
2238 *fLog << err << "Relative number of unreliable pixels: " << 100.*unrelcnt/fGeom->GetNumPixels()
2239 << "% exceeds limit of " << fUnreliablesLimit*100 << "%" << endl;
2240 return kFALSE;
2241 }
2242 return kTRUE;
2243}
2244
2245// -----------------------------------------------------------------------------------------------
2246//
2247// Print out statistics about BadPixels of type UncalibratedType_t
2248//
2249void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
2250{
2251
2252 UInt_t countinner = 0;
2253 UInt_t countouter = 0;
2254
2255 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
2256 {
2257 if ((*fBadPixels)[i].IsUncalibrated(typ))
2258 {
2259 if (fGeom->GetPixRatio(i) == 1.)
2260 countinner++;
2261 else
2262 countouter++;
2263 }
2264 }
2265
2266 *fLog << " " << setw(7) << text << "Inner: " << Form("%3i",countinner)
2267 << " Outer: " << Form("%3i", countouter) << endl;
2268}
2269
2270// --------------------------------------------------------------------------
2271//
2272// Read the environment for the following data members:
2273// - fChargeLimit
2274// - fChargeErrLimit
2275// - fChargeRelErrLimit
2276// - fDebug
2277// - fFFactorErrLimit
2278// - fLambdaErrLimit
2279// - fLambdaCheckErrLimit
2280// - fPheErrLimit
2281//
2282Int_t MCalibrationChargeCalc::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
2283{
2284
2285 Bool_t rc = kFALSE;
2286 if (IsEnvDefined(env, prefix, "ChargeLimit", print))
2287 {
2288 SetChargeLimit(GetEnvValue(env, prefix, "ChargeLimit", fChargeLimit));
2289 rc = kTRUE;
2290 }
2291 if (IsEnvDefined(env, prefix, "ChargeErrLimit", print))
2292 {
2293 SetChargeErrLimit(GetEnvValue(env, prefix, "ChargeErrLimit", fChargeErrLimit));
2294 rc = kTRUE;
2295 }
2296 if (IsEnvDefined(env, prefix, "ChargeRelErrLimit", print))
2297 {
2298 SetChargeRelErrLimit(GetEnvValue(env, prefix, "ChargeRelErrLimit", fChargeRelErrLimit));
2299 rc = kTRUE;
2300 }
2301 if (IsEnvDefined(env, prefix, "Debug", print))
2302 {
2303 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
2304 rc = kTRUE;
2305 }
2306 if (IsEnvDefined(env, prefix, "ArrTimeRmsLimit", print))
2307 {
2308 SetArrTimeRmsLimit(GetEnvValue(env, prefix, "ArrTimeRmsLimit", fArrTimeRmsLimit));
2309 rc = kTRUE;
2310 }
2311 if (IsEnvDefined(env, prefix, "FFactorErrLimit", print))
2312 {
2313 SetFFactorErrLimit(GetEnvValue(env, prefix, "FFactorErrLimit", fFFactorErrLimit));
2314 rc = kTRUE;
2315 }
2316 if (IsEnvDefined(env, prefix, "LambdaErrLimit", print))
2317 {
2318 SetLambdaErrLimit(GetEnvValue(env, prefix, "LambdaErrLimit", fLambdaErrLimit));
2319 rc = kTRUE;
2320 }
2321 if (IsEnvDefined(env, prefix, "LambdaCheckLimit", print))
2322 {
2323 SetLambdaCheckLimit(GetEnvValue(env, prefix, "LambdaCheckLimit", fLambdaCheckLimit));
2324 rc = kTRUE;
2325 }
2326 if (IsEnvDefined(env, prefix, "CheckDeadPixels", print))
2327 {
2328 SetCheckDeadPixels(GetEnvValue(env, prefix, "CheckDeadPixels", IsCheckDeadPixels()));
2329 rc = kTRUE;
2330 }
2331 if (IsEnvDefined(env, prefix, "CheckDeviatingBehavior", print))
2332 {
2333 SetCheckDeviatingBehavior(GetEnvValue(env, prefix, "CheckDeviatingBehavior", IsCheckDeviatingBehavior()));
2334 rc = kTRUE;
2335 }
2336 if (IsEnvDefined(env, prefix, "CheckExtractionWindow", print))
2337 {
2338 SetCheckExtractionWindow(GetEnvValue(env, prefix, "CheckExtractionWindow", IsCheckExtractionWindow()));
2339 rc = kTRUE;
2340 }
2341 if (IsEnvDefined(env, prefix, "CheckHistOverflow", print))
2342 {
2343 SetCheckHistOverflow(GetEnvValue(env, prefix, "CheckHistOverflow", IsCheckHistOverflow()));
2344 rc = kTRUE;
2345 }
2346 if (IsEnvDefined(env, prefix, "CheckOscillations", print))
2347 {
2348 SetCheckOscillations(GetEnvValue(env, prefix, "CheckOscillations", IsCheckOscillations()));
2349 rc = kTRUE;
2350 }
2351 if (IsEnvDefined(env, prefix, "CheckArrivalTimes", print))
2352 {
2353 SetCheckArrivalTimes(GetEnvValue(env, prefix, "CheckArrivalTimes", IsCheckArrivalTimes()));
2354 rc = kTRUE;
2355 }
2356 if (IsEnvDefined(env, prefix, "PheErrLowerLimit", print))
2357 {
2358 SetPheErrLowerLimit(GetEnvValue(env, prefix, "PheErrLowerLimit", fPheErrLowerLimit));
2359 rc = kTRUE;
2360 }
2361 if (IsEnvDefined(env, prefix, "PheErrUpperLimit", print))
2362 {
2363 SetPheErrUpperLimit(GetEnvValue(env, prefix, "PheErrUpperLimit", fPheErrUpperLimit));
2364 rc = kTRUE;
2365 }
2366 if (IsEnvDefined(env, prefix, "UseExtractorRes", print))
2367 {
2368 SetUseExtractorRes(GetEnvValue(env, prefix, "UseExtractorRes", IsUseExtractorRes()));
2369 rc = kTRUE;
2370 }
2371 if (IsEnvDefined(env, prefix, "UseUnreliables", print))
2372 {
2373 SetUseUnreliables(GetEnvValue(env, prefix, "UseUnreliables", IsUseUnreliables()));
2374 rc = kTRUE;
2375 }
2376
2377 if (IsEnvDefined(env, prefix, "UseExternalNumPhes", print))
2378 {
2379 SetUseExternalNumPhes(GetEnvValue(env, prefix, "UseExternalNumPhes", IsUseExternalNumPhes()));
2380 rc = kTRUE;
2381 }
2382
2383 if (IsEnvDefined(env, prefix, "UnsuitablesLimit", print))
2384 {
2385 SetUnsuitablesLimit(GetEnvValue(env, prefix, "UnsuitablesLimit", fUnsuitablesLimit));
2386 rc = kTRUE;
2387 }
2388
2389 if (IsEnvDefined(env, prefix, "UnreliablesLimit", print))
2390 {
2391 SetUnreliablesLimit(GetEnvValue(env, prefix, "UnreliablesLimit", fUnreliablesLimit));
2392 rc = kTRUE;
2393 }
2394
2395
2396 return rc;
2397}
Note: See TracBrowser for help on using the repository browser.