source: trunk/MagicSoft/Mars/mcalib/MCalibrationChargeCalc.cc@ 9423

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