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

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