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

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