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

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