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

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