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

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