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

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