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

Last change on this file since 6013 was 6013, checked in by gaug, 20 years ago
*** empty log message ***
File size: 75.3 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24//////////////////////////////////////////////////////////////////////////////
25//
26// MCalibrationChargeCalc
27//
28// Task to calculate the calibration conversion factors and quantum efficiencies
29// from the fit results to the summed FADC slice distributions delivered by
30// MCalibrationChargeCam, MCalibrationChargePix, MCalibrationChargeBlindPix and
31// MCalibrationChargePINDiode, calculated and filled by MHCalibrationChargeCam,
32// MHCalibrationChargePix, MHCalibrationChargeBlindPix and MHCalibrationChargePINDiode.
33//
34// PreProcess(): Initialize pointers to MCalibrationChargeCam, MCalibrationChargeBlindPix
35// MCalibrationChargePINDiode and MCalibrationQECam
36//
37// Initialize pulser light wavelength
38//
39// ReInit(): MCalibrationCam::InitSize(NumPixels) is called from MGeomApply (which allocates
40// memory in a TClonesArray of type MCalibrationChargePix)
41// Initializes pointer to MBadPixelsCam
42//
43// Process(): Nothing to be done, histograms getting filled by MHCalibrationChargeCam
44//
45// PostProcess(): - FinalizePedestals()
46// - FinalizeCharges()
47// - FinalizeFFactorMethod()
48// - FinalizeBadPixels()
49// - FinalizeBlindCam()
50// - FinalizePINDiode()
51// - FinalizeFFactorQECam()
52// - FinalizeBlindPixelQECam()
53// - FinalizePINDiodeQECam()
54//
55// Input Containers:
56// MCalibrationChargeCam
57// MCalibrationChargeBlindPix
58// MCalibrationChargePINDiode
59// MCalibrationQECam
60// MPedestalCam
61// MBadPixelsCam
62// MGeomCam
63// MTime
64//
65// Output Containers:
66// MCalibrationChargeCam
67// MCalibrationChargeBlindPix
68// MCalibrationChargePINDiode
69// MCalibrationQECam
70// MBadPixelsCam
71//
72//
73// Preliminary description of the calibration in photons (email from 12/02/04)
74//
75// Why calibrating in photons:
76// ===========================
77//
78// At the Barcelona meeting in 2002, we decided to calibrate the camera in
79// photons. This for the following reasons:
80//
81// * The physical quantity arriving at the camera are photons. This is
82// the direct physical information from the air shower. The photons
83// have a flux and a spectrum.
84//
85// * The photon fluxes depend mostly on the shower energy (with
86// corrections deriving from the observation conditions), while the photon
87// spectra depend mostly on the observation conditions: zenith angle,
88// quality of the air, also the impact parameter of the shower.
89//
90// * The photomultiplier, in turn, has different response properties
91// (quantum efficiencies) for photons of different colour. (Moreover,
92// different pixels have slightly different quantum efficiencies).
93// The resulting number of photo-electrons is then amplified (linearly)
94// with respect to the photo-electron flux.
95//
96// * In the ideal case, one would like to disentagle the effects
97// of the observation conditions from the primary particle energy (which
98// one likes to measure). To do so, one needs:
99//
100// 1) A reliable calibration relating the FADC counts to the photo-electron
101// flux -> This is accomplished with the F-Factor method.
102//
103// 2) A reliable calibration of the wavelength-dependent quantum efficiency
104// -> This is accomplished with the combination of the three methods,
105// together with QE-measurements performed by David in order to do
106// the interpolation.
107//
108// 3) A reliable calibration of the observation conditions. This means:
109// - Tracing the atmospheric conditions -> LIDAR
110// - Tracing the observation zenith angle -> Drive System
111//
112// 4) Some knowlegde about the impact parameter:
113// - This is the only part which cannot be accomplished well with a
114// single telescope. We would thus need to convolute the spectrum
115// over the distribution of impact parameters.
116//
117//
118// How an ideal calibration would look like:
119// =========================================
120//
121// We know from the combined PIN-Diode and Blind-Pixel Method the response of
122// each pixel to well-measured light fluxes in three representative
123// wavelengths (green, blue, UV). We also know the response to these light
124// fluxes in photo-electrons. Thus, we can derive:
125//
126// - conversion factors to photo-electrons
127// - conversion factors to photons in three wavelengths.
128//
129// Together with David's measurements and some MC-simulation, we should be
130// able to derive tables for typical Cherenkov-photon spectra - convoluted
131// with the impact parameters and depending on the athmospheric conditions
132// and the zenith angle (the "outer parameters").
133//
134// From these tables we can create "calibration tables" containing some
135// effective quantum efficiency depending on these outer parameters and which
136// are different for each pixel.
137//
138// In an ideal MCalibrate, one would thus have to convert first the FADC
139// slices to Photo-electrons and then, depending on the outer parameters,
140// look up the effective quantum efficiency and get the mean number of
141// photons which is then used for the further analysis.
142//
143// How the (first) MAGIC calibration should look like:
144// ===================================================
145//
146// For the moment, we have only one reliable calibration method, although
147// with very large systematic errors. This is the F-Factor method. Knowing
148// that the light is uniform over the whole camera (which I would not at all
149// guarantee in the case of the CT1 pulser), one could in principle already
150// perform a relative calibration of the quantum efficiencies in the UV.
151// However, the spread in QE at UV is about 10-15% (according to the plot
152// that Abelardo sent around last time. The spread in photo-electrons is 15%
153// for the inner pixels, but much larger (40%) for the outer ones.
154//
155// I'm not sure if we can already say that we have measured the relative
156// difference in quantum efficiency for the inner pixels and produce a first
157// QE-table for each pixel. To so, I would rather check in other wavelengths
158// (which we can do in about one-two weeks when the optical transmission of
159// the calibration trigger is installed).
160//
161// Thus, for the moment being, I would join Thomas proposal to calibrate in
162// photo-electrons and apply one stupid average quantum efficiency for all
163// pixels. This keeping in mind that we will have much preciser information
164// in about one to two weeks.
165//
166//
167// What MCalibrate should calculate and what should be stored:
168// ===========================================================
169//
170// It is clear that in the end, MCerPhotEvt will store photons.
171// MCalibrationCam stores the conversionfactors to photo-electrons and also
172// some tables of how to apply the conversion to photons, given the outer
173// parameters. This is not yet implemented and not even discussed.
174//
175// To start, I would suggest that we define the "average quantum efficiency"
176// (maybe something like 25+-3%) and apply them equally to all
177// photo-electrons. Later, this average factor can be easily replaced by a
178// pixel-dependent factor and later by a (pixel-dependent) table.
179//
180//
181//
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 fBlindPixelFlags.Set(0);
343 fPINDiodeFlags .Set(0);
344 fResultFlags .Set(0);
345}
346
347
348// -----------------------------------------------------------------------------------
349//
350// The following container are searched for and execution aborted if not in MParList:
351// - MPedestalCam
352// - MCalibrationPattern
353// - MExtractedSignalCam
354//
355Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
356{
357
358 /*
359 if (IsInterlaced())
360 {
361 fTrigPattern = (MTriggerPattern*)pList->FindObject("MTriggerPattern");
362 if (!fTrigPattern)
363 {
364 *fLog << err << "MTriggerPattern not found... abort." << endl;
365 return kFALSE;
366 }
367 }
368 */
369
370 fCalibPattern = (MCalibrationPattern*)pList->FindObject("MCalibrationPattern");
371 if (!fCalibPattern)
372 {
373 *fLog << err << "MCalibrationPattern not found... abort." << endl;
374 return kFALSE;
375 }
376
377 //
378 // Containers that have to be there.
379 //
380 fSignal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
381 if (!fSignal)
382 {
383 *fLog << err << "MExtractedSignalCam not found... aborting" << endl;
384 return kFALSE;
385 }
386
387 if (fPedestals)
388 return kTRUE;
389
390 fPedestals = (MPedestalCam*)pList->FindObject(AddSerialNumber(fNamePedestalCam), "MPedestalCam");
391 if (!fPedestals)
392 {
393 *fLog << err << GetDescriptor() << ": " << fNamePedestalCam << " not found... aborting" << endl;
394 return kFALSE;
395 }
396
397 return kTRUE;
398}
399
400
401// --------------------------------------------------------------------------
402//
403// Search for the following input containers and abort if not existing:
404// - MGeomCam
405// - MCalibrationIntensityChargeCam or MCalibrationChargeCam
406// - MCalibrationIntensityQECam or MCalibrationQECam
407// - MBadPixelsIntensityCam or MBadPixelsCam
408//
409// Search for the following input containers and give a warning if not existing:
410// - MCalibrationBlindPix
411// - MCalibrationChargePINDiode
412//
413// It retrieves the following variables from MCalibrationChargeCam:
414//
415// - fNumHiGainSamples
416// - fNumLoGainSamples
417//
418// It defines the PixId of every pixel in:
419//
420// - MCalibrationIntensityChargeCam
421// - MCalibrationChargeCam
422// - MCalibrationQECam
423//
424// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
425//
426// - MCalibrationChargePix
427// - MCalibrationQEPix
428//
429// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
430//
431// - MCalibrationChargeCam
432// - MCalibrationBlindPix (if existing)
433// - MCalibrationChargePINDiode (if existing)
434//
435Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
436{
437
438 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
439 if (!fGeom)
440 {
441 *fLog << err << "No MGeomCam found... aborting." << endl;
442 return kFALSE;
443 }
444
445 fIntensCam = (MCalibrationIntensityChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityChargeCam"));
446 if (fIntensCam)
447 *fLog << inf << "Found MCalibrationIntensityChargeCam ... " << endl;
448 else
449 {
450 fCam = (MCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
451 if (!fCam)
452 {
453 *fLog << err << "Cannot find MCalibrationChargeCam ... abort." << endl;
454 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
455 return kFALSE;
456 }
457 }
458
459 fHCam = (MHCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeCam"));
460 if (!fHCam)
461 {
462 *fLog << err << "Cannot find MHCalibrationChargeCam ... abort." << endl;
463 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
464 return kFALSE;
465 }
466
467 fIntensQE = (MCalibrationIntensityQECam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityQECam"));
468 if (fIntensQE)
469 *fLog << inf << "Found MCalibrationIntensityQECam ... " << endl;
470 else
471 {
472 fQECam = (MCalibrationQECam*)pList->FindObject(AddSerialNumber("MCalibrationQECam"));
473 if (!fQECam)
474 {
475 *fLog << err << "Cannot find MCalibrationQECam ... abort." << endl;
476 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationQECam before..." << endl;
477 return kFALSE;
478 }
479 }
480
481 fIntensBlind = (MCalibrationIntensityBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityBlindCam"));
482 if (fIntensBlind)
483 *fLog << inf << "Found MCalibrationIntensityBlindCam ... " << endl;
484 else
485 {
486 fBlindCam = (MCalibrationBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationBlindCam"));
487 if (!fBlindCam)
488 {
489 *fLog << endl;
490 *fLog << warn << GetDescriptor()
491 << ": No MCalibrationBlindCam found... no Blind Pixel method! " << endl;
492 }
493 }
494
495 fHBlindCam = (MHCalibrationChargeBlindCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeBlindCam"));
496 if (!fHBlindCam)
497 {
498 *fLog << endl;
499 *fLog << warn << GetDescriptor()
500 << ": No MHCalibrationChargeBlindCam found... no Blind Pixel method! " << endl;
501 }
502
503 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
504 if (fIntensBad)
505 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
506 else
507 {
508 fBadPixels = (MBadPixelsCam*)pList->FindObject(AddSerialNumber("MBadPixelsCam"));
509 if (!fBadPixels)
510 {
511 *fLog << err << "Cannot find MBadPixelsCam ... abort." << endl;
512 return kFALSE;
513 }
514 }
515
516 //
517 // Optional Containers
518 //
519 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
520 if (!fPINDiode)
521 {
522 *fLog << endl;
523 *fLog << warn << GetDescriptor()
524 << ": MCalibrationChargePINDiode not found... no PIN Diode method! " << endl;
525 }
526
527 MCalibrationQECam *qecam = fIntensQE
528 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
529 MCalibrationChargeCam *chargecam = fIntensCam
530 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
531 MBadPixelsCam *badcam = fIntensBad
532 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
533
534 UInt_t npixels = fGeom->GetNumPixels();
535
536 for (UInt_t i=0; i<npixels; i++)
537 {
538
539 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
540 MCalibrationQEPix &pqe = (MCalibrationQEPix&) (*qecam) [i];
541 MBadPixelsPix &bad = (*badcam) [i];
542
543 if (bad.IsBad())
544 {
545 pix.SetExcluded();
546 pqe.SetExcluded();
547 continue;
548 }
549
550 if (IsDebug())
551 pix.SetDebug();
552 }
553
554 return kTRUE;
555}
556
557// ----------------------------------------------------------------------------------
558//
559// Set the correct colour to the charge containers
560//
561Int_t MCalibrationChargeCalc::Process()
562{
563
564 /*
565 if (IsInterlaced())
566 {
567 if (!(fTrigPattern->GetPrescaled() & MTriggerPattern::kCalibration))
568 {
569 if (IsDebug())
570 {
571 for (Int_t i=16; i>= 0; i--)
572 *fLog << err << (fTrigPattern->GetPrescaled() >> i & 1);
573 *fLog << endl;
574 }
575 return kTRUE;
576 }
577 else
578 {
579 if (IsDebug())
580 {
581 for (Int_t i=16; i>= 0; i--)
582 *fLog << inf << (fTrigPattern->GetPrescaled() >> i & 1);
583 *fLog << endl;
584 }
585 }
586 }
587 */
588 const MCalibrationCam::PulserColor_t col = fCalibPattern->GetPulserColor();
589
590 if (col == fPulserColor)
591 {
592 fNumProcessed++;
593 return kTRUE;
594 }
595
596 if (col == MCalibrationCam::kNONE)
597 return kTRUE;
598
599
600 //
601 // Now retrieve the colour and check if not various colours have been used
602 //
603 if (!fIntensCam)
604 {
605 if (fPulserColor != MCalibrationCam::kNONE)
606 {
607 *fLog << warn << GetDescriptor()
608 << "Multiple colours used simultaneously!" ;
609 fHCam->Finalize();
610 fHBlindCam->Finalize();
611
612 Finalize();
613
614 fHCam->ResetHists();
615 fHBlindCam->ResetHists();
616
617 *fLog << inf << "Starting next calibration... " << flush;
618
619 fHCam->SetColor(col);
620 fHBlindCam->SetColor(col);
621
622 fCam->SetPulserColor(col);
623 if (fBlindCam)
624 fBlindCam->SetPulserColor(col);
625 }
626 }
627
628 fPulserColor = col;
629
630 *fLog << inf << "Found new colour ... " << flush;
631
632 switch (col)
633 {
634 case MCalibrationCam::kGREEN: *fLog << "Green"; break;
635 case MCalibrationCam::kBLUE: *fLog << "Blue"; break;
636 case MCalibrationCam::kUV: *fLog << "UV"; break;
637 case MCalibrationCam::kCT1: *fLog << "CT1"; break;
638 default: break;
639 }
640
641 *fLog << inf << " with strength: " << fCalibPattern->GetPulserStrength() << endl;
642
643 fHCam->SetColor(col);
644 fHBlindCam->SetColor(col);
645
646 fCam->SetPulserColor(col);
647
648 if (fBlindCam)
649 fBlindCam->SetPulserColor(col);
650 if (fPINDiode)
651 fPINDiode->SetColor(col);
652
653 fNumProcessed = 0;
654
655 return kTRUE;
656}
657
658// -----------------------------------------------------------------------
659//
660// Return if number of executions is null.
661//
662Int_t MCalibrationChargeCalc::PostProcess()
663{
664
665 if (GetNumExecutions()==0)
666 return kFALSE;
667
668 if (fPulserColor == MCalibrationCam::kNONE)
669 return kTRUE;
670
671 if (fNumProcessed == 0)
672 return kTRUE;
673
674 *fLog << endl;
675
676 return Finalize();
677}
678
679// -----------------------------------------------------------------------
680//
681// Return kTRUE if fPulserColor is kNONE
682//
683// First loop over pixels, average areas and sectors, call:
684// - FinalizePedestals()
685// - FinalizeCharges()
686// for every entry. Count number of valid pixels in loop and return kFALSE
687// if there are none (the "Michele check").
688//
689// Call FinalizeBadPixels()
690//
691// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
692//
693// Call FinalizeBlindCam()
694// Call FinalizePINDiode()
695//
696// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
697// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
698// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
699//
700// Call FinalizeUnsuitablePixels()
701//
702// Call MParContainer::SetReadyToSave() for fIntensCam, fCam, fQECam, fBadPixels and
703// fBlindCam and fPINDiode if they exist
704//
705// Print out some statistics
706//
707Int_t MCalibrationChargeCalc::Finalize()
708{
709
710 fNumHiGainSamples = fSignal->GetNumUsedHiGainFADCSlices();
711 fNumLoGainSamples = fSignal->GetNumUsedLoGainFADCSlices();
712
713 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
714 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
715
716 if (fPINDiode)
717 if (!fPINDiode->IsValid())
718 {
719 *fLog << warn << GetDescriptor()
720 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
721 fPINDiode = NULL;
722 }
723
724 MCalibrationBlindCam *blindcam = fIntensBlind
725 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
726 MCalibrationQECam *qecam = fIntensQE
727 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
728 MCalibrationChargeCam *chargecam = fIntensCam
729 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
730 MBadPixelsCam *badcam = fIntensBad
731 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
732
733 //
734 // First loop over pixels, call FinalizePedestals and FinalizeCharges
735 //
736 Int_t nvalid = 0;
737
738 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
739 {
740
741 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[pixid];
742 //
743 // Check if the pixel has been excluded from the fits
744 //
745 if (pix.IsExcluded())
746 continue;
747
748 MPedestalPix &ped = (*fPedestals)[pixid];
749 MBadPixelsPix &bad = (*badcam) [pixid];
750
751 const Int_t aidx = (*fGeom)[pixid].GetAidx();
752
753 FinalizePedestals(ped,pix,aidx);
754
755 if (FinalizeCharges(pix,bad,"pixel "))
756 nvalid++;
757 }
758
759 *fLog << endl;
760 //
761 // The Michele check ...
762 //
763 if (nvalid == 0)
764 {
765 if (!fIntensCam)
766 {
767 *fLog << warn << GetDescriptor() << ": All pixels have non-valid calibration. "
768 << "Did you forget to fill the histograms "
769 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
770 *fLog << warn << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
771 << "instead of a calibration run " << endl;
772 return kFALSE;
773 }
774 }
775
776 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
777 {
778
779 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
780 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
781
782 FinalizePedestals(ped,pix,aidx);
783 FinalizeCharges(pix,
784 fIntensCam ? fIntensCam->GetAverageBadArea(aidx) : fCam->GetAverageBadArea(aidx),
785 "area id");
786 }
787
788 *fLog << endl;
789
790 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
791 {
792
793 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
794
795 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
796 FinalizePedestals(ped,pix, 0);
797 }
798
799 *fLog << endl;
800
801 //
802 // Finalize Bad Pixels
803 //
804 FinalizeBadPixels();
805
806 //
807 // Finalize F-Factor method
808 //
809 if (FinalizeFFactorMethod())
810 chargecam->SetFFactorMethodValid(kTRUE);
811 else
812 {
813 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
814 chargecam->SetFFactorMethodValid(kFALSE);
815 return kFALSE;
816 }
817
818 *fLog << endl;
819 //
820 // Finalize Blind Pixel
821 //
822 qecam->SetBlindPixelMethodValid(FinalizeBlindCam());
823
824 //
825 // Finalize PIN Diode
826 //
827 qecam->SetBlindPixelMethodValid(FinalizePINDiode());
828
829 //
830 // Finalize QE Cam
831 //
832 FinalizeFFactorQECam();
833 FinalizeBlindPixelQECam();
834 FinalizePINDiodeQECam();
835 FinalizeCombinedQECam();
836
837 //
838 // Re-direct the output to an ascii-file from now on:
839 //
840 MLog *oldlog = fLog;
841 MLog asciilog;
842 if (!fOutputFile.IsNull())
843 {
844 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
845 SetLogStream(&asciilog);
846 }
847
848 //
849 // Finalize calibration statistics
850 //
851 FinalizeUnsuitablePixels();
852
853 chargecam->SetReadyToSave();
854 qecam ->SetReadyToSave();
855 badcam ->SetReadyToSave();
856
857 if (blindcam)
858 blindcam->SetReadyToSave();
859 if (fPINDiode)
860 fPINDiode->SetReadyToSave();
861
862 *fLog << inf << endl;
863 *fLog << GetDescriptor() << ": Fatal errors statistics:" << endl;
864
865 PrintUncalibrated(MBadPixelsPix::kChargeIsPedestal,
866 Form("%s%2.1f%s","Signal less than ",fChargeLimit," Pedestal RMS: "));
867 PrintUncalibrated(MBadPixelsPix::kChargeRelErrNotValid,
868 Form("%s%2.1f%s","Signal Error bigger than ",fChargeRelErrLimit," times Mean Signal: "));
869 PrintUncalibrated(MBadPixelsPix::kLoGainSaturation,
870 "Low Gain Saturation: ");
871 PrintUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin,
872 Form("%s%2.1f%s","Mean Abs. Arr. Time in First ",1.," Bin(s): "));
873 PrintUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins,
874 Form("%s%2.1f%s","Mean Abs. Arr. Time in Last ",2.," Bin(s): "));
875 PrintUncalibrated(MBadPixelsPix::kHiGainOverFlow,
876 "Pixels with High Gain Overflow: ");
877 PrintUncalibrated(MBadPixelsPix::kLoGainOverFlow,
878 "Pixels with Low Gain Overflow : ");
879
880 *fLog << inf << endl;
881 *fLog << GetDescriptor() << ": Unreliable errors statistics:" << endl;
882
883 PrintUncalibrated(MBadPixelsPix::kChargeSigmaNotValid,
884 "Signal Sigma smaller than Pedestal RMS: ");
885 PrintUncalibrated(MBadPixelsPix::kHiGainOscillating,
886 "Changing Hi Gain signal over time: ");
887 PrintUncalibrated(MBadPixelsPix::kLoGainOscillating,
888 "Changing Lo Gain signal over time: ");
889 PrintUncalibrated(MBadPixelsPix::kHiGainNotFitted,
890 "Unsuccesful Gauss fit to the Hi Gain: ");
891 PrintUncalibrated(MBadPixelsPix::kLoGainNotFitted,
892 "Unsuccesful Gauss fit to the Lo Gain: ");
893 PrintUncalibrated(MBadPixelsPix::kDeviatingNumPhes,
894 "Deviating number of phes: ");
895 PrintUncalibrated(MBadPixelsPix::kDeviatingFFactor,
896 "Deviating F-Factor: ");
897
898 if (!fOutputFile.IsNull())
899 SetLogStream(oldlog);
900
901 return kTRUE;
902}
903
904// ----------------------------------------------------------------------------------
905//
906// Retrieves pedestal and pedestal RMS from MPedestalPix
907// Retrieves total entries from MPedestalCam
908// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
909// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
910//
911// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
912// - MCalibrationChargePix::CalcLoGainPedestal()
913//
914void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
915{
916
917 //
918 // get the pedestals
919 //
920 const Float_t pedes = ped.GetPedestal();
921 const Float_t prms = ped.GetPedestalRms();
922 const Int_t num = fPedestals->GetTotalEntries();
923
924 //
925 // RMS error set by PedCalcFromLoGain, 0 in case MPedCalcPedRun was used.
926 //
927 const Float_t prmserr = num>0 ? prms/TMath::Sqrt(2.*num) : ped.GetPedestalRmsError();
928
929 //
930 // set them in the calibration camera
931 //
932 if (cal.IsHiGainSaturation())
933 {
934 cal.SetPedestal(pedes * fNumLoGainSamples,
935 prms * fSqrtLoGainSamples,
936 prmserr * fSqrtLoGainSamples);
937 cal.CalcLoGainPedestal(fNumLoGainSamples);
938 }
939 else
940 {
941
942 cal.SetPedestal(pedes * fNumHiGainSamples,
943 prms * fSqrtHiGainSamples,
944 prmserr * fSqrtHiGainSamples);
945 }
946
947}
948
949// ----------------------------------------------------------------------------------------------------
950//
951// Check fit results validity. Bad Pixels flags are set if:
952//
953// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
954// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
955// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
956// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
957// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
958//
959// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
960//
961// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
962// and returns kFALSE if not succesful.
963//
964// Calls MCalibrationChargePix::CalcFFactor() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
965// and returns kFALSE if not succesful.
966//
967// Calls MCalibrationChargePix::CalcConvFFactor()and sets flag: MBadPixelsPix::kDeviatingNumPhes)
968// and returns kFALSE if not succesful.
969//
970Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad, const char* what)
971{
972
973 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
974 return kFALSE;
975
976 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
977 {
978 *fLog << warn
979 << Form("Fitted Charge: %5.2f < %2.1f",cal.GetMean(),fChargeLimit)
980 << Form(" * Pedestal RMS %5.2f in %s%3i",cal.GetPedRms(),what,cal.GetPixId()) << endl;
981 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
982 }
983
984 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
985 {
986 *fLog << warn
987 << Form("Fitted Charge: %4.2f < %2.1f",cal.GetMean(),fChargeRelErrLimit)
988 << Form(" * its error %4.2f in %s%3i",cal.GetMeanErr(),what,cal.GetPixId()) << endl;
989 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
990 }
991
992 if (cal.GetSigma() < cal.GetPedRms())
993 {
994 *fLog << warn
995 << Form("Sigma of Fitted Charge: %6.2f <",cal.GetSigma())
996 << Form(" Ped. RMS=%5.2f in %s%3i",cal.GetPedRms(),what,cal.GetPixId()) << endl;
997 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
998 return kFALSE;
999 }
1000
1001 if (!cal.CalcReducedSigma())
1002 {
1003 *fLog << warn
1004 << Form("Could not calculate the reduced sigma in %s: ",what)
1005 << Form(" %4i",cal.GetPixId())
1006 << endl;
1007 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
1008 return kFALSE;
1009 }
1010
1011 if (!cal.CalcFFactor())
1012 {
1013 *fLog << warn
1014 << Form("Could not calculate the F-Factor in %s: ",what)
1015 << Form(" %4i",cal.GetPixId())
1016 << endl;
1017 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1018 return kFALSE;
1019 }
1020
1021 if (cal.GetPheFFactorMethod() < 0.)
1022 {
1023 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1024 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1025 cal.SetFFactorMethodValid(kFALSE);
1026 return kFALSE;
1027 }
1028
1029 if (!cal.CalcConvFFactor())
1030 {
1031 *fLog << warn
1032 << Form("Could not calculate the Conv. FADC counts to Phes in %s: ",what)
1033 << Form(" %4i",cal.GetPixId())
1034 << endl;
1035 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
1036 return kFALSE;
1037 }
1038
1039 return kTRUE;
1040}
1041
1042
1043
1044// -----------------------------------------------------------------------------------
1045//
1046// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
1047// - MBadPixelsPix::kChargeIsPedestal
1048// - MBadPixelsPix::kChargeRelErrNotValid
1049// - MBadPixelsPix::kMeanTimeInFirstBin
1050// - MBadPixelsPix::kMeanTimeInLast2Bins
1051// - MBadPixelsPix::kDeviatingNumPhes
1052// - MBadPixelsPix::kHiGainOverFlow
1053// - MBadPixelsPix::kLoGainOverFlow
1054//
1055// - Call MCalibrationPix::SetExcluded() for the bad pixels
1056//
1057// Sets pixel to MBadPixelsPix::kUnreliableRun, if one of the following flags is set:
1058// - MBadPixelsPix::kChargeSigmaNotValid
1059//
1060void MCalibrationChargeCalc::FinalizeBadPixels()
1061{
1062
1063 MBadPixelsCam *badcam = fIntensBad
1064 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1065 MCalibrationChargeCam *chargecam = fIntensCam
1066 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1067
1068 for (Int_t i=0; i<badcam->GetSize(); i++)
1069 {
1070
1071 MBadPixelsPix &bad = (*badcam) [i];
1072 MCalibrationPix &pix = (*chargecam)[i];
1073
1074 if (IsCheckDeadPixels())
1075 {
1076 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
1077 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1078
1079 if (bad.IsUncalibrated( MBadPixelsPix::kChargeErrNotValid ))
1080 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1081
1082 if (bad.IsUncalibrated( MBadPixelsPix::kChargeRelErrNotValid ))
1083 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1084 }
1085
1086 if (IsCheckExtractionWindow())
1087 {
1088 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
1089 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1090
1091 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
1092 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1093 }
1094
1095 if (IsCheckDeviatingBehavior())
1096 {
1097 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
1098 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
1099 }
1100
1101 if (IsCheckHistOverflow())
1102 {
1103 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOverFlow ))
1104 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1105
1106 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOverFlow ))
1107 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
1108 }
1109
1110 if (bad.IsUnsuitable( MBadPixelsPix::kUnsuitableRun ))
1111 pix.SetExcluded();
1112
1113 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
1114 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
1115 }
1116}
1117
1118// ------------------------------------------------------------------------
1119//
1120//
1121// First loop: Calculate a mean and mean RMS of photo-electrons per area index
1122// Include only pixels which are not MBadPixelsPix::kUnsuitableRun nor
1123// MBadPixelsPix::kChargeSigmaNotValid (see FinalizeBadPixels()) and set
1124// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
1125//
1126// Second loop: Get mean number of photo-electrons and its RMS including
1127// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
1128// and further exclude those deviating by more than fPheErrLimit mean
1129// sigmas from the mean (obtained in first loop). Set
1130// MBadPixelsPix::kDeviatingNumPhes if excluded.
1131//
1132// For the suitable pixels with flag MBadPixelsPix::kChargeSigmaNotValid
1133// set the number of photo-electrons as the mean number of photo-electrons
1134// calculated in that area index.
1135//
1136// Set weighted mean and variance of photo-electrons per area index in:
1137// average area pixels of MCalibrationChargeCam (obtained from:
1138// MCalibrationChargeCam::GetAverageArea() )
1139//
1140// Set weighted mean and variance of photo-electrons per sector in:
1141// average sector pixels of MCalibrationChargeCam (obtained from:
1142// MCalibrationChargeCam::GetAverageSector() )
1143//
1144//
1145// Third loop: Set mean number of photo-electrons and its RMS in the pixels
1146// only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1147//
1148Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
1149{
1150
1151 MBadPixelsCam *badcam = fIntensBad
1152 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1153 MCalibrationChargeCam *chargecam = fIntensCam
1154 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1155
1156 const Int_t npixels = fGeom->GetNumPixels();
1157 const Int_t nareas = fGeom->GetNumAreas();
1158 const Int_t nsectors = fGeom->GetNumSectors();
1159
1160 TArrayF lowlim (nareas);
1161 TArrayF upplim (nareas);
1162 TArrayD areavars (nareas);
1163 TArrayD areaweights (nareas);
1164 TArrayD sectorweights (nsectors);
1165 TArrayD areaphes (nareas);
1166 TArrayD sectorphes (nsectors);
1167 TArrayI numareavalid (nareas);
1168 TArrayI numsectorvalid(nsectors);
1169
1170 //
1171 // First loop: Get mean number of photo-electrons and the RMS
1172 // The loop is only to recognize later pixels with very deviating numbers
1173 //
1174 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
1175
1176 for (Int_t i=0; i<npixels; i++)
1177 {
1178
1179 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1180 MBadPixelsPix &bad = (*badcam)[i];
1181
1182 if (!pix.IsFFactorMethodValid())
1183 continue;
1184
1185 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1186 {
1187 pix.SetFFactorMethodValid(kFALSE);
1188 continue;
1189 }
1190
1191 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1192 continue;
1193
1194 const Float_t nphe = pix.GetPheFFactorMethod();
1195 const Int_t aidx = (*fGeom)[i].GetAidx();
1196
1197 camphes.Fill(i,nphe);
1198 camphes.SetUsed(i);
1199
1200 areaphes [aidx] += nphe;
1201 areavars [aidx] += nphe*nphe;
1202 numareavalid[aidx] ++;
1203 }
1204
1205 for (Int_t i=0; i<nareas; i++)
1206 {
1207 if (numareavalid[i] == 0)
1208 {
1209 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
1210 << "in area index: " << i << endl;
1211 continue;
1212 }
1213
1214 if (numareavalid[i] == 1)
1215 areavars[i] = 0.;
1216 else if (numareavalid[i] == 0)
1217 {
1218 areaphes[i] = -1.;
1219 areaweights[i] = -1.;
1220 }
1221 else
1222 {
1223 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1);
1224 areaphes[i] = areaphes[i] / numareavalid[i];
1225 }
1226
1227 if (areavars[i] < 0.)
1228 {
1229 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of photo-electrons found "
1230 << "in area index: " << i << endl;
1231 continue;
1232 }
1233
1234 lowlim [i] = areaphes[i] - fPheErrLimit*TMath::Sqrt(areavars[i]);
1235 upplim [i] = areaphes[i] + fPheErrLimit*TMath::Sqrt(areavars[i]);
1236
1237 TH1D *hist = camphes.ProjectionS(TArrayI(),TArrayI(1,&i),"_py",100);
1238 hist->Fit("gaus","Q");
1239 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1240 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1241 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1242
1243 if (IsDebug())
1244 hist->DrawClone();
1245
1246 if (ndf < 2)
1247 {
1248 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1249 << "in the camera with area index: " << i << endl;
1250 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1251 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1252 delete hist;
1253 continue;
1254 }
1255
1256 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1257
1258 if (prob < 0.001)
1259 {
1260 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1261 << "in the camera with area index: " << i << endl;
1262 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1263 << " is smaller than 0.001 " << endl;
1264 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1265 delete hist;
1266 continue;
1267 }
1268
1269 if (mean < 0.)
1270 {
1271 *fLog << inf << GetDescriptor() << ": Fitted mean number of photo-electrons "
1272 << "with area idx " << i << ": " << mean << " is smaller than 0. " << endl;
1273 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1274 delete hist;
1275 continue;
1276 }
1277
1278 *fLog << inf << GetDescriptor() << ": Mean number of phes with area idx " << i << ": "
1279 << Form("%7.2f+-%6.2f",mean,sigma) << endl;
1280
1281 lowlim [i] = mean - fPheErrLimit*sigma;
1282 upplim [i] = mean + fPheErrLimit*sigma;
1283
1284 delete hist;
1285 }
1286
1287 *fLog << endl;
1288
1289 numareavalid.Reset();
1290 areaphes .Reset();
1291 areavars .Reset();
1292 //
1293 // Second loop: Get mean number of photo-electrons and its RMS excluding
1294 // pixels deviating by more than fPheErrLimit sigma.
1295 // Set the conversion factor FADC counts to photo-electrons
1296 //
1297 for (Int_t i=0; i<npixels; i++)
1298 {
1299
1300 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1301
1302 if (!pix.IsFFactorMethodValid())
1303 continue;
1304
1305 MBadPixelsPix &bad = (*badcam)[i];
1306
1307 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1308 continue;
1309
1310 const Float_t nvar = pix.GetPheFFactorMethodVar();
1311 if (nvar <= 0.)
1312 {
1313 pix.SetFFactorMethodValid(kFALSE);
1314 continue;
1315 }
1316
1317 const Int_t aidx = (*fGeom)[i].GetAidx();
1318 const Int_t sector = (*fGeom)[i].GetSector();
1319 const Float_t area = (*fGeom)[i].GetA();
1320 const Float_t nphe = pix.GetPheFFactorMethod();
1321
1322 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1323 {
1324 *fLog << warn << "Number of phes: "
1325 << Form("%7.2f out of %3.1f sigma limit: ",nphe,fPheErrLimit)
1326 << Form("[%7.2f,%7.2f] pixel%4i",lowlim[aidx],upplim[aidx],i) << endl;
1327 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1328 if (IsCheckDeviatingBehavior())
1329 {
1330 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1331 // pix.SetFFactorMethodValid(kFALSE);
1332 }
1333 continue;
1334 }
1335
1336 areaweights [aidx] += nphe*nphe;
1337 areaphes [aidx] += nphe;
1338 numareavalid [aidx] ++;
1339
1340 if (aidx == 0)
1341 fNumInnerFFactorMethodUsed++;
1342
1343 sectorweights [sector] += nphe*nphe/area/area;
1344 sectorphes [sector] += nphe/area;
1345 numsectorvalid[sector] ++;
1346 }
1347
1348 *fLog << endl;
1349
1350 for (Int_t aidx=0; aidx<nareas; aidx++)
1351 {
1352
1353 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1354
1355 if (numareavalid[aidx] == 1)
1356 areaweights[aidx] = 0.;
1357 else if (numareavalid[aidx] == 0)
1358 {
1359 areaphes[aidx] = -1.;
1360 areaweights[aidx] = -1.;
1361 }
1362 else
1363 {
1364 areaweights[aidx] = (areaweights[aidx] - areaphes[aidx]*areaphes[aidx]/numareavalid[aidx])
1365 / (numareavalid[aidx]-1);
1366 areaphes[aidx] /= numareavalid[aidx];
1367 }
1368
1369 if (areaweights[aidx] < 0. || areaphes[aidx] <= 0.)
1370 {
1371 *fLog << warn << GetDescriptor()
1372 << ": Mean number phes from area index " << aidx << " could not be calculated: "
1373 << " Mean: " << areaphes[aidx]
1374 << " Variance: " << areaweights[aidx] << endl;
1375 apix.SetFFactorMethodValid(kFALSE);
1376 continue;
1377 }
1378
1379 *fLog << inf << GetDescriptor()
1380 << ": Average total phes for area idx " << aidx << ": "
1381 << Form("%7.2f +- 6.2f",areaphes[aidx],TMath::Sqrt(areaweights[aidx])) << endl;
1382
1383 apix.SetPheFFactorMethod ( areaphes[aidx] );
1384 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1385 apix.SetFFactorMethodValid ( kTRUE );
1386
1387 }
1388
1389 *fLog << endl;
1390
1391 for (Int_t sector=0; sector<nsectors; sector++)
1392 {
1393
1394 if (numsectorvalid[sector] == 1)
1395 sectorweights[sector] = 0.;
1396 else if (numsectorvalid[sector] == 0)
1397 {
1398 sectorphes[sector] = -1.;
1399 sectorweights[sector] = -1.;
1400 }
1401 else
1402 {
1403 sectorweights[sector] = (sectorweights[sector]
1404 - sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector]
1405 )
1406 / (numsectorvalid[sector]-1.);
1407 sectorphes[sector] /= numsectorvalid[sector];
1408 }
1409
1410 MCalibrationChargePix &spix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
1411
1412 if (sectorweights[sector] < 0. || sectorphes[sector] <= 0.)
1413 {
1414 *fLog << warn << GetDescriptor()
1415 <<": Mean number phes/area for sector " << sector << " could not be calculated: "
1416 << " Mean: " << sectorphes[sector]
1417 << " Variance: " << sectorweights[sector] << endl;
1418 spix.SetFFactorMethodValid(kFALSE);
1419 continue;
1420 }
1421
1422 *fLog << inf << GetDescriptor()
1423 << ": Avg number phes/mm^2 in sector " << sector << ": "
1424 << Form("%5.3f+-%4.3f",sectorphes[sector],TMath::Sqrt(sectorweights[sector]))
1425 << endl;
1426
1427 spix.SetPheFFactorMethod ( sectorphes[sector] );
1428 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1429 spix.SetFFactorMethodValid ( kTRUE );
1430
1431 }
1432
1433 //
1434 // Third loop: Set mean number of photo-electrons and its RMS in the pixels
1435 // only excluded as: MBadPixelsPix::kChargeSigmaNotValid
1436 //
1437 for (Int_t i=0; i<npixels; i++)
1438 {
1439
1440 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1441 MBadPixelsPix &bad = (*badcam)[i];
1442
1443 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1444 continue;
1445
1446 if (bad.IsUncalibrated(MBadPixelsPix::kChargeSigmaNotValid))
1447 {
1448 const Int_t aidx = (*fGeom)[i].GetAidx();
1449 MCalibrationChargePix &apix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
1450
1451 pix.SetPheFFactorMethod ( apix.GetPheFFactorMethod() );
1452 pix.SetPheFFactorMethodVar( apix.GetPheFFactorMethodVar() );
1453 if (!pix.CalcConvFFactor())
1454 {
1455 *fLog << warn << GetDescriptor()
1456 << ": Could not calculate the Conv. FADC counts to Phes in pixel: "
1457 << Form(" %4i",pix.GetPixId())
1458 << endl;
1459 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1460 if (IsCheckDeviatingBehavior())
1461 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1462 }
1463 }
1464 }
1465
1466 return kTRUE;
1467}
1468
1469
1470
1471// ------------------------------------------------------------------------
1472//
1473// Returns kFALSE if pointer to MCalibrationBlindCam is NULL
1474//
1475// The check returns kFALSE if:
1476//
1477// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1478// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1479//
1480// Calls:
1481// - MCalibrationBlindPix::CalcFluxInsidePlexiglass()
1482//
1483Bool_t MCalibrationChargeCalc::FinalizeBlindCam()
1484{
1485
1486 MCalibrationBlindCam *blindcam = fIntensBlind
1487 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
1488
1489 if (!blindcam)
1490 return kFALSE;
1491
1492 Int_t nvalid = 0;
1493
1494 for (Int_t i=0; i<blindcam->GetSize(); i++)
1495 {
1496
1497 MCalibrationBlindPix &blindpix = (MCalibrationBlindPix&)(*blindcam)[i];
1498
1499 if (!blindpix.IsValid())
1500 continue;
1501
1502 const Float_t lambda = blindpix.GetLambda();
1503 const Float_t lambdaerr = blindpix.GetLambdaErr();
1504 const Float_t lambdacheck = blindpix.GetLambdaCheck();
1505
1506 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1507 {
1508 *fLog << warn << GetDescriptor()
1509 << Form("%s%4.2f%s%4.2f%s%4.2f%s%2i",": Lambda: ",lambda," and Lambda-Check: ",
1510 lambdacheck," differ by more than ",fLambdaCheckLimit," in the Blind Pixel Nr.",i)
1511 << endl;
1512 blindpix.SetValid(kFALSE);
1513 continue;
1514 }
1515
1516 if (lambdaerr > fLambdaErrLimit)
1517 {
1518 *fLog << warn << GetDescriptor()
1519 << Form("%s%4.2f%s%4.2f%s%2i",": Error of Fitted Lambda: ",lambdaerr," is greater than ",
1520 fLambdaErrLimit," in Blind Pixel Nr.",i) << endl;
1521 blindpix.SetValid(kFALSE);
1522 continue;
1523 }
1524
1525 if (!blindpix.CalcFluxInsidePlexiglass())
1526 {
1527 *fLog << warn << "Could not calculate the flux of photons from Blind Pixel Nr." << i << endl;
1528 blindpix.SetValid(kFALSE);
1529 continue;
1530 }
1531
1532 nvalid++;
1533 }
1534
1535 if (!nvalid)
1536 return kFALSE;
1537
1538 return kTRUE;
1539}
1540
1541// ------------------------------------------------------------------------
1542//
1543// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1544//
1545// The check returns kFALSE if:
1546//
1547// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1548// 2) PINDiode has a fit error smaller than fChargeErrLimit
1549// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1550// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1551//
1552// Calls:
1553// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1554//
1555Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1556{
1557
1558 if (!fPINDiode)
1559 return kFALSE;
1560
1561 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1562 {
1563 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1564 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1565 return kFALSE;
1566 }
1567
1568 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1569 {
1570 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1571 << fChargeErrLimit << " in PINDiode " << endl;
1572 return kFALSE;
1573 }
1574
1575 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1576 {
1577 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1578 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1579 return kFALSE;
1580 }
1581
1582 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1583 {
1584 *fLog << warn << GetDescriptor()
1585 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1586 return kFALSE;
1587 }
1588
1589
1590 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1591 {
1592 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1593 << "will skip PIN Diode Calibration " << endl;
1594 return kFALSE;
1595 }
1596
1597 return kTRUE;
1598}
1599
1600// ------------------------------------------------------------------------
1601//
1602// Calculate the average number of photons outside the plexiglass with the
1603// formula:
1604//
1605// av.Num.photons(area index) = av.Num.Phes(area index)
1606// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1607// / MCalibrationQEPix::GetPMTCollectionEff()
1608// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1609// / MCalibrationQECam::GetPlexiglassQE()
1610//
1611// Calculate the variance on the average number of photons assuming that the error on the
1612// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1613// values keeps it ordinary error since it is systematic.
1614//
1615// Loop over pixels:
1616//
1617// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1618// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1619//
1620// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1621// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1622//
1623// - Calculate the quantum efficiency with the formula:
1624//
1625// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1626//
1627// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1628//
1629// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1630// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1631//
1632// - Call MCalibrationQEPix::UpdateFFactorMethod()
1633//
1634void MCalibrationChargeCalc::FinalizeFFactorQECam()
1635{
1636
1637 if (fNumInnerFFactorMethodUsed < 2)
1638 {
1639 *fLog << warn << GetDescriptor()
1640 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1641 return;
1642 }
1643
1644 MCalibrationQECam *qecam = fIntensQE
1645 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1646 MCalibrationChargeCam *chargecam = fIntensCam
1647 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1648 MBadPixelsCam *badcam = fIntensBad
1649 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1650
1651 MCalibrationChargePix &avpix = (MCalibrationChargePix&)chargecam->GetAverageArea(0);
1652 MCalibrationQEPix &qepix = (MCalibrationQEPix&) qecam->GetAverageArea(0);
1653
1654 const Float_t avphotons = avpix.GetPheFFactorMethod()
1655 / qepix.GetDefaultQE(fPulserColor)
1656 / qepix.GetPMTCollectionEff()
1657 / qepix.GetLightGuidesEff(fPulserColor)
1658 / qecam->GetPlexiglassQE();
1659
1660 const Float_t avphotrelvar = avpix.GetPheFFactorMethodRelVar()
1661 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1662 + qepix.GetPMTCollectionEffRelVar()
1663 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1664 + qecam->GetPlexiglassQERelVar();
1665
1666 const UInt_t nareas = fGeom->GetNumAreas();
1667
1668 //
1669 // Set the results in the MCalibrationChargeCam
1670 //
1671 chargecam->SetNumPhotonsFFactorMethod (avphotons);
1672
1673 if (avphotrelvar > 0.)
1674 chargecam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1675
1676 TArrayF lowlim (nareas);
1677 TArrayF upplim (nareas);
1678 TArrayD avffactorphotons (nareas);
1679 TArrayD avffactorphotvar (nareas);
1680 TArrayI numffactor (nareas);
1681
1682 const UInt_t npixels = fGeom->GetNumPixels();
1683
1684 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1685
1686 for (UInt_t i=0; i<npixels; i++)
1687 {
1688
1689 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1690 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1691 MBadPixelsPix &bad = (*badcam) [i];
1692
1693 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1694 continue;
1695
1696 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1697 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1698
1699 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1700
1701 qepix.SetQEFFactor ( qe , fPulserColor );
1702 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1703 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1704
1705 if (!qepix.UpdateFFactorMethod( qecam->GetPlexiglassQE() ))
1706 *fLog << warn << GetDescriptor()
1707 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1708
1709 //
1710 // The following pixels are those with deviating sigma, but otherwise OK,
1711 // probably those with stars during the pedestal run, but not the cal. run.
1712 //
1713 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1714 {
1715 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1716 if (IsCheckDeviatingBehavior())
1717 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1718 continue;
1719 }
1720
1721 const Int_t aidx = (*fGeom)[i].GetAidx();
1722 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1723
1724 camffactor.Fill(i,ffactor);
1725 camffactor.SetUsed(i);
1726
1727 avffactorphotons[aidx] += ffactor;
1728 avffactorphotvar[aidx] += ffactor*ffactor;
1729 numffactor[aidx]++;
1730 }
1731
1732 for (UInt_t i=0; i<nareas; i++)
1733 {
1734
1735 if (numffactor[i] == 0)
1736 {
1737 *fLog << warn << GetDescriptor() << ": No pixels with valid total F-Factor found "
1738 << "in area index: " << i << endl;
1739 continue;
1740 }
1741
1742 avffactorphotvar[i] = (avffactorphotvar[i] - avffactorphotons[i]*avffactorphotons[i]/numffactor[i])
1743 / (numffactor[i]-1.);
1744 avffactorphotons[i] = avffactorphotons[i] / numffactor[i];
1745
1746 if (avffactorphotvar[i] < 0.)
1747 {
1748 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of total F-Factor found "
1749 << "in area index: " << i << endl;
1750 continue;
1751 }
1752
1753 lowlim [i] = 1.; // Lowest known F-Factor of a PMT
1754 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1755
1756 TArrayI area(1);
1757 area[0] = i;
1758
1759 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1760 hist->Fit("gaus","Q");
1761 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1762 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1763 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1764
1765 if (IsDebug())
1766 camffactor.DrawClone();
1767
1768 if (ndf < 2)
1769 {
1770 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1771 << "in the camera with area index: " << i << endl;
1772 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1773 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1774 delete hist;
1775 continue;
1776 }
1777
1778 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1779
1780 if (prob < 0.001)
1781 {
1782 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1783 << "in the camera with area index: " << i << endl;
1784 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1785 << " is smaller than 0.001 " << endl;
1786 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1787 delete hist;
1788 continue;
1789 }
1790
1791 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1792 << "with area index #" << i << ": "
1793 << Form("%4.2f+-%4.2f",mean,sigma) << endl;
1794
1795 lowlim [i] = 1.;
1796 upplim [i] = mean + fFFactorErrLimit*sigma;
1797
1798 delete hist;
1799 }
1800
1801 *fLog << endl;
1802
1803 for (UInt_t i=0; i<npixels; i++)
1804 {
1805
1806 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1807 MBadPixelsPix &bad = (*badcam) [i];
1808
1809 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1810 continue;
1811
1812 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1813 const Int_t aidx = (*fGeom)[i].GetAidx();
1814
1815 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1816 {
1817 *fLog << warn << "Overall F-Factor "
1818 << Form("%5.2f",ffactor) << " out of range ["
1819 << Form("%5.2f,%5.2f",lowlim[aidx],upplim[aidx]) << "] Pixel " << i << endl;
1820
1821 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1822 if (IsCheckDeviatingBehavior())
1823 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1824 }
1825 }
1826
1827 for (UInt_t i=0; i<npixels; i++)
1828 {
1829
1830 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1831 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1832 MBadPixelsPix &bad = (*badcam) [i];
1833
1834 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1835 {
1836 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1837 pix.SetFFactorMethodValid(kFALSE);
1838 pix.SetExcluded();
1839 continue;
1840 }
1841 }
1842}
1843
1844
1845// ------------------------------------------------------------------------
1846//
1847// Loop over pixels:
1848//
1849// - Continue, if not MCalibrationBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1850// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1851//
1852// - Calculate the quantum efficiency with the formula:
1853//
1854// QE = Num.Phes / MCalibrationBlindPix::GetFluxInsidePlexiglass()
1855// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1856//
1857// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1858// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1859// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1860//
1861// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1862//
1863void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1864{
1865
1866
1867 MCalibrationBlindCam *blindcam = fIntensBlind
1868 ? (MCalibrationBlindCam*) fIntensBlind->GetCam(): fBlindCam;
1869 MBadPixelsCam *badcam = fIntensBad
1870 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1871 MCalibrationQECam *qecam = fIntensQE
1872 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1873 MCalibrationChargeCam *chargecam = fIntensCam
1874 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1875
1876 //
1877 // Set the results in the MCalibrationChargeCam
1878 //
1879 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
1880 {
1881
1882 const Float_t photons = blindcam->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
1883 / qecam->GetPlexiglassQE();
1884 chargecam->SetNumPhotonsBlindPixelMethod(photons);
1885
1886 const Float_t photrelvar = blindcam->GetFluxInsidePlexiglassRelVar()
1887 + qecam->GetPlexiglassQERelVar();
1888
1889 if (photrelvar > 0.)
1890 chargecam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
1891 }
1892
1893 //
1894 // With the knowledge of the overall photon flux, calculate the
1895 // quantum efficiencies after the Blind Pixel and PIN Diode method
1896 //
1897 const UInt_t npixels = fGeom->GetNumPixels();
1898 for (UInt_t i=0; i<npixels; i++)
1899 {
1900
1901 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1902
1903 if (!blindcam || !(blindcam->IsFluxInsidePlexiglassAvailable()))
1904 {
1905 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1906 continue;
1907 }
1908
1909 MBadPixelsPix &bad = (*badcam) [i];
1910
1911 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1912 {
1913 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1914 continue;
1915 }
1916
1917 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1918 MGeomPix &geo = (*fGeom) [i];
1919
1920 const Float_t qe = pix.GetPheFFactorMethod()
1921 / blindcam->GetFluxInsidePlexiglass()
1922 / geo.GetA()
1923 * qecam->GetPlexiglassQE();
1924
1925 const Float_t qerelvar = blindcam->GetFluxInsidePlexiglassRelVar()
1926 + qecam->GetPlexiglassQERelVar()
1927 + pix.GetPheFFactorMethodRelVar();
1928
1929 qepix.SetQEBlindPixel ( qe , fPulserColor );
1930 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
1931 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
1932
1933 if (!qepix.UpdateBlindPixelMethod( qecam->GetPlexiglassQE()))
1934 *fLog << warn << GetDescriptor()
1935 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
1936 }
1937}
1938
1939// ------------------------------------------------------------------------
1940//
1941// Loop over pixels:
1942//
1943// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
1944// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
1945//
1946// - Calculate the quantum efficiency with the formula:
1947//
1948// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
1949//
1950// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
1951// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
1952// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
1953//
1954// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
1955//
1956void MCalibrationChargeCalc::FinalizePINDiodeQECam()
1957{
1958
1959 const UInt_t npixels = fGeom->GetNumPixels();
1960
1961 MCalibrationQECam *qecam = fIntensQE
1962 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1963 MCalibrationChargeCam *chargecam = fIntensCam
1964 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1965 MBadPixelsCam *badcam = fIntensBad
1966 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1967 //
1968 // With the knowledge of the overall photon flux, calculate the
1969 // quantum efficiencies after the PIN Diode method
1970 //
1971 for (UInt_t i=0; i<npixels; i++)
1972 {
1973
1974 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1975
1976 if (!fPINDiode)
1977 {
1978 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1979 continue;
1980 }
1981
1982 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
1983 {
1984 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1985 continue;
1986 }
1987
1988 MBadPixelsPix &bad = (*badcam) [i];
1989
1990 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1991 {
1992 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1993 continue;
1994 }
1995
1996 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1997 MGeomPix &geo = (*fGeom) [i];
1998
1999 const Float_t qe = pix.GetPheFFactorMethod()
2000 / fPINDiode->GetFluxOutsidePlexiglass()
2001 / geo.GetA();
2002
2003 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
2004
2005 qepix.SetQEPINDiode ( qe , fPulserColor );
2006 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
2007 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
2008
2009 if (!qepix.UpdatePINDiodeMethod())
2010 *fLog << warn << GetDescriptor()
2011 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
2012 }
2013}
2014
2015// ------------------------------------------------------------------------
2016//
2017// Loop over pixels:
2018//
2019// - Call MCalibrationQEPix::UpdateCombinedMethod()
2020//
2021void MCalibrationChargeCalc::FinalizeCombinedQECam()
2022{
2023
2024 const UInt_t npixels = fGeom->GetNumPixels();
2025
2026 MCalibrationQECam *qecam = fIntensQE
2027 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
2028 MBadPixelsCam *badcam = fIntensBad
2029 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
2030
2031 for (UInt_t i=0; i<npixels; i++)
2032 {
2033
2034 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2035 MBadPixelsPix &bad = (*badcam) [i];
2036
2037 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2038 {
2039 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2040 continue;
2041 }
2042
2043 qepix.UpdateCombinedMethod();
2044 }
2045}
2046
2047// -----------------------------------------------------------------------------------------------
2048//
2049// - Print out statistics about BadPixels of type UnsuitableType_t
2050// - store numbers of bad pixels of each type in fCam or fIntensCam
2051//
2052void MCalibrationChargeCalc::FinalizeUnsuitablePixels()
2053{
2054
2055 *fLog << inf << endl;
2056 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
2057 *fLog << dec << setfill(' ');
2058
2059 const Int_t nareas = fGeom->GetNumAreas();
2060
2061 TArrayI counts(nareas);
2062
2063 MBadPixelsCam *badcam = fIntensBad
2064 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2065 MCalibrationChargeCam *chargecam = fIntensCam
2066 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
2067
2068 for (Int_t i=0; i<badcam->GetSize(); i++)
2069 {
2070 MBadPixelsPix &bad = (*badcam)[i];
2071 if (!bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
2072 {
2073 const Int_t aidx = (*fGeom)[i].GetAidx();
2074 counts[aidx]++;
2075 }
2076 }
2077
2078 if (fGeom->InheritsFrom("MGeomCamMagic"))
2079 *fLog << " " << setw(7) << "Successfully calibrated Pixels: "
2080 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2081
2082 counts.Reset();
2083
2084 for (Int_t i=0; i<badcam->GetSize(); i++)
2085 {
2086 MBadPixelsPix &bad = (*badcam)[i];
2087
2088 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
2089 {
2090 const Int_t aidx = (*fGeom)[i].GetAidx();
2091 counts[aidx]++;
2092 }
2093 }
2094
2095 for (Int_t aidx=0; aidx<nareas; aidx++)
2096 chargecam->SetNumUnsuitable(counts[aidx], aidx);
2097
2098 if (fGeom->InheritsFrom("MGeomCamMagic"))
2099 *fLog << " " << setw(7) << "Uncalibrated Pixels: "
2100 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2101
2102 counts.Reset();
2103
2104 for (Int_t i=0; i<badcam->GetSize(); i++)
2105 {
2106
2107 MBadPixelsPix &bad = (*badcam)[i];
2108
2109 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
2110 {
2111 const Int_t aidx = (*fGeom)[i].GetAidx();
2112 counts[aidx]++;
2113 }
2114 }
2115
2116 for (Int_t aidx=0; aidx<nareas; aidx++)
2117 chargecam->SetNumUnreliable(counts[aidx], aidx);
2118
2119 *fLog << " " << setw(7) << "Unreliable Pixels: "
2120 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2121
2122}
2123
2124// -----------------------------------------------------------------------------------------------
2125//
2126// Print out statistics about BadPixels of type UncalibratedType_t
2127//
2128void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
2129{
2130
2131 UInt_t countinner = 0;
2132 UInt_t countouter = 0;
2133
2134 MBadPixelsCam *badcam = fIntensBad
2135 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2136
2137 for (Int_t i=0; i<badcam->GetSize(); i++)
2138 {
2139 MBadPixelsPix &bad = (*badcam)[i];
2140
2141 if (bad.IsUncalibrated(typ))
2142 {
2143 if (fGeom->GetPixRatio(i) == 1.)
2144 countinner++;
2145 else
2146 countouter++;
2147 }
2148 }
2149
2150 *fLog << " " << setw(7) << text
2151 << Form("%s%3i%s%3i","Inner: ",countinner," Outer: ",countouter) << endl;
2152}
2153
2154// --------------------------------------------------------------------------
2155//
2156// Set the path for output file
2157//
2158void MCalibrationChargeCalc::SetOutputPath(TString path)
2159{
2160 fOutputPath = path;
2161 if (fOutputPath.EndsWith("/"))
2162 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
2163}
2164
2165// --------------------------------------------------------------------------
2166//
2167// Set the output file
2168//
2169void MCalibrationChargeCalc::SetOutputFile(TString file)
2170{
2171 fOutputFile = file;
2172}
2173
2174// --------------------------------------------------------------------------
2175//
2176// Get the output file
2177//
2178const char* MCalibrationChargeCalc::GetOutputFile()
2179{
2180 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
2181}
2182
2183// --------------------------------------------------------------------------
2184//
2185// Read the environment for the following data members:
2186// - fChargeLimit
2187// - fChargeErrLimit
2188// - fChargeRelErrLimit
2189// - fDebug
2190// - fFFactorErrLimit
2191// - fLambdaErrLimit
2192// - fLambdaCheckErrLimit
2193// - fPheErrLimit
2194//
2195Int_t MCalibrationChargeCalc::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
2196{
2197
2198 Bool_t rc = kFALSE;
2199 if (IsEnvDefined(env, prefix, "ChargeLimit", print))
2200 {
2201 SetChargeLimit(GetEnvValue(env, prefix, "ChargeLimit", fChargeLimit));
2202 rc = kTRUE;
2203 }
2204 if (IsEnvDefined(env, prefix, "ChargeErrLimit", print))
2205 {
2206 SetChargeErrLimit(GetEnvValue(env, prefix, "ChargeErrLimit", fChargeErrLimit));
2207 rc = kTRUE;
2208 }
2209 if (IsEnvDefined(env, prefix, "ChargeRelErrLimit", print))
2210 {
2211 SetChargeRelErrLimit(GetEnvValue(env, prefix, "ChargeRelErrLimit", fChargeRelErrLimit));
2212 rc = kTRUE;
2213 }
2214 if (IsEnvDefined(env, prefix, "Debug", print))
2215 {
2216 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
2217 rc = kTRUE;
2218 }
2219 if (IsEnvDefined(env, prefix, "FFactorErrLimit", print))
2220 {
2221 SetFFactorErrLimit(GetEnvValue(env, prefix, "FFactorErrLimit", fFFactorErrLimit));
2222 rc = kTRUE;
2223 }
2224 if (IsEnvDefined(env, prefix, "LambdaErrLimit", print))
2225 {
2226 SetLambdaErrLimit(GetEnvValue(env, prefix, "LambdaErrLimit", fLambdaErrLimit));
2227 rc = kTRUE;
2228 }
2229 if (IsEnvDefined(env, prefix, "LambdaCheckLimit", print))
2230 {
2231 SetLambdaCheckLimit(GetEnvValue(env, prefix, "LambdaCheckLimit", fLambdaCheckLimit));
2232 rc = kTRUE;
2233 }
2234 if (IsEnvDefined(env, prefix, "PheErrLimit", print))
2235 {
2236 SetPheErrLimit(GetEnvValue(env, prefix, "PheErrLimit", fPheErrLimit));
2237 rc = kTRUE;
2238 }
2239 if (IsEnvDefined(env, prefix, "CheckDeadPixels", print))
2240 {
2241 SetCheckDeadPixels(GetEnvValue(env, prefix, "CheckDeadPixels", IsCheckDeadPixels()));
2242 rc = kTRUE;
2243 }
2244 if (IsEnvDefined(env, prefix, "CheckDeviatingBehavior", print))
2245 {
2246 SetCheckDeviatingBehavior(GetEnvValue(env, prefix, "CheckDeviatingBehavior", IsCheckDeviatingBehavior()));
2247 rc = kTRUE;
2248 }
2249 if (IsEnvDefined(env, prefix, "CheckExtractionWindow", print))
2250 {
2251 SetCheckExtractionWindow(GetEnvValue(env, prefix, "CheckExtractionWindow", IsCheckExtractionWindow()));
2252 rc = kTRUE;
2253 }
2254 if (IsEnvDefined(env, prefix, "CheckHistOverflow", print))
2255 {
2256 SetCheckHistOverflow(GetEnvValue(env, prefix, "CheckHistOverflow", IsCheckHistOverflow()));
2257 rc = kTRUE;
2258 }
2259 if (IsEnvDefined(env, prefix, "CheckOscillations", print))
2260 {
2261 SetCheckOscillations(GetEnvValue(env, prefix, "CheckOscillations", IsCheckOscillations()));
2262 rc = kTRUE;
2263 }
2264
2265 return rc;
2266}
2267
Note: See TracBrowser for help on using the repository browser.