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

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