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

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