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

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