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

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