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

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