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

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