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

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