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

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