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

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