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

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