source: tags/Mars-V0.9/mcalib/MCalibrationChargeCalc.cc

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