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

Last change on this file since 7669 was 7292, 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 ? (MCalibrationBlindCam*) fIntensBlind->GetCam() : fBlindCam;
742 MCalibrationQECam *qecam = fIntensQE ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
743 MCalibrationChargeCam *chargecam = fIntensCam ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
744 MBadPixelsCam *badcam = fIntensBad ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
745
746 //
747 // First loop over pixels, call FinalizePedestals and FinalizeCharges
748 //
749 Int_t nvalid = 0;
750
751 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
752 {
753
754 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[pixid];
755 //
756 // Check if the pixel has been excluded from the fits
757 //
758 if (pix.IsExcluded())
759 continue;
760
761 MPedestalPix &ped = (*fPedestals)[pixid];
762 MBadPixelsPix &bad = (*badcam) [pixid];
763
764 const Int_t aidx = (*fGeom)[pixid].GetAidx();
765
766 FinalizePedestals(ped,pix,aidx);
767
768 if (FinalizeCharges(pix,bad,"pixel "))
769 nvalid++;
770
771 FinalizeArrivalTimes(pix,bad,"pixel ");
772 }
773
774 *fLog << endl;
775
776 //
777 // The Michele check ...
778 //
779 if (nvalid == 0)
780 {
781 if (!fIntensCam)
782 {
783 *fLog << warn << GetDescriptor() << ": All pixels have non-valid calibration. "
784 << "Did you forget to fill the histograms "
785 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
786 *fLog << warn << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
787 << "instead of a calibration run " << endl;
788 return kFALSE;
789 }
790 }
791
792 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
793 {
794
795 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
796 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
797 const MArrayI &arr = fHCam->GetAverageAreaNum();
798
799 FinalizePedestals(ped,pix,aidx);
800
801 //
802 // Correct for sqrt(number of valid pixels) in the pedestal RMS
803 // (already done for calibration sigma in MHCalibrationCam::CalcAverageSigma()
804 //
805 pix.SetPedRMS(pix.GetPedRms()*TMath::Sqrt((Float_t)arr[aidx]),
806 pix.GetPedRmsErr()*TMath::Sqrt((Float_t)arr[aidx]));
807 pix.SetSigma (pix.GetSigma()/pix.GetFFactorFADC2Phe());
808
809 FinalizeCharges(pix, chargecam->GetAverageBadArea(aidx),"area id");
810 FinalizeArrivalTimes(pix, chargecam->GetAverageBadArea(aidx), "area id");
811 }
812
813 *fLog << endl;
814
815 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
816 {
817
818 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
819
820 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
821 FinalizePedestals(ped,pix, 0);
822 }
823
824 *fLog << endl;
825
826 //
827 // Finalize Bad Pixels
828 //
829 FinalizeBadPixels();
830
831 //
832 // Finalize F-Factor method
833 //
834 if (FinalizeFFactorMethod())
835 chargecam->SetFFactorMethodValid(kTRUE);
836 else
837 {
838 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
839 chargecam->SetFFactorMethodValid(kFALSE);
840 if (!fIntensCam)
841 return kFALSE;
842 }
843
844 *fLog << endl;
845
846 //
847 // Finalize Blind Pixel
848 //
849 qecam->SetBlindPixelMethodValid(FinalizeBlindCam());
850
851 //
852 // Finalize PIN Diode
853 //
854 qecam->SetBlindPixelMethodValid(FinalizePINDiode());
855
856 //
857 // Finalize QE Cam
858 //
859 FinalizeFFactorQECam();
860 FinalizeBlindPixelQECam();
861 FinalizePINDiodeQECam();
862 FinalizeCombinedQECam();
863
864 //
865 // Re-direct the output to an ascii-file from now on:
866 //
867 MLog *oldlog = fLog;
868 MLog asciilog;
869 if (!fOutputFile.IsNull())
870 {
871 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
872 SetLogStream(&asciilog);
873 }
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 if (!fOutputFile.IsNull())
920 SetLogStream(oldlog);
921
922 chargecam->SetReadyToSave();
923 qecam ->SetReadyToSave();
924 badcam ->SetReadyToSave();
925
926 if (blindcam)
927 blindcam->SetReadyToSave();
928 if (fPINDiode)
929 fPINDiode->SetReadyToSave();
930
931 //
932 // Finalize calibration statistics
933 //
934 return FinalizeUnsuitablePixels();
935}
936
937// ----------------------------------------------------------------------------------
938//
939// Retrieves pedestal and pedestal RMS from MPedestalPix
940// Retrieves total entries from MPedestalCam
941// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
942// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
943//
944// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
945// - MCalibrationChargePix::CalcLoGainPedestal()
946//
947void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
948{
949
950 //
951 // get the pedestals
952 //
953 const Float_t pedes = ped.GetPedestal();
954 const Float_t prms = ped.GetPedestalRms();
955 const Int_t num = fPedestals->GetTotalEntries();
956
957 //
958 // RMS error set by PedCalcFromLoGain, 0 in case MPedCalcPedRun was used.
959 //
960 const Float_t prmserr = num>0 ? prms/TMath::Sqrt(2.*num) : ped.GetPedestalRmsError();
961
962 //
963 // set them in the calibration camera
964 //
965 if (cal.IsHiGainSaturation())
966 {
967 cal.SetPedestal(pedes * fNumLoGainSamples,
968 prms * fSqrtLoGainSamples,
969 prmserr * fSqrtLoGainSamples);
970 cal.CalcLoGainPedestal(fNumLoGainSamples);
971 }
972 else
973 {
974
975 cal.SetPedestal(pedes * fNumHiGainSamples,
976 prms * fSqrtHiGainSamples,
977 prmserr * fSqrtHiGainSamples);
978 }
979
980}
981
982// ----------------------------------------------------------------------------------------------------
983//
984// Check fit results validity. Bad Pixels flags are set if:
985//
986// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
987// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
988// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
989// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
990// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
991//
992// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
993//
994// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
995// and returns kFALSE if not succesful.
996//
997// Calls MCalibrationChargePix::CalcFFactor() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
998// and returns kFALSE if not succesful.
999//
1000// Calls MCalibrationChargePix::CalcConvFFactor()and sets flag: MBadPixelsPix::kDeviatingNumPhes)
1001// and returns kFALSE if not succesful.
1002//
1003Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad, const char* what)
1004{
1005
1006 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1007 return kFALSE;
1008
1009 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
1010 {
1011 *fLog << warn
1012 << Form("Fitted Charge: %5.2f < %2.1f",cal.GetMean(),fChargeLimit)
1013 << Form(" * Pedestal RMS %5.2f in %s%3i",cal.GetPedRms(),what,cal.GetPixId()) << endl;
1014 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
1015 }
1016
1017 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
1018 {
1019 *fLog << warn
1020 << Form("Fitted Charge: %4.2f < %2.1f",cal.GetMean(),fChargeRelErrLimit)
1021 << Form(" * its error %4.2f in %s%3i",cal.GetMeanErr(),what,cal.GetPixId()) << endl;
1022 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
1023 }
1024
1025 if (cal.GetSigma() < cal.GetPedRms())
1026 {
1027 *fLog << warn << "Sigma of Fitted Charge: "
1028 << Form("%6.2f <",cal.GetSigma()) << " Ped. RMS="
1029 << Form("%5.2f", cal.GetPedRms()) << " in " << what
1030 << Form("%3i",cal.GetPixId()) << endl;
1031 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1032 return kFALSE;
1033 }
1034
1035 if (!cal.CalcReducedSigma())
1036 {
1037 *fLog << warn << "Could not calculate the reduced sigma in " << what
1038 << ": " << Form("%4i",cal.GetPixId())
1039 << endl;
1040 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1041 return kFALSE;
1042 }
1043
1044 if (!cal.CalcFFactor())
1045 {
1046 *fLog << warn << "Could not calculate the F-Factor in " << what
1047 << ": " << Form("%4i",cal.GetPixId())
1048 << endl;
1049 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1050 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1051 return kFALSE;
1052 }
1053
1054 if (cal.GetPheFFactorMethod() < 0.)
1055 {
1056 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1057 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1058 cal.SetFFactorMethodValid(kFALSE);
1059 return kFALSE;
1060 }
1061
1062 if (!cal.CalcConvFFactor())
1063 {
1064 *fLog << warn << "Could not calculate the Conv. FADC counts to Phes in ";
1065 *fLog << what << ": " << Form("%4i", cal.GetPixId()) << endl;
1066 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1067 return kFALSE;
1068 }
1069
1070 if (!IsUseExtractorRes())
1071 return kTRUE;
1072
1073 if (!fExtractor)
1074 {
1075 *fLog << err << "Extractor resolution has been chosen, but not extractor is set. Cannot calibrate" << endl;
1076 return kFALSE;
1077 }
1078
1079 const Float_t resinphes = cal.IsHiGainSaturation()
1080 ? cal.GetPheFFactorMethod()*fExtractor->GetResolutionPerPheLoGain()
1081 : cal.GetPheFFactorMethod()*fExtractor->GetResolutionPerPheHiGain();
1082
1083 Float_t resinfadc = cal.IsHiGainSaturation()
1084 ? resinphes/cal.GetMeanConvFADC2Phe()/cal.GetConversionHiLo()
1085 : resinphes/cal.GetMeanConvFADC2Phe();
1086
1087 if (resinfadc > 3.0*cal.GetPedRms() )
1088 {
1089 *fLog << warn << " Extractor Resolution: " << resinfadc << " bigger than 3 Pedestal RMS "
1090 << cal.GetPedRms() << endl;
1091 resinfadc = 3.0*cal.GetPedRms();
1092 }
1093
1094 if (!cal.CalcReducedSigma(resinfadc))
1095 {
1096 *fLog << warn << "Could not calculate the reduced sigma in " << what;
1097 *fLog << ": " << Form("%4i",cal.GetPixId()) << endl;
1098 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1099 return kFALSE;
1100 }
1101
1102 if (!cal.CalcFFactor())
1103 {
1104 *fLog << warn << "Could not calculate the F-Factor in " << what;
1105 *fLog << ": " << Form("%4i",cal.GetPixId()) << endl;
1106 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1107 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1108 return kFALSE;
1109 }
1110
1111 if (cal.GetPheFFactorMethod() < 0.)
1112 {
1113 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1114 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1115 cal.SetFFactorMethodValid(kFALSE);
1116 return kFALSE;
1117 }
1118
1119 if (!cal.CalcConvFFactor())
1120 {
1121 *fLog << warn << "Could not calculate the Conv. FADC counts to Phes in "
1122 << what << ": " << Form("%4i",cal.GetPixId())
1123 << endl;
1124 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1125 return kFALSE;
1126 }
1127
1128 return kTRUE;
1129}
1130
1131// -----------------------------------------------------------------------------------
1132//
1133// Test the arrival Times RMS of the pixel and set the bit
1134// - MBadPixelsPix::kFluctuatingArrivalTimes
1135//
1136void MCalibrationChargeCalc::FinalizeArrivalTimes(MCalibrationChargePix &cal, MBadPixelsPix &bad, const char* what)
1137{
1138 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1139 return;
1140
1141 if (cal.GetAbsTimeRms() > fArrTimeRmsLimit)
1142 {
1143 *fLog << warn;
1144 *fLog << "RMS of pulse arrival times: " << Form("%2.1f", cal.GetAbsTimeRms());
1145 *fLog << " FADC sl. < " << Form("%2.1f", fArrTimeRmsLimit);
1146 *fLog << " in " << what << Form("%3i", cal.GetPixId()) << endl;
1147 bad.SetUncalibrated( MBadPixelsPix::kFluctuatingArrivalTimes);
1148 }
1149}
1150
1151// -----------------------------------------------------------------------------------
1152//
1153// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
1154// - MBadPixelsPix::kChargeIsPedestal
1155// - MBadPixelsPix::kChargeRelErrNotValid
1156// - MBadPixelsPix::kMeanTimeInFirstBin
1157// - MBadPixelsPix::kMeanTimeInLast2Bins
1158// - MBadPixelsPix::kDeviatingNumPhes
1159// - MBadPixelsPix::kHiGainOverFlow
1160// - MBadPixelsPix::kLoGainOverFlow
1161//
1162// - Call MCalibrationPix::SetExcluded() for the bad pixels
1163//
1164// Sets pixel to MBadPixelsPix::kUnreliableRun, if one of the following flags is set:
1165// - MBadPixelsPix::kChargeSigmaNotValid
1166//
1167void MCalibrationChargeCalc::FinalizeBadPixels()
1168{
1169
1170 MBadPixelsCam *badcam = fIntensBad ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
1171
1172 for (Int_t i=0; i<badcam->GetSize(); i++)
1173 {
1174
1175 MBadPixelsPix &bad = (*badcam)[i];
1176
1177 if (IsCheckDeadPixels())
1178 {
1179 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
1180 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1181 }
1182
1183 if (IsCheckExtractionWindow())
1184 {
1185 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
1186 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1187
1188 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
1189 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1190 }
1191
1192 if (IsCheckDeviatingBehavior())
1193 {
1194 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
1195 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1196 }
1197
1198 if (IsCheckHistOverflow())
1199 {
1200 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOverFlow ))
1201 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1202
1203 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOverFlow ))
1204 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1205 }
1206
1207 if (IsCheckArrivalTimes())
1208 {
1209 if (bad.IsUncalibrated( MBadPixelsPix::kFluctuatingArrivalTimes ))
1210 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1211 }
1212
1213 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
1214 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
1215 }
1216}
1217
1218// ------------------------------------------------------------------------
1219//
1220//
1221// First loop: Calculate a mean and mean RMS of photo-electrons per area index
1222// Include only pixels which are not MBadPixelsPix::kUnsuitableRun nor
1223// MBadPixelsPix::kChargeSigmaNotValid (see FinalizeBadPixels()) and set
1224// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
1225//
1226// Second loop: Get mean number of photo-electrons and its RMS including
1227// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
1228// and further exclude those deviating by more than fPheErrLimit mean
1229// sigmas from the mean (obtained in first loop). Set
1230// MBadPixelsPix::kDeviatingNumPhes if excluded.
1231//
1232// For the suitable pixels with flag MBadPixelsPix::kChargeSigmaNotValid
1233// set the number of photo-electrons as the mean number of photo-electrons
1234// calculated in that area index.
1235//
1236// Set weighted mean and variance of photo-electrons per area index in:
1237// average area pixels of MCalibrationChargeCam (obtained from:
1238// MCalibrationChargeCam::GetAverageArea() )
1239//
1240// Set weighted mean and variance of photo-electrons per sector in:
1241// average sector pixels of MCalibrationChargeCam (obtained from:
1242// MCalibrationChargeCam::GetAverageSector() )
1243//
1244//
1245// Third loop: Set mean number of photo-electrons and its RMS in the pixels
1246// only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1247//
1248Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
1249{
1250 MBadPixelsCam *badcam = fIntensBad
1251 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1252 MCalibrationChargeCam *chargecam = fIntensCam
1253 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1254
1255 const Int_t npixels = fGeom->GetNumPixels();
1256 const Int_t nareas = fGeom->GetNumAreas();
1257 const Int_t nsectors = fGeom->GetNumSectors();
1258
1259 TArrayF lowlim (nareas);
1260 TArrayF upplim (nareas);
1261 TArrayD areavars (nareas);
1262 TArrayD areaweights (nareas);
1263 TArrayD sectorweights (nsectors);
1264 TArrayD areaphes (nareas);
1265 TArrayD sectorphes (nsectors);
1266 TArrayI numareavalid (nareas);
1267 TArrayI numsectorvalid(nsectors);
1268
1269 //
1270 // First loop: Get mean number of photo-electrons and the RMS
1271 // The loop is only to recognize later pixels with very deviating numbers
1272 //
1273 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
1274
1275 for (Int_t i=0; i<npixels; i++)
1276 {
1277 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1278 MBadPixelsPix &bad = (*badcam)[i];
1279
1280 if (!pix.IsFFactorMethodValid())
1281 continue;
1282
1283 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1284 {
1285 pix.SetFFactorMethodValid(kFALSE);
1286 continue;
1287 }
1288
1289 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1290 continue;
1291
1292 if (!IsUseUnreliables())
1293 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1294 continue;
1295
1296 const Float_t nphe = pix.GetPheFFactorMethod();
1297 const Int_t aidx = (*fGeom)[i].GetAidx();
1298 camphes.Fill(i,nphe);
1299 camphes.SetUsed(i);
1300 areaphes [aidx] += nphe;
1301 areavars [aidx] += nphe*nphe;
1302 numareavalid[aidx] ++;
1303 }
1304
1305 for (Int_t i=0; i<nareas; i++)
1306 {
1307 if (numareavalid[i] == 0)
1308 {
1309 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
1310 << "in area index: " << i << endl;
1311 continue;
1312 }
1313
1314 if (numareavalid[i] == 1)
1315 areavars[i] = 0.;
1316 else if (numareavalid[i] == 0)
1317 {
1318 areaphes[i] = -1.;
1319 areaweights[i] = -1.;
1320 }
1321 else
1322 {
1323 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1);
1324 areaphes[i] = areaphes[i] / numareavalid[i];
1325 }
1326
1327 if (areavars[i] < 0.)
1328 {
1329 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of photo-electrons found "
1330 << "in area index: " << i << endl;
1331 continue;
1332 }
1333
1334 lowlim[i] = areaphes[i] - fPheErrLowerLimit*TMath::Sqrt(areavars[i]);
1335 upplim[i] = areaphes[i] + fPheErrUpperLimit*TMath::Sqrt(areavars[i]);
1336
1337
1338 TH1D *hist = camphes.ProjectionS(TArrayI(),TArrayI(1,&i),"_py",100);
1339 hist->Fit("gaus","Q0");
1340 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1341 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1342 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1343
1344 if (IsDebug())
1345 hist->DrawClone();
1346
1347 if (ndf < 5)
1348 {
1349 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons " << endl;
1350 *fLog << " in the camera with area index: " << i << endl;
1351 *fLog << " Number of dof.: " << ndf << " is smaller than 5 " << endl;
1352 *fLog << " Will use the simple mean and rms " << endl;
1353 delete hist;
1354 SetPheFitOK(i,kFALSE);
1355 continue;
1356 }
1357
1358 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1359
1360 if (prob < 0.001)
1361 {
1362 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons " << endl;
1363 *fLog << " in the camera with area index: " << i << endl;
1364 *fLog << " Fit probability " << prob << " is smaller than 0.001 " << endl;
1365 *fLog << " Will use the simple mean and rms " << endl;
1366 delete hist;
1367 SetPheFitOK(i,kFALSE);
1368 continue;
1369 }
1370
1371 if (mean < 0.)
1372 {
1373 *fLog << inf << GetDescriptor() << ": Fitted mean number of photo-electrons " << endl;
1374 *fLog << " with area idx " << i << ": " << mean << " is smaller than 0. " << endl;
1375 *fLog << warn << " Will use the simple mean and rms " << endl;
1376 SetPheFitOK(i,kFALSE);
1377 delete hist;
1378 continue;
1379 }
1380
1381 *fLog << inf << GetDescriptor() << ": Mean number of phes with area idx " << i << ": "
1382 << Form("%7.2f+-%6.2f",mean,sigma) << endl;
1383
1384 lowlim[i] = mean - fPheErrLowerLimit*sigma;
1385 upplim[i] = mean + fPheErrUpperLimit*sigma;
1386
1387 delete hist;
1388
1389 SetPheFitOK(i,kTRUE);
1390 }
1391
1392 *fLog << endl;
1393
1394 numareavalid.Reset();
1395 areaphes .Reset();
1396 areavars .Reset();
1397 //
1398 // Second loop: Get mean number of photo-electrons and its RMS excluding
1399 // pixels deviating by more than fPheErrLimit sigma.
1400 // Set the conversion factor FADC counts to photo-electrons
1401 //
1402 for (Int_t i=0; i<npixels; i++)
1403 {
1404
1405 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1406
1407 if (!pix.IsFFactorMethodValid())
1408 continue;
1409
1410 MBadPixelsPix &bad = (*badcam)[i];
1411
1412 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1413 continue;
1414
1415 if (!IsUseUnreliables())
1416 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1417 continue;
1418
1419 const Float_t nvar = pix.GetPheFFactorMethodVar();
1420 if (nvar <= 0.)
1421 {
1422 pix.SetFFactorMethodValid(kFALSE);
1423 continue;
1424 }
1425
1426 const Int_t aidx = (*fGeom)[i].GetAidx();
1427 const Int_t sector = (*fGeom)[i].GetSector();
1428 const Float_t area = (*fGeom)[i].GetA();
1429 const Float_t nphe = pix.GetPheFFactorMethod();
1430
1431 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1432 {
1433 *fLog << warn << "Number of phes: "
1434 << Form("%7.2f out of +%3.1f-%3.1f sigma limit: ",nphe,fPheErrUpperLimit,fPheErrLowerLimit)
1435 << Form("[%7.2f,%7.2f] pixel%4i",lowlim[aidx],upplim[aidx],i) << endl;
1436 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1437
1438 if (IsCheckDeviatingBehavior())
1439 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1440 continue;
1441 }
1442
1443 areaweights [aidx] += nphe*nphe;
1444 areaphes [aidx] += nphe;
1445 numareavalid [aidx] ++;
1446
1447 if (aidx == 0)
1448 fNumInnerFFactorMethodUsed++;
1449
1450 sectorweights [sector] += nphe*nphe/area/area;
1451 sectorphes [sector] += nphe/area;
1452 numsectorvalid[sector] ++;
1453 }
1454
1455 *fLog << endl;
1456
1457 for (Int_t aidx=0; aidx<nareas; aidx++)
1458 {
1459
1460 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1461
1462 if (numareavalid[aidx] == 1)
1463 areaweights[aidx] = 0.;
1464 else if (numareavalid[aidx] == 0)
1465 {
1466 areaphes[aidx] = -1.;
1467 areaweights[aidx] = -1.;
1468 }
1469 else
1470 {
1471 areaweights[aidx] = (areaweights[aidx] - areaphes[aidx]*areaphes[aidx]/numareavalid[aidx])
1472 / (numareavalid[aidx]-1);
1473 areaphes[aidx] /= numareavalid[aidx];
1474 }
1475
1476 if (areaweights[aidx] < 0. || areaphes[aidx] <= 0.)
1477 {
1478 *fLog << warn << GetDescriptor()
1479 << ": Mean number phes from area index " << aidx << " could not be calculated: "
1480 << " Mean: " << areaphes[aidx]
1481 << " Variance: " << areaweights[aidx] << endl;
1482 apix.SetFFactorMethodValid(kFALSE);
1483 continue;
1484 }
1485
1486 *fLog << inf << GetDescriptor()
1487 << ": Average total phes for area idx " << aidx << ": "
1488 << Form("%7.2f +- %6.2f",areaphes[aidx],TMath::Sqrt(areaweights[aidx])) << endl;
1489
1490 apix.SetPheFFactorMethod ( areaphes[aidx] );
1491 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1492 apix.SetFFactorMethodValid ( kTRUE );
1493
1494 }
1495
1496 *fLog << endl;
1497
1498 for (Int_t sector=0; sector<nsectors; sector++)
1499 {
1500
1501 if (numsectorvalid[sector] == 1)
1502 sectorweights[sector] = 0.;
1503 else if (numsectorvalid[sector] == 0)
1504 {
1505 sectorphes[sector] = -1.;
1506 sectorweights[sector] = -1.;
1507 }
1508 else
1509 {
1510 sectorweights[sector] = (sectorweights[sector]
1511 - sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector]
1512 )
1513 / (numsectorvalid[sector]-1.);
1514 sectorphes[sector] /= numsectorvalid[sector];
1515 }
1516
1517 MCalibrationChargePix &spix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
1518
1519 if (sectorweights[sector] < 0. || sectorphes[sector] <= 0.)
1520 {
1521 *fLog << warn << GetDescriptor()
1522 <<": Mean number phes/area for sector " << sector << " could not be calculated: "
1523 << " Mean: " << sectorphes[sector]
1524 << " Variance: " << sectorweights[sector] << endl;
1525 spix.SetFFactorMethodValid(kFALSE);
1526 continue;
1527 }
1528
1529 *fLog << inf << GetDescriptor()
1530 << ": Avg number phes/mm^2 in sector " << sector << ": "
1531 << Form("%5.3f+-%4.3f",sectorphes[sector],TMath::Sqrt(sectorweights[sector]))
1532 << endl;
1533
1534 spix.SetPheFFactorMethod ( sectorphes[sector] );
1535 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1536 spix.SetFFactorMethodValid ( kTRUE );
1537
1538 }
1539
1540 //
1541 // Third loop: Set mean number of photo-electrons and its RMS in the pixels
1542 // only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1543 //
1544 for (Int_t i=0; i<npixels; i++)
1545 {
1546
1547 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1548 MBadPixelsPix &bad = (*badcam)[i];
1549
1550 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1551 continue;
1552
1553 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1554 {
1555 const Int_t aidx = (*fGeom)[i].GetAidx();
1556 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1557
1558 pix.SetPheFFactorMethod ( apix.GetPheFFactorMethod() );
1559 pix.SetPheFFactorMethodVar( apix.GetPheFFactorMethodVar() );
1560
1561 if (!pix.CalcConvFFactor())
1562 {
1563 *fLog << warn << GetDescriptor()
1564 << ": Could not calculate the Conv. FADC counts to Phes in pixel: "
1565 << Form(" %4i",pix.GetPixId())
1566 << endl;
1567 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1568 if (IsCheckDeviatingBehavior())
1569 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1570 }
1571
1572 }
1573 }
1574
1575 return kTRUE;
1576}
1577
1578
1579
1580// ------------------------------------------------------------------------
1581//
1582// Returns kFALSE if pointer to MCalibrationBlindCam is NULL
1583//
1584// The check returns kFALSE if:
1585//
1586// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1587// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1588//
1589// Calls:
1590// - MCalibrationBlindPix::CalcFluxInsidePlexiglass()
1591//
1592Bool_t MCalibrationChargeCalc::FinalizeBlindCam()
1593{
1594
1595 MCalibrationBlindCam *blindcam = fIntensBlind
1596 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
1597
1598 if (!blindcam)
1599 return kFALSE;
1600
1601 Int_t nvalid = 0;
1602
1603 for (Int_t i=0; i<blindcam->GetSize(); i++)
1604 {
1605
1606 MCalibrationBlindPix &blindpix = (MCalibrationBlindPix&)(*blindcam)[i];
1607
1608 if (!blindpix.IsValid())
1609 continue;
1610
1611 const Float_t lambda = blindpix.GetLambda();
1612 const Float_t lambdaerr = blindpix.GetLambdaErr();
1613 const Float_t lambdacheck = blindpix.GetLambdaCheck();
1614
1615 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1616 {
1617 *fLog << warn << GetDescriptor() << ": Lambda="
1618 << Form("%4.2f", lambda) << " and Lambda-Check="
1619 << Form("%4.2f", lambdacheck) << " differ by more than "
1620 << Form("%4.2f", fLambdaCheckLimit) << " in the Blind Pixel Nr."
1621 << Form("%2i", i) << endl;
1622 blindpix.SetValid(kFALSE);
1623 continue;
1624 }
1625
1626 if (lambdaerr > fLambdaErrLimit)
1627 {
1628 *fLog << warn << GetDescriptor() << ": Error of Fitted Lambda="
1629 << Form("%4.2f", lambdaerr) << " is greater than "
1630 << Form("%4.2f", fLambdaErrLimit)
1631 << " in Blind Pixel Nr." << Form("%2d", i) << endl;
1632 blindpix.SetValid(kFALSE);
1633 continue;
1634 }
1635
1636 if (!blindpix.CalcFluxInsidePlexiglass())
1637 {
1638 *fLog << warn << "Could not calculate the flux of photons from Blind Pixel Nr." << i << endl;
1639 blindpix.SetValid(kFALSE);
1640 continue;
1641 }
1642
1643 nvalid++;
1644 }
1645
1646 if (!nvalid)
1647 return kFALSE;
1648
1649 return kTRUE;
1650}
1651
1652// ------------------------------------------------------------------------
1653//
1654// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1655//
1656// The check returns kFALSE if:
1657//
1658// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1659// 2) PINDiode has a fit error smaller than fChargeErrLimit
1660// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1661// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1662//
1663// Calls:
1664// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1665//
1666Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1667{
1668
1669 if (!fPINDiode)
1670 return kFALSE;
1671
1672 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1673 {
1674 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1675 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1676 return kFALSE;
1677 }
1678
1679 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1680 {
1681 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1682 << fChargeErrLimit << " in PINDiode " << endl;
1683 return kFALSE;
1684 }
1685
1686 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1687 {
1688 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1689 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1690 return kFALSE;
1691 }
1692
1693 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1694 {
1695 *fLog << warn << GetDescriptor()
1696 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1697 return kFALSE;
1698 }
1699
1700
1701 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1702 {
1703 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1704 << "will skip PIN Diode Calibration " << endl;
1705 return kFALSE;
1706 }
1707
1708 return kTRUE;
1709}
1710
1711// ------------------------------------------------------------------------
1712//
1713// Calculate the average number of photons outside the plexiglass with the
1714// formula:
1715//
1716// av.Num.photons(area index) = av.Num.Phes(area index)
1717// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1718// / MCalibrationQEPix::GetPMTCollectionEff()
1719// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1720// / MCalibrationQECam::GetPlexiglassQE()
1721//
1722// Calculate the variance on the average number of photons assuming that the error on the
1723// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1724// values keeps it ordinary error since it is systematic.
1725//
1726// Loop over pixels:
1727//
1728// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1729// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1730//
1731// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1732// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1733//
1734// - Calculate the quantum efficiency with the formula:
1735//
1736// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1737//
1738// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1739//
1740// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1741// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1742//
1743// - Call MCalibrationQEPix::UpdateFFactorMethod()
1744//
1745void MCalibrationChargeCalc::FinalizeFFactorQECam()
1746{
1747
1748 if (fNumInnerFFactorMethodUsed < 2)
1749 {
1750 *fLog << warn << GetDescriptor()
1751 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1752 return;
1753 }
1754
1755 MCalibrationQECam *qecam = fIntensQE
1756 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1757 MCalibrationChargeCam *chargecam = fIntensCam
1758 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1759 MBadPixelsCam *badcam = fIntensBad
1760 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1761
1762 MCalibrationChargePix &avpix = (MCalibrationChargePix&)chargecam->GetAverageArea(0);
1763 MCalibrationQEPix &qepix = (MCalibrationQEPix&) qecam->GetAverageArea(0);
1764
1765 if (IsDebug())
1766 *fLog << dbginf << "External Phes: " << fExternalNumPhes
1767 << " Internal Phes: " << avpix.GetPheFFactorMethod() << endl;
1768
1769 const Float_t avphotons = ( IsUseExternalNumPhes()
1770 ? fExternalNumPhes
1771 : avpix.GetPheFFactorMethod() )
1772 / qepix.GetDefaultQE(fPulserColor)
1773 / qepix.GetPMTCollectionEff()
1774 / qepix.GetLightGuidesEff(fPulserColor)
1775 / qecam->GetPlexiglassQE();
1776
1777 const Float_t avphotrelvar = ( IsUseExternalNumPhes()
1778 ? fExternalNumPhesRelVar
1779 : avpix.GetPheFFactorMethodRelVar() )
1780 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1781 + qepix.GetPMTCollectionEffRelVar()
1782 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1783 + qecam->GetPlexiglassQERelVar();
1784
1785 const UInt_t nareas = fGeom->GetNumAreas();
1786
1787 //
1788 // Set the results in the MCalibrationChargeCam
1789 //
1790 chargecam->SetNumPhotonsFFactorMethod (avphotons);
1791
1792 if (avphotrelvar > 0.)
1793 chargecam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1794
1795 TArrayF lowlim (nareas);
1796 TArrayF upplim (nareas);
1797 TArrayD avffactorphotons (nareas);
1798 TArrayD avffactorphotvar (nareas);
1799 TArrayI numffactor (nareas);
1800
1801 const UInt_t npixels = fGeom->GetNumPixels();
1802
1803 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1804
1805 for (UInt_t i=0; i<npixels; i++)
1806 {
1807
1808 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1809 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1810 MBadPixelsPix &bad = (*badcam) [i];
1811
1812 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1813 continue;
1814
1815 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1816 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1817
1818 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1819
1820 qepix.SetQEFFactor ( qe , fPulserColor );
1821 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1822 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1823
1824 if (!qepix.UpdateFFactorMethod( qecam->GetPlexiglassQE() ))
1825 *fLog << warn << GetDescriptor()
1826 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1827
1828 //
1829 // The following pixels are those with deviating sigma, but otherwise OK,
1830 // probably those with stars during the pedestal run, but not the cal. run.
1831 //
1832 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1833 {
1834 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1835 if (IsCheckDeviatingBehavior())
1836 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1837 continue;
1838 }
1839
1840 const Int_t aidx = (*fGeom)[i].GetAidx();
1841 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1842
1843 camffactor.Fill(i,ffactor);
1844 camffactor.SetUsed(i);
1845
1846 avffactorphotons[aidx] += ffactor;
1847 avffactorphotvar[aidx] += ffactor*ffactor;
1848 numffactor[aidx]++;
1849 }
1850
1851 for (UInt_t i=0; i<nareas; i++)
1852 {
1853
1854 if (numffactor[i] == 0)
1855 {
1856 *fLog << warn << GetDescriptor() << ": No pixels with valid total F-Factor found "
1857 << "in area index: " << i << endl;
1858 continue;
1859 }
1860
1861 avffactorphotvar[i] = (avffactorphotvar[i] - avffactorphotons[i]*avffactorphotons[i]/numffactor[i])
1862 / (numffactor[i]-1.);
1863 avffactorphotons[i] = avffactorphotons[i] / numffactor[i];
1864
1865 if (avffactorphotvar[i] < 0.)
1866 {
1867 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of total F-Factor found "
1868 << "in area index: " << i << endl;
1869 continue;
1870 }
1871
1872 lowlim [i] = 1.; // Lowest known F-Factor of a PMT
1873 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1874
1875 TArrayI area(1);
1876 area[0] = i;
1877
1878 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1879 hist->Fit("gaus","Q0");
1880 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1881 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1882 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1883
1884 if (IsDebug())
1885 camffactor.DrawClone();
1886
1887 if (ndf < 2)
1888 {
1889 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor " << endl;
1890 *fLog << " in the camera with area index: " << i << endl;
1891 *fLog << " Number of dof.: " << ndf << " is smaller than 2 " << endl;
1892 *fLog << " Will use the simple mean and rms." << endl;
1893 delete hist;
1894 SetFFactorFitOK(i,kFALSE);
1895 continue;
1896 }
1897
1898 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1899
1900 if (prob < 0.001)
1901 {
1902 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor " << endl;
1903 *fLog << " in the camera with area index: " << i << endl;
1904 *fLog << " Fit probability " << prob << " is smaller than 0.001 " << endl;
1905 *fLog << " Will use the simple mean and rms." << endl;
1906 delete hist;
1907 SetFFactorFitOK(i,kFALSE);
1908 continue;
1909 }
1910
1911 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1912 << "with area index #" << i << ": " << Form("%4.2f+-%4.2f",mean,sigma) << endl;
1913
1914 lowlim [i] = 1.;
1915 upplim [i] = mean + fFFactorErrLimit*sigma;
1916
1917 delete hist;
1918
1919 SetFFactorFitOK(i,kTRUE);
1920 }
1921
1922 *fLog << endl;
1923
1924 for (UInt_t i=0; i<npixels; i++)
1925 {
1926
1927 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1928 MBadPixelsPix &bad = (*badcam) [i];
1929
1930 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1931 continue;
1932
1933 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1934 const Int_t aidx = (*fGeom)[i].GetAidx();
1935
1936 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1937 {
1938 *fLog << warn << "Overall F-Factor "
1939 << Form("%5.2f",ffactor) << " out of range ["
1940 << Form("%5.2f,%5.2f",lowlim[aidx],upplim[aidx]) << "] Pixel " << i << endl;
1941
1942 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1943 if (IsCheckDeviatingBehavior())
1944 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1945 }
1946 }
1947
1948 for (UInt_t i=0; i<npixels; i++)
1949 {
1950
1951 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1952 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1953 MBadPixelsPix &bad = (*badcam) [i];
1954
1955 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1956 {
1957 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1958 pix.SetFFactorMethodValid(kFALSE);
1959 continue;
1960 }
1961 }
1962}
1963
1964
1965// ------------------------------------------------------------------------
1966//
1967// Loop over pixels:
1968//
1969// - Continue, if not MCalibrationBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1970// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1971//
1972// - Calculate the quantum efficiency with the formula:
1973//
1974// QE = Num.Phes / MCalibrationBlindPix::GetFluxInsidePlexiglass()
1975// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1976//
1977// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1978// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1979// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1980//
1981// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1982//
1983void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1984{
1985
1986
1987 MCalibrationBlindCam *blindcam = fIntensBlind
1988 ? (MCalibrationBlindCam*) fIntensBlind->GetCam(): fBlindCam;
1989 MBadPixelsCam *badcam = fIntensBad
1990 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1991 MCalibrationQECam *qecam = fIntensQE
1992 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1993 MCalibrationChargeCam *chargecam = fIntensCam
1994 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1995
1996 if (!blindcam)
1997 return;
1998
1999 //
2000 // Set the results in the MCalibrationChargeCam
2001 //
2002 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
2003 {
2004
2005 const Float_t photons = blindcam->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
2006 / qecam->GetPlexiglassQE();
2007 chargecam->SetNumPhotonsBlindPixelMethod(photons);
2008
2009 const Float_t photrelvar = blindcam->GetFluxInsidePlexiglassRelVar()
2010 + qecam->GetPlexiglassQERelVar();
2011
2012 if (photrelvar > 0.)
2013 chargecam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
2014 }
2015
2016 //
2017 // With the knowledge of the overall photon flux, calculate the
2018 // quantum efficiencies after the Blind Pixel and PIN Diode method
2019 //
2020 const UInt_t npixels = fGeom->GetNumPixels();
2021 for (UInt_t i=0; i<npixels; i++)
2022 {
2023
2024 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2025
2026 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
2027 {
2028 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
2029 continue;
2030 }
2031
2032 MBadPixelsPix &bad = (*badcam) [i];
2033
2034 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2035 {
2036 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
2037 continue;
2038 }
2039
2040 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
2041 MGeomPix &geo = (*fGeom) [i];
2042
2043 const Float_t qe = pix.GetPheFFactorMethod()
2044 / blindcam->GetFluxInsidePlexiglass()
2045 / geo.GetA()
2046 * qecam->GetPlexiglassQE();
2047
2048 const Float_t qerelvar = blindcam->GetFluxInsidePlexiglassRelVar()
2049 + qecam->GetPlexiglassQERelVar()
2050 + pix.GetPheFFactorMethodRelVar();
2051
2052 qepix.SetQEBlindPixel ( qe , fPulserColor );
2053 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
2054 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
2055
2056 if (!qepix.UpdateBlindPixelMethod( qecam->GetPlexiglassQE()))
2057 *fLog << warn << GetDescriptor()
2058 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
2059 }
2060}
2061
2062// ------------------------------------------------------------------------
2063//
2064// Loop over pixels:
2065//
2066// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
2067// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
2068//
2069// - Calculate the quantum efficiency with the formula:
2070//
2071// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
2072//
2073// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
2074// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
2075// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
2076//
2077// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
2078//
2079void MCalibrationChargeCalc::FinalizePINDiodeQECam()
2080{
2081
2082 const UInt_t npixels = fGeom->GetNumPixels();
2083
2084 MCalibrationQECam *qecam = fIntensQE
2085 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
2086 MCalibrationChargeCam *chargecam = fIntensCam
2087 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
2088 MBadPixelsCam *badcam = fIntensBad
2089 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
2090
2091 if (!fPINDiode)
2092 return;
2093
2094 //
2095 // With the knowledge of the overall photon flux, calculate the
2096 // quantum efficiencies after the PIN Diode method
2097 //
2098 for (UInt_t i=0; i<npixels; i++)
2099 {
2100
2101 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2102
2103 if (!fPINDiode)
2104 {
2105 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2106 continue;
2107 }
2108
2109 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
2110 {
2111 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2112 continue;
2113 }
2114
2115 MBadPixelsPix &bad = (*badcam) [i];
2116
2117 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2118 {
2119 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2120 continue;
2121 }
2122
2123 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
2124 MGeomPix &geo = (*fGeom) [i];
2125
2126 const Float_t qe = pix.GetPheFFactorMethod()
2127 / fPINDiode->GetFluxOutsidePlexiglass()
2128 / geo.GetA();
2129
2130 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
2131
2132 qepix.SetQEPINDiode ( qe , fPulserColor );
2133 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
2134 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
2135
2136 if (!qepix.UpdatePINDiodeMethod())
2137 *fLog << warn << GetDescriptor()
2138 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
2139 }
2140}
2141
2142// ------------------------------------------------------------------------
2143//
2144// Loop over pixels:
2145//
2146// - Call MCalibrationQEPix::UpdateCombinedMethod()
2147//
2148void MCalibrationChargeCalc::FinalizeCombinedQECam()
2149{
2150
2151 const UInt_t npixels = fGeom->GetNumPixels();
2152
2153 MCalibrationQECam *qecam = fIntensQE
2154 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
2155 MBadPixelsCam *badcam = fIntensBad
2156 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
2157
2158 for (UInt_t i=0; i<npixels; i++)
2159 {
2160
2161 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2162 MBadPixelsPix &bad = (*badcam) [i];
2163
2164 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2165 {
2166 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2167 continue;
2168 }
2169
2170 qepix.UpdateCombinedMethod();
2171 }
2172}
2173
2174// -----------------------------------------------------------------------------------------------
2175//
2176// - Print out statistics about BadPixels of type UnsuitableType_t
2177// - store numbers of bad pixels of each type in fCam or fIntensCam
2178//
2179Bool_t MCalibrationChargeCalc::FinalizeUnsuitablePixels()
2180{
2181
2182 *fLog << inf << endl;
2183 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
2184 *fLog << dec << setfill(' ');
2185
2186 const Int_t nareas = fGeom->GetNumAreas();
2187
2188 TArrayI suit(nareas);
2189 TArrayI unsuit(nareas);
2190 TArrayI unrel(nareas);
2191
2192 const MBadPixelsCam *badcam = fIntensBad ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2193 MCalibrationChargeCam *chargecam = fIntensCam ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
2194
2195 Int_t unsuitcnt=0;
2196 Int_t unrelcnt =0;
2197
2198 // Count number of succesfully calibrated pixels
2199 for (Int_t aidx=0; aidx<nareas; aidx++)
2200 {
2201 suit[aidx] = badcam->GetNumSuitable(MBadPixelsPix::kUnsuitableRun, fGeom, aidx);
2202 unsuit[aidx] = badcam->GetNumUnsuitable(MBadPixelsPix::kUnsuitableRun, fGeom, aidx);
2203 unrel[aidx] = badcam->GetNumUnsuitable(MBadPixelsPix::kUnreliableRun, fGeom, aidx);
2204
2205 unsuitcnt += unsuit[aidx];
2206 unrelcnt += unrel[aidx];
2207
2208 chargecam->SetNumUnsuitable(unsuit[aidx], aidx);
2209 chargecam->SetNumUnreliable(unrel[aidx], aidx);
2210 }
2211
2212 TArrayI counts(nareas);
2213 for (Int_t i=0; i<chargecam->GetSize(); i++)
2214 {
2215 MCalibrationPix &pix = (*chargecam)[i];
2216 if (pix.IsHiGainSaturation())
2217 {
2218 const Int_t aidx = (*fGeom)[i].GetAidx();
2219 counts[aidx]++;
2220 }
2221 }
2222
2223 if (fGeom->InheritsFrom("MGeomCamMagic"))
2224 {
2225 *fLog << " Successfully calibrated Pixels: Inner: "
2226 << Form("%3i",suit[0]) << " Outer: " << Form("%3i",suit[1]) << endl;
2227 *fLog << " Uncalibrated Pixels: Inner: "
2228 << Form("%3i",unsuit[0]) << " Outer: " << Form("%3i",unsuit[1]) << endl;
2229 *fLog << " Unreliable Pixels: Inner: "
2230 << Form("%3i",unrel[0]) << " Outer: " << Form("%3i",unrel[1]) << endl;
2231 *fLog << " High-gain saturated Pixels: Inner: "
2232 << Form("%3i",counts[0]) << " Outer: " << Form("%3i",counts[1]) << endl;
2233 *fLog << endl;
2234 }
2235
2236 if (unsuitcnt > fUnsuitablesLimit*fGeom->GetNumPixels())
2237 {
2238 *fLog << err << "Number of unsuitable pixels: " << 100.*unsuitcnt/fGeom->GetNumPixels()
2239 << "% exceeds limit of " << fUnsuitablesLimit*100 << "%" << endl;
2240 return kFALSE;
2241 }
2242
2243 if (unrelcnt > fUnreliablesLimit*fGeom->GetNumPixels())
2244 {
2245 *fLog << err << "Relative number of unreliable pixels: " << 100.*unrelcnt/fGeom->GetNumPixels()
2246 << "% exceeds limit of " << fUnreliablesLimit*100 << "%" << endl;
2247 return kFALSE;
2248 }
2249 return kTRUE;
2250}
2251
2252// -----------------------------------------------------------------------------------------------
2253//
2254// Print out statistics about BadPixels of type UncalibratedType_t
2255//
2256void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
2257{
2258
2259 UInt_t countinner = 0;
2260 UInt_t countouter = 0;
2261
2262 MBadPixelsCam *badcam = fIntensBad
2263 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2264
2265 for (Int_t i=0; i<badcam->GetSize(); i++)
2266 {
2267 MBadPixelsPix &bad = (*badcam)[i];
2268
2269 if (bad.IsUncalibrated(typ))
2270 {
2271 if (fGeom->GetPixRatio(i) == 1.)
2272 countinner++;
2273 else
2274 countouter++;
2275 }
2276 }
2277
2278 *fLog << " " << setw(7) << text << "Inner: " << Form("%3i",countinner)
2279 << " Outer: " << Form("%3i", countouter) << endl;
2280}
2281
2282// --------------------------------------------------------------------------
2283//
2284// Set the path for output file
2285//
2286void MCalibrationChargeCalc::SetOutputPath(TString path)
2287{
2288 fOutputPath = path;
2289 if (fOutputPath.EndsWith("/"))
2290 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
2291}
2292
2293// --------------------------------------------------------------------------
2294//
2295// Set the output file
2296//
2297void MCalibrationChargeCalc::SetOutputFile(TString file)
2298{
2299 fOutputFile = file;
2300}
2301
2302// --------------------------------------------------------------------------
2303//
2304// Get the output file
2305//
2306const char* MCalibrationChargeCalc::GetOutputFile()
2307{
2308 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
2309}
2310
2311// --------------------------------------------------------------------------
2312//
2313// Read the environment for the following data members:
2314// - fChargeLimit
2315// - fChargeErrLimit
2316// - fChargeRelErrLimit
2317// - fDebug
2318// - fFFactorErrLimit
2319// - fLambdaErrLimit
2320// - fLambdaCheckErrLimit
2321// - fPheErrLimit
2322//
2323Int_t MCalibrationChargeCalc::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
2324{
2325
2326 Bool_t rc = kFALSE;
2327 if (IsEnvDefined(env, prefix, "ChargeLimit", print))
2328 {
2329 SetChargeLimit(GetEnvValue(env, prefix, "ChargeLimit", fChargeLimit));
2330 rc = kTRUE;
2331 }
2332 if (IsEnvDefined(env, prefix, "ChargeErrLimit", print))
2333 {
2334 SetChargeErrLimit(GetEnvValue(env, prefix, "ChargeErrLimit", fChargeErrLimit));
2335 rc = kTRUE;
2336 }
2337 if (IsEnvDefined(env, prefix, "ChargeRelErrLimit", print))
2338 {
2339 SetChargeRelErrLimit(GetEnvValue(env, prefix, "ChargeRelErrLimit", fChargeRelErrLimit));
2340 rc = kTRUE;
2341 }
2342 if (IsEnvDefined(env, prefix, "Debug", print))
2343 {
2344 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
2345 rc = kTRUE;
2346 }
2347 if (IsEnvDefined(env, prefix, "ArrTimeRmsLimit", print))
2348 {
2349 SetArrTimeRmsLimit(GetEnvValue(env, prefix, "ArrTimeRmsLimit", fArrTimeRmsLimit));
2350 rc = kTRUE;
2351 }
2352 if (IsEnvDefined(env, prefix, "FFactorErrLimit", print))
2353 {
2354 SetFFactorErrLimit(GetEnvValue(env, prefix, "FFactorErrLimit", fFFactorErrLimit));
2355 rc = kTRUE;
2356 }
2357 if (IsEnvDefined(env, prefix, "LambdaErrLimit", print))
2358 {
2359 SetLambdaErrLimit(GetEnvValue(env, prefix, "LambdaErrLimit", fLambdaErrLimit));
2360 rc = kTRUE;
2361 }
2362 if (IsEnvDefined(env, prefix, "LambdaCheckLimit", print))
2363 {
2364 SetLambdaCheckLimit(GetEnvValue(env, prefix, "LambdaCheckLimit", fLambdaCheckLimit));
2365 rc = kTRUE;
2366 }
2367 if (IsEnvDefined(env, prefix, "CheckDeadPixels", print))
2368 {
2369 SetCheckDeadPixels(GetEnvValue(env, prefix, "CheckDeadPixels", IsCheckDeadPixels()));
2370 rc = kTRUE;
2371 }
2372 if (IsEnvDefined(env, prefix, "CheckDeviatingBehavior", print))
2373 {
2374 SetCheckDeviatingBehavior(GetEnvValue(env, prefix, "CheckDeviatingBehavior", IsCheckDeviatingBehavior()));
2375 rc = kTRUE;
2376 }
2377 if (IsEnvDefined(env, prefix, "CheckExtractionWindow", print))
2378 {
2379 SetCheckExtractionWindow(GetEnvValue(env, prefix, "CheckExtractionWindow", IsCheckExtractionWindow()));
2380 rc = kTRUE;
2381 }
2382 if (IsEnvDefined(env, prefix, "CheckHistOverflow", print))
2383 {
2384 SetCheckHistOverflow(GetEnvValue(env, prefix, "CheckHistOverflow", IsCheckHistOverflow()));
2385 rc = kTRUE;
2386 }
2387 if (IsEnvDefined(env, prefix, "CheckOscillations", print))
2388 {
2389 SetCheckOscillations(GetEnvValue(env, prefix, "CheckOscillations", IsCheckOscillations()));
2390 rc = kTRUE;
2391 }
2392 if (IsEnvDefined(env, prefix, "CheckArrivalTimes", print))
2393 {
2394 SetCheckArrivalTimes(GetEnvValue(env, prefix, "CheckArrivalTimes", IsCheckArrivalTimes()));
2395 rc = kTRUE;
2396 }
2397 if (IsEnvDefined(env, prefix, "PheErrLowerLimit", print))
2398 {
2399 SetPheErrLowerLimit(GetEnvValue(env, prefix, "PheErrLowerLimit", fPheErrLowerLimit));
2400 rc = kTRUE;
2401 }
2402 if (IsEnvDefined(env, prefix, "PheErrUpperLimit", print))
2403 {
2404 SetPheErrUpperLimit(GetEnvValue(env, prefix, "PheErrUpperLimit", fPheErrUpperLimit));
2405 rc = kTRUE;
2406 }
2407 if (IsEnvDefined(env, prefix, "UseExtractorRes", print))
2408 {
2409 SetUseExtractorRes(GetEnvValue(env, prefix, "UseExtractorRes", IsUseExtractorRes()));
2410 rc = kTRUE;
2411 }
2412 if (IsEnvDefined(env, prefix, "UseUnreliables", print))
2413 {
2414 SetUseUnreliables(GetEnvValue(env, prefix, "UseUnreliables", IsUseUnreliables()));
2415 rc = kTRUE;
2416 }
2417
2418 if (IsEnvDefined(env, prefix, "UseExternalNumPhes", print))
2419 {
2420 SetUseExternalNumPhes(GetEnvValue(env, prefix, "UseExternalNumPhes", IsUseExternalNumPhes()));
2421 rc = kTRUE;
2422 }
2423
2424 if (IsEnvDefined(env, prefix, "UnsuitablesLimit", print))
2425 {
2426 SetUnsuitablesLimit(GetEnvValue(env, prefix, "UnsuitablesLimit", fUnsuitablesLimit));
2427 rc = kTRUE;
2428 }
2429
2430 if (IsEnvDefined(env, prefix, "UnreliablesLimit", print))
2431 {
2432 SetUnreliablesLimit(GetEnvValue(env, prefix, "UnreliablesLimit", fUnreliablesLimit));
2433 rc = kTRUE;
2434 }
2435
2436
2437 return rc;
2438}
2439
2440
2441
Note: See TracBrowser for help on using the repository browser.