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

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