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

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