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

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