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

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