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

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