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

Last change on this file since 5863 was 5852, 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 "MRawEvtHeader.h"
197
198#include "MGeomCam.h"
199#include "MGeomPix.h"
200#include "MHCamera.h"
201
202#include "MPedestalCam.h"
203#include "MPedestalPix.h"
204
205#include "MCalibrationIntensityChargeCam.h"
206#include "MCalibrationIntensityQECam.h"
207#include "MCalibrationIntensityBlindCam.h"
208
209#include "MHCalibrationChargeCam.h"
210#include "MHCalibrationChargeBlindCam.h"
211
212#include "MCalibrationChargeCam.h"
213#include "MCalibrationChargePix.h"
214#include "MCalibrationChargePINDiode.h"
215#include "MCalibrationBlindPix.h"
216#include "MCalibrationBlindCam.h"
217
218#include "MExtractedSignalCam.h"
219#include "MExtractedSignalPix.h"
220#include "MExtractedSignalBlindPixel.h"
221#include "MExtractedSignalPINDiode.h"
222
223#include "MBadPixelsIntensityCam.h"
224#include "MBadPixelsCam.h"
225
226#include "MCalibrationQECam.h"
227#include "MCalibrationQEPix.h"
228
229ClassImp(MCalibrationChargeCalc);
230
231using namespace std;
232
233const Float_t MCalibrationChargeCalc::fgChargeLimit = 2.5;
234const Float_t MCalibrationChargeCalc::fgChargeErrLimit = 0.;
235const Float_t MCalibrationChargeCalc::fgChargeRelErrLimit = 1.;
236const Float_t MCalibrationChargeCalc::fgLambdaErrLimit = 0.2;
237const Float_t MCalibrationChargeCalc::fgLambdaCheckLimit = 0.5;
238const Float_t MCalibrationChargeCalc::fgPheErrLimit = 4.5;
239const Float_t MCalibrationChargeCalc::fgFFactorErrLimit = 4.5;
240const TString MCalibrationChargeCalc::fgNamePedestalCam = "MPedestalCam";
241// --------------------------------------------------------------------------
242//
243// 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), fEvtHeader(NULL)
274{
275
276 fName = name ? name : "MCalibrationChargeCalc";
277 fTitle = title ? title : "Task to calculate the calibration constants and MCalibrationCam ";
278
279 AddToBranchList("MRawEvtData.fHiGainPixId");
280 AddToBranchList("MRawEvtData.fLoGainPixId");
281 AddToBranchList("MRawEvtData.fHiGainFadcSamples");
282 AddToBranchList("MRawEvtData.fLoGainFadcSamples");
283
284 SetChargeLimit ();
285 SetChargeErrLimit ();
286 SetChargeRelErrLimit ();
287 SetFFactorErrLimit ();
288 SetLambdaCheckLimit ();
289 SetLambdaErrLimit ();
290 SetNamePedestalCam ();
291 SetPheErrLimit ();
292 SetOutputPath ();
293 SetOutputFile ();
294 SetDebug ( kFALSE );
295
296 SetCheckDeadPixels ();
297 SetCheckDeviatingBehavior();
298 SetCheckExtractionWindow ();
299 SetCheckHistOverflow ();
300 SetCheckOscillations ();
301
302 Clear();
303
304}
305
306// --------------------------------------------------------------------------
307//
308// Sets:
309// - all variables to 0.,
310// - all flags to kFALSE
311// - all pointers to NULL
312// - the pulser colour to kNONE
313// - fBlindPixelFlags to 0
314// - fPINDiodeFlags to 0
315//
316void MCalibrationChargeCalc::Clear(const Option_t *o)
317{
318
319 fNumHiGainSamples = 0.;
320 fNumLoGainSamples = 0.;
321 fSqrtHiGainSamples = 0.;
322 fSqrtLoGainSamples = 0.;
323 fNumInnerFFactorMethodUsed = 0;
324
325 fNumProcessed = 0;
326
327 fIntensBad = NULL;
328 fBadPixels = NULL;
329 fIntensCam = NULL;
330 fCam = NULL;
331 fHCam = NULL;
332 fIntensQE = NULL;
333 fQECam = NULL;
334 fIntensBlind = NULL;
335 fBlindCam = NULL;
336 fHBlindCam = NULL;
337 fPINDiode = NULL;
338 fPedestals = NULL;
339
340 SetPulserColor ( MCalibrationCam::kNONE );
341
342 fBlindPixelFlags.Set(0);
343 fPINDiodeFlags .Set(0);
344 fResultFlags .Set(0);
345}
346
347
348// -----------------------------------------------------------------------------------
349//
350// The following container are searched for and execution aborted if not in MParList:
351// - MPedestalCam
352// - MRawEvtHeader
353// - MExtractedSignalCam
354//
355Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
356{
357
358 /*
359 if (IsInterlaced())
360 {
361 fTrigPattern = (MTriggerPattern*)pList->FindObject("MTriggerPattern");
362 if (!fTrigPattern)
363 {
364 *fLog << err << "MTriggerPattern not found... abort." << endl;
365 return kFALSE;
366 }
367 }
368 */
369
370 fEvtHeader = (MRawEvtHeader*)pList->FindObject("MRawEvtHeader");
371 if (!fEvtHeader)
372 {
373 *fLog << err << "MRawEvtHeader not found... abort." << endl;
374 return kFALSE;
375 }
376
377 //
378 // Containers that have to be there.
379 //
380 fSignal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
381 if (!fSignal)
382 {
383 *fLog << err << "MExtractedSignalCam not found... aborting" << endl;
384 return kFALSE;
385 }
386
387 if (fPedestals)
388 return kTRUE;
389
390 fPedestals = (MPedestalCam*)pList->FindObject(AddSerialNumber(fNamePedestalCam), "MPedestalCam");
391 if (!fPedestals)
392 {
393 *fLog << err << GetDescriptor() << ": " << fNamePedestalCam << " not found... aborting" << endl;
394 return kFALSE;
395 }
396
397 return kTRUE;
398}
399
400
401// --------------------------------------------------------------------------
402//
403// Search for the following input containers and abort if not existing:
404// - MGeomCam
405// - MCalibrationIntensityChargeCam or MCalibrationChargeCam
406// - MCalibrationIntensityQECam or MCalibrationQECam
407// - MBadPixelsIntensityCam or MBadPixelsCam
408//
409// Search for the following input containers and give a warning if not existing:
410// - MCalibrationBlindPix
411// - MCalibrationChargePINDiode
412//
413// It retrieves the following variables from MCalibrationChargeCam:
414//
415// - fNumHiGainSamples
416// - fNumLoGainSamples
417//
418// It defines the PixId of every pixel in:
419//
420// - MCalibrationIntensityChargeCam
421// - MCalibrationChargeCam
422// - MCalibrationQECam
423//
424// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
425//
426// - MCalibrationChargePix
427// - MCalibrationQEPix
428//
429// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
430//
431// - MCalibrationChargeCam
432// - MCalibrationBlindPix (if existing)
433// - MCalibrationChargePINDiode (if existing)
434//
435Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
436{
437
438 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
439 if (!fGeom)
440 {
441 *fLog << err << "No MGeomCam found... aborting." << endl;
442 return kFALSE;
443 }
444
445 fIntensCam = (MCalibrationIntensityChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityChargeCam"));
446 if (fIntensCam)
447 *fLog << inf << "Found MCalibrationIntensityChargeCam ... " << endl;
448 else
449 {
450 fCam = (MCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
451 if (!fCam)
452 {
453 *fLog << err << "Cannot find MCalibrationChargeCam ... abort." << endl;
454 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
455 return kFALSE;
456 }
457 }
458
459 fHCam = (MHCalibrationChargeCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeCam"));
460 if (!fHCam)
461 {
462 *fLog << err << "Cannot find MHCalibrationChargeCam ... abort." << endl;
463 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
464 return kFALSE;
465 }
466
467 fIntensQE = (MCalibrationIntensityQECam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityQECam"));
468 if (fIntensQE)
469 *fLog << inf << "Found MCalibrationIntensityQECam ... " << endl;
470 else
471 {
472 fQECam = (MCalibrationQECam*)pList->FindObject(AddSerialNumber("MCalibrationQECam"));
473 if (!fQECam)
474 {
475 *fLog << err << "Cannot find MCalibrationQECam ... abort." << endl;
476 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationQECam before..." << endl;
477 return kFALSE;
478 }
479 }
480
481 fIntensBlind = (MCalibrationIntensityBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityBlindCam"));
482 if (fIntensBlind)
483 *fLog << inf << "Found MCalibrationIntensityBlindCam ... " << endl;
484 else
485 {
486 fBlindCam = (MCalibrationBlindCam*)pList->FindObject(AddSerialNumber("MCalibrationBlindCam"));
487 if (!fBlindCam)
488 {
489 *fLog << endl;
490 *fLog << warn << GetDescriptor()
491 << ": No MCalibrationBlindCam found... no Blind Pixel method! " << endl;
492 }
493 }
494
495 fHBlindCam = (MHCalibrationChargeBlindCam*)pList->FindObject(AddSerialNumber("MHCalibrationChargeBlindCam"));
496 if (!fHBlindCam)
497 {
498 *fLog << endl;
499 *fLog << warn << GetDescriptor()
500 << ": No MHCalibrationChargeBlindCam found... no Blind Pixel method! " << endl;
501 }
502
503 fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
504 if (fIntensBad)
505 *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
506 else
507 {
508 fBadPixels = (MBadPixelsCam*)pList->FindObject(AddSerialNumber("MBadPixelsCam"));
509 if (!fBadPixels)
510 {
511 *fLog << err << "Cannot find MBadPixelsCam ... abort." << endl;
512 return kFALSE;
513 }
514 }
515
516 //
517 // Optional Containers
518 //
519 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
520 if (!fPINDiode)
521 {
522 *fLog << endl;
523 *fLog << warn << GetDescriptor()
524 << ": MCalibrationChargePINDiode not found... no PIN Diode method! " << endl;
525 }
526
527 MCalibrationQECam *qecam = fIntensQE
528 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
529 MCalibrationChargeCam *chargecam = fIntensCam
530 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
531 MBadPixelsCam *badcam = fIntensBad
532 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
533
534 UInt_t npixels = fGeom->GetNumPixels();
535
536 for (UInt_t i=0; i<npixels; i++)
537 {
538
539 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
540 MCalibrationQEPix &pqe = (MCalibrationQEPix&) (*qecam) [i];
541 MBadPixelsPix &bad = (*badcam) [i];
542
543 if (bad.IsBad())
544 {
545 pix.SetExcluded();
546 pqe.SetExcluded();
547 continue;
548 }
549
550 if (IsDebug())
551 pix.SetDebug();
552 }
553
554 return kTRUE;
555}
556
557// ----------------------------------------------------------------------------------
558//
559// Set the correct colour to the charge containers
560//
561Int_t MCalibrationChargeCalc::Process()
562{
563
564 /*
565 if (IsInterlaced())
566 {
567 if (!(fTrigPattern->GetPrescaled() & MTriggerPattern::kCalibration))
568 {
569 if (IsDebug())
570 {
571 for (Int_t i=16; i>= 0; i--)
572 *fLog << err << (fTrigPattern->GetPrescaled() >> i & 1);
573 *fLog << endl;
574 }
575 return kTRUE;
576 }
577 else
578 {
579 if (IsDebug())
580 {
581 for (Int_t i=16; i>= 0; i--)
582 *fLog << inf << (fTrigPattern->GetPrescaled() >> i & 1);
583 *fLog << endl;
584 }
585 }
586 }
587 */
588 const MCalibrationCam::PulserColor_t col = fEvtHeader->GetPulserColor();
589
590 if (col == fPulserColor)
591 {
592 fNumProcessed++;
593 return kTRUE;
594 }
595
596 if (col == MCalibrationCam::kNONE)
597 return kTRUE;
598
599
600 //
601 // Now retrieve the colour and check if not various colours have been used
602 //
603 if (!fIntensCam)
604 {
605 if (fPulserColor != MCalibrationCam::kNONE)
606 {
607 *fLog << err << "Multiple colours used simultaneously in calibration file.";
608 *fLog << "Use Intenstiy Calibration for this case! " << endl;
609 return kFALSE;
610 /*
611 fHCam->Finalize();
612 fHBlindCam->Finalize();
613 fHCam->ResetHists();
614 fHBlindCam->ResetHists();
615 Finalize();
616 */
617 }
618 }
619
620 fPulserColor = col;
621
622 *fLog << inf << "Found new colour ... " << flush;
623
624 switch (col)
625 {
626 case MCalibrationCam::kGREEN: *fLog << "Green"; break;
627 case MCalibrationCam::kBLUE: *fLog << "Blue"; break;
628 case MCalibrationCam::kUV: *fLog << "UV"; break;
629 case MCalibrationCam::kCT1: *fLog << "CT1"; break;
630 default: break;
631 }
632
633 *fLog << inf << " with strength: " << fEvtHeader->GetPulserStrength() << endl;
634
635 if (fPINDiode)
636 fPINDiode->SetColor( fPulserColor );
637
638 fNumProcessed = 0;
639
640 return kTRUE;
641}
642
643// -----------------------------------------------------------------------
644//
645// Return if number of executions is null.
646//
647Int_t MCalibrationChargeCalc::PostProcess()
648{
649
650 if (GetNumExecutions()==0)
651 return kFALSE;
652
653 if (fPulserColor == MCalibrationCam::kNONE)
654 return kTRUE;
655
656 if (fNumProcessed == 0)
657 return kTRUE;
658
659 *fLog << endl;
660
661 return Finalize();
662}
663
664// -----------------------------------------------------------------------
665//
666// Return kTRUE if fPulserColor is kNONE
667//
668// First loop over pixels, average areas and sectors, call:
669// - FinalizePedestals()
670// - FinalizeCharges()
671// for every entry. Count number of valid pixels in loop and return kFALSE
672// if there are none (the "Michele check").
673//
674// Call FinalizeBadPixels()
675//
676// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
677//
678// Call FinalizeBlindCam()
679// Call FinalizePINDiode()
680//
681// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
682// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
683// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
684//
685// Call FinalizeUnsuitablePixels()
686//
687// Call MParContainer::SetReadyToSave() for fIntensCam, fCam, fQECam, fBadPixels and
688// fBlindCam and fPINDiode if they exist
689//
690// Print out some statistics
691//
692Int_t MCalibrationChargeCalc::Finalize()
693{
694
695 fNumHiGainSamples = fSignal->GetNumUsedHiGainFADCSlices();
696 fNumLoGainSamples = fSignal->GetNumUsedLoGainFADCSlices();
697
698 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
699 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
700
701 if (fPINDiode)
702 if (!fPINDiode->IsValid())
703 {
704 *fLog << warn << GetDescriptor()
705 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
706 fPINDiode = NULL;
707 }
708
709 MCalibrationBlindCam *blindcam = fIntensBlind
710 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
711 MCalibrationQECam *qecam = fIntensQE
712 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
713 MCalibrationChargeCam *chargecam = fIntensCam
714 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
715 MBadPixelsCam *badcam = fIntensBad
716 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
717
718 //
719 // First loop over pixels, call FinalizePedestals and FinalizeCharges
720 //
721 Int_t nvalid = 0;
722
723 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
724 {
725
726 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[pixid];
727 //
728 // Check if the pixel has been excluded from the fits
729 //
730 if (pix.IsExcluded())
731 continue;
732
733 MPedestalPix &ped = (*fPedestals)[pixid];
734 MBadPixelsPix &bad = (*badcam) [pixid];
735
736 const Int_t aidx = (*fGeom)[pixid].GetAidx();
737
738 FinalizePedestals(ped,pix,aidx);
739
740 if (FinalizeCharges(pix,bad,"pixel "))
741 nvalid++;
742 }
743
744 *fLog << endl;
745 //
746 // The Michele check ...
747 //
748 if (nvalid == 0)
749 {
750 if (!fIntensCam)
751 {
752 *fLog << warn << GetDescriptor() << ": All pixels have non-valid calibration. "
753 << "Did you forget to fill the histograms "
754 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
755 *fLog << warn << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
756 << "instead of a calibration run " << endl;
757 return kFALSE;
758 }
759 }
760
761 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
762 {
763
764 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
765 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(aidx);
766
767 FinalizePedestals(ped,pix,aidx);
768 FinalizeCharges(pix,
769 fIntensCam ? fIntensCam->GetAverageBadArea(aidx) : fCam->GetAverageBadArea(aidx),
770 "area id");
771 }
772
773 *fLog << endl;
774
775 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
776 {
777
778 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
779
780 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(sector);
781 FinalizePedestals(ped,pix, 0);
782 }
783
784 *fLog << endl;
785
786 //
787 // Finalize Bad Pixels
788 //
789 FinalizeBadPixels();
790
791 //
792 // Finalize F-Factor method
793 //
794 if (FinalizeFFactorMethod())
795 chargecam->SetFFactorMethodValid(kTRUE);
796 else
797 {
798 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
799 chargecam->SetFFactorMethodValid(kFALSE);
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 if (!pix.CalcConvFFactor())
1439 {
1440 *fLog << warn << GetDescriptor()
1441 << ": Could not calculate the Conv. FADC counts to Phes in pixel: "
1442 << Form(" %4i",pix.GetPixId())
1443 << endl;
1444 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1445 if (IsCheckDeviatingBehavior())
1446 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1447 }
1448 }
1449 }
1450
1451 return kTRUE;
1452}
1453
1454
1455
1456// ------------------------------------------------------------------------
1457//
1458// Returns kFALSE if pointer to MCalibrationBlindCam is NULL
1459//
1460// The check returns kFALSE if:
1461//
1462// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1463// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1464//
1465// Calls:
1466// - MCalibrationBlindPix::CalcFluxInsidePlexiglass()
1467//
1468Bool_t MCalibrationChargeCalc::FinalizeBlindCam()
1469{
1470
1471 MCalibrationBlindCam *blindcam = fIntensBlind
1472 ? (MCalibrationBlindCam*)fIntensBlind->GetCam() : fBlindCam;
1473
1474 if (!blindcam)
1475 return kFALSE;
1476
1477 Int_t nvalid = 0;
1478
1479 for (Int_t i=0; i<blindcam->GetSize(); i++)
1480 {
1481
1482 MCalibrationBlindPix &blindpix = (MCalibrationBlindPix&)(*blindcam)[i];
1483
1484 if (!blindpix.IsValid())
1485 continue;
1486
1487 const Float_t lambda = blindpix.GetLambda();
1488 const Float_t lambdaerr = blindpix.GetLambdaErr();
1489 const Float_t lambdacheck = blindpix.GetLambdaCheck();
1490
1491 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1492 {
1493 *fLog << warn << GetDescriptor()
1494 << Form("%s%4.2f%s%4.2f%s%4.2f%s%2i",": Lambda: ",lambda," and Lambda-Check: ",
1495 lambdacheck," differ by more than ",fLambdaCheckLimit," in the Blind Pixel Nr.",i)
1496 << endl;
1497 blindpix.SetValid(kFALSE);
1498 continue;
1499 }
1500
1501 if (lambdaerr > fLambdaErrLimit)
1502 {
1503 *fLog << warn << GetDescriptor()
1504 << Form("%s%4.2f%s%4.2f%s%2i",": Error of Fitted Lambda: ",lambdaerr," is greater than ",
1505 fLambdaErrLimit," in Blind Pixel Nr.",i) << endl;
1506 blindpix.SetValid(kFALSE);
1507 continue;
1508 }
1509
1510 if (!blindpix.CalcFluxInsidePlexiglass())
1511 {
1512 *fLog << warn << "Could not calculate the flux of photons from Blind Pixel Nr." << i << endl;
1513 blindpix.SetValid(kFALSE);
1514 continue;
1515 }
1516
1517 nvalid++;
1518 }
1519
1520 if (!nvalid)
1521 return kFALSE;
1522
1523 return kTRUE;
1524}
1525
1526// ------------------------------------------------------------------------
1527//
1528// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1529//
1530// The check returns kFALSE if:
1531//
1532// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1533// 2) PINDiode has a fit error smaller than fChargeErrLimit
1534// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1535// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1536//
1537// Calls:
1538// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1539//
1540Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1541{
1542
1543 if (!fPINDiode)
1544 return kFALSE;
1545
1546 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1547 {
1548 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1549 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1550 return kFALSE;
1551 }
1552
1553 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1554 {
1555 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1556 << fChargeErrLimit << " in PINDiode " << endl;
1557 return kFALSE;
1558 }
1559
1560 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1561 {
1562 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1563 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1564 return kFALSE;
1565 }
1566
1567 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1568 {
1569 *fLog << warn << GetDescriptor()
1570 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1571 return kFALSE;
1572 }
1573
1574
1575 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1576 {
1577 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1578 << "will skip PIN Diode Calibration " << endl;
1579 return kFALSE;
1580 }
1581
1582 return kTRUE;
1583}
1584
1585// ------------------------------------------------------------------------
1586//
1587// Calculate the average number of photons outside the plexiglass with the
1588// formula:
1589//
1590// av.Num.photons(area index) = av.Num.Phes(area index)
1591// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1592// / MCalibrationQEPix::GetPMTCollectionEff()
1593// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1594// / MCalibrationQECam::GetPlexiglassQE()
1595//
1596// Calculate the variance on the average number of photons assuming that the error on the
1597// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1598// values keeps it ordinary error since it is systematic.
1599//
1600// Loop over pixels:
1601//
1602// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1603// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1604//
1605// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1606// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1607//
1608// - Calculate the quantum efficiency with the formula:
1609//
1610// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1611//
1612// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1613//
1614// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1615// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1616//
1617// - Call MCalibrationQEPix::UpdateFFactorMethod()
1618//
1619void MCalibrationChargeCalc::FinalizeFFactorQECam()
1620{
1621
1622 if (fNumInnerFFactorMethodUsed < 2)
1623 {
1624 *fLog << warn << GetDescriptor()
1625 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1626 return;
1627 }
1628
1629 MCalibrationQECam *qecam = fIntensQE
1630 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1631 MCalibrationChargeCam *chargecam = fIntensCam
1632 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
1633 MBadPixelsCam *badcam = fIntensBad
1634 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1635
1636 MCalibrationChargePix &avpix = (MCalibrationChargePix&)chargecam->GetAverageArea(0);
1637 MCalibrationQEPix &qepix = (MCalibrationQEPix&) qecam->GetAverageArea(0);
1638
1639 const Float_t avphotons = avpix.GetPheFFactorMethod()
1640 / qepix.GetDefaultQE(fPulserColor)
1641 / qepix.GetPMTCollectionEff()
1642 / qepix.GetLightGuidesEff(fPulserColor)
1643 / qecam->GetPlexiglassQE();
1644
1645 const Float_t avphotrelvar = avpix.GetPheFFactorMethodRelVar()
1646 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1647 + qepix.GetPMTCollectionEffRelVar()
1648 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1649 + qecam->GetPlexiglassQERelVar();
1650
1651 const UInt_t nareas = fGeom->GetNumAreas();
1652
1653 //
1654 // Set the results in the MCalibrationChargeCam
1655 //
1656 chargecam->SetNumPhotonsFFactorMethod (avphotons);
1657
1658 if (avphotrelvar > 0.)
1659 chargecam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1660
1661 TArrayF lowlim (nareas);
1662 TArrayF upplim (nareas);
1663 TArrayD avffactorphotons (nareas);
1664 TArrayD avffactorphotvar (nareas);
1665 TArrayI numffactor (nareas);
1666
1667 const UInt_t npixels = fGeom->GetNumPixels();
1668
1669 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1670
1671 for (UInt_t i=0; i<npixels; i++)
1672 {
1673
1674 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1675 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1676 MBadPixelsPix &bad = (*badcam) [i];
1677
1678 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1679 continue;
1680
1681 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1682 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1683
1684 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1685
1686 qepix.SetQEFFactor ( qe , fPulserColor );
1687 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1688 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1689
1690 if (!qepix.UpdateFFactorMethod( qecam->GetPlexiglassQE() ))
1691 *fLog << warn << GetDescriptor()
1692 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1693
1694 //
1695 // The following pixels are those with deviating sigma, but otherwise OK,
1696 // probably those with stars during the pedestal run, but not the cal. run.
1697 //
1698 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1699 {
1700 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1701 if (IsCheckDeviatingBehavior())
1702 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1703 continue;
1704 }
1705
1706 const Int_t aidx = (*fGeom)[i].GetAidx();
1707 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1708
1709 camffactor.Fill(i,ffactor);
1710 camffactor.SetUsed(i);
1711
1712 avffactorphotons[aidx] += ffactor;
1713 avffactorphotvar[aidx] += ffactor*ffactor;
1714 numffactor[aidx]++;
1715 }
1716
1717 for (UInt_t i=0; i<nareas; i++)
1718 {
1719
1720 if (numffactor[i] == 0)
1721 {
1722 *fLog << warn << GetDescriptor() << ": No pixels with valid total F-Factor found "
1723 << "in area index: " << i << endl;
1724 continue;
1725 }
1726
1727 avffactorphotvar[i] = (avffactorphotvar[i] - avffactorphotons[i]*avffactorphotons[i]/numffactor[i])
1728 / (numffactor[i]-1.);
1729 avffactorphotons[i] = avffactorphotons[i] / numffactor[i];
1730
1731 if (avffactorphotvar[i] < 0.)
1732 {
1733 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of total F-Factor found "
1734 << "in area index: " << i << endl;
1735 continue;
1736 }
1737
1738 lowlim [i] = 1.; // Lowest known F-Factor of a PMT
1739 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1740
1741 TArrayI area(1);
1742 area[0] = i;
1743
1744 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1745 hist->Fit("gaus","Q");
1746 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1747 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1748 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1749
1750 if (IsDebug())
1751 camffactor.DrawClone();
1752
1753 if (ndf < 2)
1754 {
1755 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1756 << "in the camera with area index: " << i << endl;
1757 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1758 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1759 delete hist;
1760 continue;
1761 }
1762
1763 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1764
1765 if (prob < 0.001)
1766 {
1767 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1768 << "in the camera with area index: " << i << endl;
1769 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1770 << " is smaller than 0.001 " << endl;
1771 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1772 delete hist;
1773 continue;
1774 }
1775
1776 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1777 << "with area index #" << i << ": "
1778 << Form("%4.2f+-%4.2f",mean,sigma) << endl;
1779
1780 lowlim [i] = 1.;
1781 upplim [i] = mean + fFFactorErrLimit*sigma;
1782
1783 delete hist;
1784 }
1785
1786 *fLog << endl;
1787
1788 for (UInt_t i=0; i<npixels; i++)
1789 {
1790
1791 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1792 MBadPixelsPix &bad = (*badcam) [i];
1793
1794 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1795 continue;
1796
1797 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1798 const Int_t aidx = (*fGeom)[i].GetAidx();
1799
1800 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1801 {
1802 *fLog << warn << "Overall F-Factor "
1803 << Form("%5.2f",ffactor) << " out of range ["
1804 << Form("%5.2f,%5.2f",lowlim[aidx],upplim[aidx]) << "] Pixel " << i << endl;
1805
1806 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1807 if (IsCheckDeviatingBehavior())
1808 bad.SetUnsuitable ( MBadPixelsPix::kUnreliableRun );
1809 }
1810 }
1811
1812 for (UInt_t i=0; i<npixels; i++)
1813 {
1814
1815 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1816 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1817 MBadPixelsPix &bad = (*badcam) [i];
1818
1819 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1820 {
1821 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1822 pix.SetFFactorMethodValid(kFALSE);
1823 pix.SetExcluded();
1824 continue;
1825 }
1826 }
1827}
1828
1829
1830// ------------------------------------------------------------------------
1831//
1832// Loop over pixels:
1833//
1834// - Continue, if not MCalibrationBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1835// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1836//
1837// - Calculate the quantum efficiency with the formula:
1838//
1839// QE = Num.Phes / MCalibrationBlindPix::GetFluxInsidePlexiglass()
1840// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1841//
1842// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1843// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1844// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1845//
1846// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1847//
1848void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1849{
1850
1851
1852 MCalibrationBlindCam *blindcam = fIntensBlind
1853 ? (MCalibrationBlindCam*) fIntensBlind->GetCam(): fBlindCam;
1854 MBadPixelsCam *badcam = fIntensBad
1855 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
1856 MCalibrationQECam *qecam = fIntensQE
1857 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
1858 MCalibrationChargeCam *chargecam = fIntensCam
1859 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
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 // With the knowledge of the overall photon flux, calculate the
1954 // quantum efficiencies after the PIN Diode method
1955 //
1956 for (UInt_t i=0; i<npixels; i++)
1957 {
1958
1959 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
1960
1961 if (!fPINDiode)
1962 {
1963 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1964 continue;
1965 }
1966
1967 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
1968 {
1969 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1970 continue;
1971 }
1972
1973 MBadPixelsPix &bad = (*badcam) [i];
1974
1975 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1976 {
1977 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1978 continue;
1979 }
1980
1981 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i];
1982 MGeomPix &geo = (*fGeom) [i];
1983
1984 const Float_t qe = pix.GetPheFFactorMethod()
1985 / fPINDiode->GetFluxOutsidePlexiglass()
1986 / geo.GetA();
1987
1988 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
1989
1990 qepix.SetQEPINDiode ( qe , fPulserColor );
1991 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
1992 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
1993
1994 if (!qepix.UpdatePINDiodeMethod())
1995 *fLog << warn << GetDescriptor()
1996 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
1997 }
1998}
1999
2000// ------------------------------------------------------------------------
2001//
2002// Loop over pixels:
2003//
2004// - Call MCalibrationQEPix::UpdateCombinedMethod()
2005//
2006void MCalibrationChargeCalc::FinalizeCombinedQECam()
2007{
2008
2009 const UInt_t npixels = fGeom->GetNumPixels();
2010
2011 MCalibrationQECam *qecam = fIntensQE
2012 ? (MCalibrationQECam*) fIntensQE->GetCam() : fQECam;
2013 MBadPixelsCam *badcam = fIntensBad
2014 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
2015
2016 for (UInt_t i=0; i<npixels; i++)
2017 {
2018
2019 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*qecam) [i];
2020 MBadPixelsPix &bad = (*badcam) [i];
2021
2022 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
2023 {
2024 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
2025 continue;
2026 }
2027
2028 qepix.UpdateCombinedMethod();
2029 }
2030}
2031
2032// -----------------------------------------------------------------------------------------------
2033//
2034// - Print out statistics about BadPixels of type UnsuitableType_t
2035// - store numbers of bad pixels of each type in fCam or fIntensCam
2036//
2037void MCalibrationChargeCalc::FinalizeUnsuitablePixels()
2038{
2039
2040 *fLog << inf << endl;
2041 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
2042 *fLog << dec << setfill(' ');
2043
2044 const Int_t nareas = fGeom->GetNumAreas();
2045
2046 TArrayI counts(nareas);
2047
2048 MBadPixelsCam *badcam = fIntensBad
2049 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2050 MCalibrationChargeCam *chargecam = fIntensCam
2051 ? (MCalibrationChargeCam*)fIntensCam->GetCam() : fCam;
2052
2053 for (Int_t i=0; i<badcam->GetSize(); i++)
2054 {
2055 MBadPixelsPix &bad = (*badcam)[i];
2056 if (!bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
2057 {
2058 const Int_t aidx = (*fGeom)[i].GetAidx();
2059 counts[aidx]++;
2060 }
2061 }
2062
2063 if (fGeom->InheritsFrom("MGeomCamMagic"))
2064 *fLog << " " << setw(7) << "Successfully calibrated Pixels: "
2065 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2066
2067 counts.Reset();
2068
2069 for (Int_t i=0; i<badcam->GetSize(); i++)
2070 {
2071 MBadPixelsPix &bad = (*badcam)[i];
2072
2073 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
2074 {
2075 const Int_t aidx = (*fGeom)[i].GetAidx();
2076 counts[aidx]++;
2077 }
2078 }
2079
2080 for (Int_t aidx=0; aidx<nareas; aidx++)
2081 chargecam->SetNumUnsuitable(counts[aidx], aidx);
2082
2083 if (fGeom->InheritsFrom("MGeomCamMagic"))
2084 *fLog << " " << setw(7) << "Uncalibrated Pixels: "
2085 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2086
2087 counts.Reset();
2088
2089 for (Int_t i=0; i<badcam->GetSize(); i++)
2090 {
2091
2092 MBadPixelsPix &bad = (*badcam)[i];
2093
2094 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
2095 {
2096 const Int_t aidx = (*fGeom)[i].GetAidx();
2097 counts[aidx]++;
2098 }
2099 }
2100
2101 for (Int_t aidx=0; aidx<nareas; aidx++)
2102 chargecam->SetNumUnreliable(counts[aidx], aidx);
2103
2104 *fLog << " " << setw(7) << "Unreliable Pixels: "
2105 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
2106
2107}
2108
2109// -----------------------------------------------------------------------------------------------
2110//
2111// Print out statistics about BadPixels of type UncalibratedType_t
2112//
2113void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
2114{
2115
2116 UInt_t countinner = 0;
2117 UInt_t countouter = 0;
2118
2119 MBadPixelsCam *badcam = fIntensBad
2120 ? (MBadPixelsCam*)fIntensBad->GetCam() : fBadPixels;
2121
2122 for (Int_t i=0; i<badcam->GetSize(); i++)
2123 {
2124 MBadPixelsPix &bad = (*badcam)[i];
2125
2126 if (bad.IsUncalibrated(typ))
2127 {
2128 if (fGeom->GetPixRatio(i) == 1.)
2129 countinner++;
2130 else
2131 countouter++;
2132 }
2133 }
2134
2135 *fLog << " " << setw(7) << text
2136 << Form("%s%3i%s%3i","Inner: ",countinner," Outer: ",countouter) << endl;
2137}
2138
2139// --------------------------------------------------------------------------
2140//
2141// Set the path for output file
2142//
2143void MCalibrationChargeCalc::SetOutputPath(TString path)
2144{
2145 fOutputPath = path;
2146 if (fOutputPath.EndsWith("/"))
2147 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
2148}
2149
2150// --------------------------------------------------------------------------
2151//
2152// Set the output file
2153//
2154void MCalibrationChargeCalc::SetOutputFile(TString file)
2155{
2156 fOutputFile = file;
2157}
2158
2159// --------------------------------------------------------------------------
2160//
2161// Get the output file
2162//
2163const char* MCalibrationChargeCalc::GetOutputFile()
2164{
2165 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
2166}
2167
2168// --------------------------------------------------------------------------
2169//
2170// Read the environment for the following data members:
2171// - fChargeLimit
2172// - fChargeErrLimit
2173// - fChargeRelErrLimit
2174// - fDebug
2175// - fFFactorErrLimit
2176// - fLambdaErrLimit
2177// - fLambdaCheckErrLimit
2178// - fPheErrLimit
2179//
2180Int_t MCalibrationChargeCalc::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
2181{
2182
2183 Bool_t rc = kFALSE;
2184 if (IsEnvDefined(env, prefix, "ChargeLimit", print))
2185 {
2186 SetChargeLimit(GetEnvValue(env, prefix, "ChargeLimit", fChargeLimit));
2187 rc = kTRUE;
2188 }
2189 if (IsEnvDefined(env, prefix, "ChargeErrLimit", print))
2190 {
2191 SetChargeErrLimit(GetEnvValue(env, prefix, "ChargeErrLimit", fChargeErrLimit));
2192 rc = kTRUE;
2193 }
2194 if (IsEnvDefined(env, prefix, "ChargeRelErrLimit", print))
2195 {
2196 SetChargeRelErrLimit(GetEnvValue(env, prefix, "ChargeRelErrLimit", fChargeRelErrLimit));
2197 rc = kTRUE;
2198 }
2199 if (IsEnvDefined(env, prefix, "Debug", print))
2200 {
2201 SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
2202 rc = kTRUE;
2203 }
2204 if (IsEnvDefined(env, prefix, "FFactorErrLimit", print))
2205 {
2206 SetFFactorErrLimit(GetEnvValue(env, prefix, "FFactorErrLimit", fFFactorErrLimit));
2207 rc = kTRUE;
2208 }
2209 if (IsEnvDefined(env, prefix, "LambdaErrLimit", print))
2210 {
2211 SetLambdaErrLimit(GetEnvValue(env, prefix, "LambdaErrLimit", fLambdaErrLimit));
2212 rc = kTRUE;
2213 }
2214 if (IsEnvDefined(env, prefix, "LambdaCheckLimit", print))
2215 {
2216 SetLambdaCheckLimit(GetEnvValue(env, prefix, "LambdaCheckLimit", fLambdaCheckLimit));
2217 rc = kTRUE;
2218 }
2219 if (IsEnvDefined(env, prefix, "PheErrLimit", print))
2220 {
2221 SetPheErrLimit(GetEnvValue(env, prefix, "PheErrLimit", fPheErrLimit));
2222 rc = kTRUE;
2223 }
2224 if (IsEnvDefined(env, prefix, "CheckDeadPixels", print))
2225 {
2226 SetCheckDeadPixels(GetEnvValue(env, prefix, "CheckDeadPixels", IsCheckDeadPixels()));
2227 rc = kTRUE;
2228 }
2229 if (IsEnvDefined(env, prefix, "CheckDeviatingBehavior", print))
2230 {
2231 SetCheckDeviatingBehavior(GetEnvValue(env, prefix, "CheckDeviatingBehavior", IsCheckDeviatingBehavior()));
2232 rc = kTRUE;
2233 }
2234 if (IsEnvDefined(env, prefix, "CheckExtractionWindow", print))
2235 {
2236 SetCheckExtractionWindow(GetEnvValue(env, prefix, "CheckExtractionWindow", IsCheckExtractionWindow()));
2237 rc = kTRUE;
2238 }
2239 if (IsEnvDefined(env, prefix, "CheckHistOverflow", print))
2240 {
2241 SetCheckHistOverflow(GetEnvValue(env, prefix, "CheckHistOverflow", IsCheckHistOverflow()));
2242 rc = kTRUE;
2243 }
2244 if (IsEnvDefined(env, prefix, "CheckOscillations", print))
2245 {
2246 SetCheckOscillations(GetEnvValue(env, prefix, "CheckOscillations", IsCheckOscillations()));
2247 rc = kTRUE;
2248 }
2249
2250 return rc;
2251}
2252
Note: See TracBrowser for help on using the repository browser.