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

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