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

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