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

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