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

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