source: tags/Mars-V0.9.3/mcalib/MCalibrationChargeCalc.cc

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