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

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