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

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