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

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