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

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