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

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