source: tags/Mars-V0.8.7pre/mcalib/MCalibrationChargeCalc.cc

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