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

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