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

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