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

Last change on this file since 4322 was 4322, checked in by gaug, 20 years ago
*** empty log message ***
File size: 60.7 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// - FinalizePINDiode()
52// - FinalizeFFactorQECam()
53// - FinalizeBlindPixelQECam()
54// - FinalizePINDiodeQECam()
55//
56// Input Containers:
57// MCalibrationChargeCam
58// MCalibrationChargeBlindPix
59// MCalibrationChargePINDiode
60// MCalibrationQECam
61// MPedestalCam
62// MBadPixelsCam
63// MGeomCam
64// MTime
65//
66// Output Containers:
67// MCalibrationChargeCam
68// MCalibrationChargeBlindPix
69// MCalibrationChargePINDiode
70// MCalibrationQECam
71// MBadPixelsCam
72//
73//
74// Preliminary description of the calibration in photons (email from 12/02/04)
75//
76// Why calibrating in photons:
77// ===========================
78//
79// At the Barcelona meeting in 2002, we decided to calibrate the camera in
80// photons. This for the following reasons:
81//
82// * The physical quantity arriving at the camera are photons. This is
83// the direct physical information from the air shower. The photons
84// have a flux and a spectrum.
85//
86// * The photon fluxes depend mostly on the shower energy (with
87// corrections deriving from the observation conditions), while the photon
88// spectra depend mostly on the observation conditions: zenith angle,
89// quality of the air, also the impact parameter of the shower.
90//
91// * The photomultiplier, in turn, has different response properties
92// (quantum efficiencies) for photons of different colour. (Moreover,
93// different pixels have slightly different quantum efficiencies).
94// The resulting number of photo-electrons is then amplified (linearly)
95// with respect to the photo-electron flux.
96//
97// * In the ideal case, one would like to disentagle the effects
98// of the observation conditions from the primary particle energy (which
99// one likes to measure). To do so, one needs:
100//
101// 1) A reliable calibration relating the FADC counts to the photo-electron
102// flux -> This is accomplished with the F-Factor method.
103//
104// 2) A reliable calibration of the wavelength-dependent quantum efficiency
105// -> This is accomplished with the combination of the three methods,
106// together with QE-measurements performed by David in order to do
107// the interpolation.
108//
109// 3) A reliable calibration of the observation conditions. This means:
110// - Tracing the atmospheric conditions -> LIDAR
111// - Tracing the observation zenith angle -> Drive System
112//
113// 4) Some knowlegde about the impact parameter:
114// - This is the only part which cannot be accomplished well with a
115// single telescope. We would thus need to convolute the spectrum
116// over the distribution of impact parameters.
117//
118//
119// How an ideal calibration would look like:
120// =========================================
121//
122// We know from the combined PIN-Diode and Blind-Pixel Method the response of
123// each pixel to well-measured light fluxes in three representative
124// wavelengths (green, blue, UV). We also know the response to these light
125// fluxes in photo-electrons. Thus, we can derive:
126//
127// - conversion factors to photo-electrons
128// - conversion factors to photons in three wavelengths.
129//
130// Together with David's measurements and some MC-simulation, we should be
131// able to derive tables for typical Cherenkov-photon spectra - convoluted
132// with the impact parameters and depending on the athmospheric conditions
133// and the zenith angle (the "outer parameters").
134//
135// From these tables we can create "calibration tables" containing some
136// effective quantum efficiency depending on these outer parameters and which
137// are different for each pixel.
138//
139// In an ideal MCalibrate, one would thus have to convert first the FADC
140// slices to Photo-electrons and then, depending on the outer parameters,
141// look up the effective quantum efficiency and get the mean number of
142// photons which is then used for the further analysis.
143//
144// How the (first) MAGIC calibration should look like:
145// ===================================================
146//
147// For the moment, we have only one reliable calibration method, although
148// with very large systematic errors. This is the F-Factor method. Knowing
149// that the light is uniform over the whole camera (which I would not at all
150// guarantee in the case of the CT1 pulser), one could in principle already
151// perform a relative calibration of the quantum efficiencies in the UV.
152// However, the spread in QE at UV is about 10-15% (according to the plot
153// that Abelardo sent around last time. The spread in photo-electrons is 15%
154// for the inner pixels, but much larger (40%) for the outer ones.
155//
156// I'm not sure if we can already say that we have measured the relative
157// difference in quantum efficiency for the inner pixels and produce a first
158// QE-table for each pixel. To so, I would rather check in other wavelengths
159// (which we can do in about one-two weeks when the optical transmission of
160// the calibration trigger is installed).
161//
162// Thus, for the moment being, I would join Thomas proposal to calibrate in
163// photo-electrons and apply one stupid average quantum efficiency for all
164// pixels. This keeping in mind that we will have much preciser information
165// in about one to two weeks.
166//
167//
168// What MCalibrate should calculate and what should be stored:
169// ===========================================================
170//
171// It is clear that in the end, MCerPhotEvt will store photons.
172// MCalibrationCam stores the conversionfactors to photo-electrons and also
173// some tables of how to apply the conversion to photons, given the outer
174// parameters. This is not yet implemented and not even discussed.
175//
176// To start, I would suggest that we define the "average quantum efficiency"
177// (maybe something like 25+-3%) and apply them equally to all
178// photo-electrons. Later, this average factor can be easily replaced by a
179// pixel-dependent factor and later by a (pixel-dependent) table.
180//
181//
182//
183//////////////////////////////////////////////////////////////////////////////
184#include "MCalibrationChargeCalc.h"
185
186#include <TSystem.h>
187#include <TH1.h>
188#include <TF1.h>
189
190#include "MLog.h"
191#include "MLogManip.h"
192
193#include "MParList.h"
194
195#include "MRawRunHeader.h"
196#include "MRawEvtPixelIter.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 "MCalibrationChargeCam.h"
206#include "MCalibrationChargePix.h"
207#include "MCalibrationChargePINDiode.h"
208#include "MCalibrationChargeBlindPix.h"
209
210#include "MExtractedSignalCam.h"
211#include "MExtractedSignalPix.h"
212#include "MExtractedSignalBlindPixel.h"
213#include "MExtractedSignalPINDiode.h"
214
215#include "MBadPixelsCam.h"
216#include "MBadPixelsPix.h"
217
218#include "MCalibrationQECam.h"
219#include "MCalibrationQEPix.h"
220
221#include "MCalibrationCam.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 = 3.5;
233const Float_t MCalibrationChargeCalc::fgFFactorErrLimit = 4.5;
234// --------------------------------------------------------------------------
235//
236// Default constructor.
237//
238// Sets all pointers 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//
258// Calls:
259// - Clear()
260//
261MCalibrationChargeCalc::MCalibrationChargeCalc(const char *name, const char *title)
262 : fBadPixels(NULL), fCam(NULL), fBlindPixel(NULL), fPINDiode(NULL),
263 fQECam(NULL), fGeom(NULL), fPedestals(NULL), fEvtTime(NULL)
264{
265
266 fName = name ? name : "MCalibrationChargeCalc";
267 fTitle = title ? title : "Task to calculate the calibration constants and MCalibrationCam ";
268
269 AddToBranchList("MRawEvtData.fHiGainPixId");
270 AddToBranchList("MRawEvtData.fLoGainPixId");
271 AddToBranchList("MRawEvtData.fHiGainFadcSamples");
272 AddToBranchList("MRawEvtData.fLoGainFadcSamples");
273
274 SetChargeLimit();
275 SetChargeErrLimit();
276 SetChargeRelErrLimit();
277 SetFFactorErrLimit();
278 SetLambdaCheckLimit();
279 SetLambdaErrLimit();
280 SetPheErrLimit();
281 SetPulserColor(MCalibrationCam::kNONE);
282 SetOutputPath();
283 SetOutputFile();
284
285 Clear();
286
287}
288
289// --------------------------------------------------------------------------
290//
291// Sets:
292// - all variables to 0.,
293// - all flags to kFALSE
294//
295void MCalibrationChargeCalc::Clear(const Option_t *o)
296{
297
298 fNumHiGainSamples = 0.;
299 fNumLoGainSamples = 0.;
300 fSqrtHiGainSamples = 0.;
301 fSqrtLoGainSamples = 0.;
302 fNumInnerFFactorMethodUsed = 0;
303 SkipHiLoGainCalibration( kFALSE );
304}
305
306
307// -----------------------------------------------------------------------------------
308//
309// The following container are searched for and execution aborted if not in MParList:
310// - MPedestalCam
311//
312// The following containers are searched and created if they were not found:
313//
314// - MCalibrationQECam
315// - MBadPixelsCam
316//
317// The following output containers are only searched, but not created. If they
318// cannot be found, the corresponding calibration part is only skipped.
319//
320// - MTime
321//
322Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
323{
324
325 //
326 // Containers that have to be there.
327 //
328 fPedestals = (MPedestalCam*)pList->FindObject("MPedestalCam");
329 if (!fPedestals)
330 {
331 *fLog << err << "MPedestalCam not found... aborting" << endl;
332 return kFALSE;
333 }
334
335 //
336 // Containers that are created in case that they are not there.
337 //
338 fQECam = (MCalibrationQECam*)pList->FindCreateObj("MCalibrationQECam");
339 if (!fQECam)
340 {
341 *fLog << err << "Cannot find nor create MCalibrationQECam... aborting" << endl;
342 return kFALSE;
343 }
344
345 fBadPixels = (MBadPixelsCam*)pList->FindCreateObj("MBadPixelsCam");
346 if (!fBadPixels)
347 {
348 *fLog << err << "Could not find or create MBadPixelsCam ... aborting." << endl;
349 return kFALSE;
350 }
351
352
353 fEvtTime = (MTime*)pList->FindObject("MTime");
354
355 //
356 // Check the pulser colour --> FIXME: this solution is only valid until the arrival of the DM's
357 //
358 if (fPulserColor == MCalibrationCam::kNONE)
359 {
360 *fLog << endl;
361 *fLog << err << GetDescriptor()
362 << ": No Pulser colour has been chosen. Since the installation of the IFAE pulser box,"
363 << " you HAVE to provide the LEDs colour, otherwise there is no calibration. " << endl;
364 *fLog << "See e.g. the macro calibration.C " << endl;
365 return kFALSE;
366 }
367
368 return kTRUE;
369}
370
371
372// --------------------------------------------------------------------------
373//
374// Search for the following input containers and abort if not existing:
375// - MGeomCam
376// - MCalibrationChargeCam
377//
378// Search for the following input containers and give a warning if not existing:
379// - MCalibrationChargeBlindPix
380// - MCalibrationChargePINDiode
381//
382// It retrieves the following variables from MCalibrationChargeCam:
383//
384// - fNumHiGainSamples
385// - fNumLoGainSamples
386//
387// It defines the PixId of every pixel in:
388//
389// - MCalibrationChargeCam
390// - MCalibrationQECam
391//
392// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
393//
394// - MCalibrationChargePix
395// - MCalibrationQEPix
396//
397// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
398//
399// - MCalibrationChargeCam
400// - MCalibrationChargeBlindPix (if existing)
401// - MCalibrationChargePINDiode (if existing)
402//
403Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
404{
405
406 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
407 if (!fGeom)
408 {
409 *fLog << err << "No MGeomCam found... aborting." << endl;
410 return kFALSE;
411 }
412
413 fCam = (MCalibrationChargeCam*)pList->FindObject("MCalibrationChargeCam");
414 if (!fCam)
415 {
416 *fLog << err << "Cannot find MCalibrationChargeCam... aborting" << endl;
417 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
418 return kFALSE;
419 }
420
421 //
422 // Optional Containers
423 //
424 fBlindPixel = (MCalibrationChargeBlindPix*)pList->FindObject("MCalibrationChargeBlindPix");
425 if (!fBlindPixel)
426 {
427 *fLog << endl;
428 *fLog << warn << GetDescriptor()
429 << ": MCalibrationChargeBlindPix not found... no Blind Pixel method! " << endl;
430 }
431
432 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
433 if (!fPINDiode)
434 {
435 *fLog << endl;
436 *fLog << warn << GetDescriptor()
437 << ": MCalibrationChargePINDiode not found... no PIN Diode method! " << endl;
438 }
439
440
441 //
442 // Initialize the pulser colours
443 //
444 if (fCam->GetPulserColor() == MCalibrationCam::kNONE)
445 {
446 fCam->SetPulserColor( fPulserColor );
447
448 if (fBlindPixel)
449 fBlindPixel->SetColor( fPulserColor );
450
451 if (fPINDiode)
452 fPINDiode->SetColor( fPulserColor );
453 }
454
455 if (fPulserColor != fCam->GetPulserColor())
456 {
457 *fLog << err << GetDescriptor()
458 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargeCam" << endl;
459 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
460 return kFALSE;
461 }
462
463 if (fBlindPixel)
464 if (fPulserColor != fBlindPixel->GetColor())
465 {
466 *fLog << err << GetDescriptor()
467 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargeBlindPix." << endl;
468 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
469 return kFALSE;
470 }
471
472 if (fPINDiode)
473 if (fPulserColor != fPINDiode->GetColor())
474 {
475 *fLog << err << GetDescriptor()
476 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargePINDiode." << endl;
477 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
478 return kFALSE;
479 }
480
481
482 fNumHiGainSamples = fCam->GetNumHiGainFADCSlices();
483 fNumLoGainSamples = fCam->GetNumLoGainFADCSlices();
484 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
485 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
486
487 UInt_t npixels = fGeom->GetNumPixels();
488
489 for (UInt_t i=0; i<npixels; i++)
490 {
491
492 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam) [i];
493 MCalibrationQEPix &pqe = (MCalibrationQEPix&) (*fQECam)[i];
494 MBadPixelsPix &bad = (*fBadPixels)[i];
495
496 pix.SetPixId(i);
497 pqe.SetPixId(i);
498
499 if (bad.IsBad())
500 {
501 pix.SetExcluded();
502 pqe.SetExcluded();
503 continue;
504 }
505
506 }
507
508 return kTRUE;
509}
510
511// ----------------------------------------------------------------------------------
512//
513// Nothing to be done in Process, but have a look at MHCalibrationChargeCam, instead
514//
515Int_t MCalibrationChargeCalc::Process()
516{
517 return kTRUE;
518}
519
520// -----------------------------------------------------------------------
521//
522// Return if number of executions is null.
523//
524// First loop over pixels, average areas and sectors, call:
525// - FinalizePedestals()
526// - FinalizeCharges()
527// for every entry. Count number of valid pixels in loop and return kFALSE
528// if there are none (the "Michele check").
529//
530// Call FinalizeBadPixels()
531//
532// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
533//
534// Call FinalizeBlindPixel()
535// Call FinalizePINDiode()
536//
537// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
538// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
539// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
540//
541// Call FinalizeUnsuitablePixels()
542//
543// Call MParContainer::SetReadyToSave() for fCam, fQECam, fBadPixels and
544// fBlindPixel and fPINDiode if they exist
545//
546// Print out some statistics
547//
548Int_t MCalibrationChargeCalc::PostProcess()
549{
550
551 if (GetNumExecutions()==0)
552 return kFALSE;
553
554 if (fPINDiode)
555 if (!fPINDiode->IsValid())
556 {
557 *fLog << warn << GetDescriptor()
558 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
559 fPINDiode = NULL;
560 }
561
562 if (fBlindPixel)
563 if (!fBlindPixel->IsValid())
564 {
565 *fLog << warn << GetDescriptor()
566 << ": MCalibrationChargeBlindPix is declared not valid... no Blind Pixel method! " << endl;
567 fBlindPixel = NULL;
568 }
569
570 //
571 // First loop over pixels, call FinalizePedestals and FinalizeCharges
572 //
573 Int_t nvalid = 0;
574
575 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
576 {
577
578 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[pixid];
579 //
580 // Check if the pixel has been excluded from the fits
581 //
582 if (pix.IsExcluded())
583 continue;
584
585 MPedestalPix &ped = (*fPedestals)[pixid];
586 MBadPixelsPix &bad = (*fBadPixels)[pixid];
587 const Int_t aidx = (*fGeom)[pixid].GetAidx();
588
589 FinalizePedestals(ped,pix,aidx);
590
591 if (FinalizeCharges(pix,bad))
592 nvalid++;
593 }
594
595 //
596 // The Michele check ...
597 //
598 if (nvalid == 0)
599 {
600 *fLog << err << GetDescriptor() << ": All pixels have non-valid calibration. "
601 << "Did you forget to fill the histograms "
602 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
603 *fLog << err << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
604 << "instead of a calibration run " << endl;
605 return kFALSE;
606 }
607
608 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
609 {
610
611 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
612 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
613
614 FinalizePedestals(ped,pix,aidx);
615 FinalizeCharges(pix, fCam->GetAverageBadArea(aidx));
616 }
617
618 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
619 {
620
621 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
622
623 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
624 FinalizePedestals(ped,pix, 0);
625 FinalizeCharges(pix, fCam->GetAverageBadSector(sector));
626 }
627
628 //
629 // Finalize Bad Pixels
630 //
631 FinalizeBadPixels();
632
633 //
634 // Finalize F-Factor method
635 //
636 if (!FinalizeFFactorMethod())
637 {
638 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
639 fCam->SetFFactorMethodValid(kFALSE);
640 return kFALSE;
641 }
642 else
643 fCam->SetFFactorMethodValid(kTRUE);
644
645 //
646 // Finalize Blind Pixel
647 //
648 if (FinalizeBlindPixel())
649 fQECam->SetBlindPixelMethodValid(kTRUE);
650 else
651 fQECam->SetBlindPixelMethodValid(kFALSE);
652
653 //
654 // Finalize PIN Diode
655 //
656 if (FinalizePINDiode())
657 fQECam->SetPINDiodeMethodValid(kTRUE);
658 else
659 fQECam->SetPINDiodeMethodValid(kFALSE);
660
661 //
662 // Finalize QE Cam
663 //
664 FinalizeFFactorQECam();
665 FinalizeBlindPixelQECam();
666 FinalizePINDiodeQECam();
667
668 //
669 // Re-direct the output to an ascii-file from now on:
670 //
671 MLog asciilog;
672 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
673 SetLogStream(&asciilog);
674 //
675 // Finalize calibration statistics
676 //
677 FinalizeUnsuitablePixels();
678
679 fCam ->SetReadyToSave();
680 fQECam ->SetReadyToSave();
681 fBadPixels->SetReadyToSave();
682
683 if (fBlindPixel)
684 fBlindPixel->SetReadyToSave();
685 if (fPINDiode)
686 fPINDiode->SetReadyToSave();
687
688 *fLog << inf << endl;
689 *fLog << GetDescriptor() << ": Fatal errors statistics:" << endl;
690
691 PrintUncalibrated(MBadPixelsPix::kChargeIsPedestal,
692 Form("%s%2.1f%s","Signal less than ",fChargeLimit," Pedestal RMS: "));
693 PrintUncalibrated(MBadPixelsPix::kChargeRelErrNotValid,
694 Form("%s%2.1f%s","Signal Error bigger than ",fChargeRelErrLimit," times Mean Signal: "));
695 PrintUncalibrated(MBadPixelsPix::kChargeSigmaNotValid,
696 "Signal Sigma smaller than Pedestal RMS: ");
697 PrintUncalibrated(MBadPixelsPix::kLoGainSaturation,
698 "Pixels with Low Gain Saturation: ");
699 PrintUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin,
700 Form("%s%2.1f%s","Mean Abs. Arr. Time in First ",1.," Bin(s): "));
701 PrintUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins,
702 Form("%s%2.1f%s","Mean Abs. Arr. Time in Last ",2.," Bin(s): "));
703 PrintUncalibrated(MBadPixelsPix::kDeviatingNumPhes,
704 "Pixels with deviating number of phes: ");
705 PrintUncalibrated(MBadPixelsPix::kDeviatingFFactor,
706 "Pixels with deviating F-Factor: ");
707
708 *fLog << inf << endl;
709 *fLog << GetDescriptor() << ": Unreliable errors statistics:" << endl;
710
711 PrintUncalibrated(MBadPixelsPix::kHiGainOscillating,
712 "Pixels with changing Hi Gain signal over time: ");
713 PrintUncalibrated(MBadPixelsPix::kLoGainOscillating,
714 "Pixels with changing Lo Gain signal over time: ");
715 PrintUncalibrated(MBadPixelsPix::kHiGainNotFitted,
716 "Pixels with unsuccesful Gauss fit to the Hi Gain: ");
717 PrintUncalibrated(MBadPixelsPix::kLoGainNotFitted,
718 "Pixels with unsuccesful Gauss fit to the Lo Gain: ");
719
720 SetLogStream(&gLog);
721
722 return kTRUE;
723}
724
725// ----------------------------------------------------------------------------------
726//
727// Retrieves pedestal and pedestal RMS from MPedestalPix
728// Retrieves total entries from MPedestalCam
729// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
730// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
731//
732// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
733// - MCalibrationChargePix::CalcLoGainPedestal()
734//
735void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
736{
737
738 //
739 // get the pedestals
740 //
741 const Float_t pedes = ped.GetPedestal();
742 const Float_t prms = ped.GetPedestalRms();
743 const Float_t num = TMath::Sqrt((Float_t)fPedestals->GetTotalEntries());
744
745 //
746 // set them in the calibration camera
747 //
748 if (cal.IsHiGainSaturation())
749 {
750 cal.SetPedestal(pedes* fNumLoGainSamples,
751 prms * fSqrtLoGainSamples,
752 prms * fNumLoGainSamples / num);
753 cal.CalcLoGainPedestal((Float_t)fNumLoGainSamples, aidx);
754 }
755 else
756 {
757 cal.SetPedestal(pedes* fNumHiGainSamples,
758 prms * fSqrtHiGainSamples,
759 prms * fNumHiGainSamples / num);
760 }
761
762}
763
764// ----------------------------------------------------------------------------------------------------
765//
766// Check fit results validity. Bad Pixels flags are set if:
767//
768// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
769// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
770// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
771// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
772// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
773//
774// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
775//
776// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
777// if not succesful.
778//
779// Calls MCalibrationChargePix::CalcFFactorMethod() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
780// if not succesful.
781//
782Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad)
783{
784
785 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
786 return kFALSE;
787
788 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
789 {
790 *fLog << warn << GetDescriptor() << ": Fitted Charge: " << cal.GetMean() << " is smaller than "
791 << fChargeLimit << " Pedestal RMS: " << cal.GetPedRms() << " in Pixel " << cal.GetPixId() << endl;
792 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
793 }
794
795 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
796 {
797 *fLog << warn << GetDescriptor() << ": Fitted Charge: " << cal.GetMean() << " is smaller than "
798 << fChargeRelErrLimit << "* its error: " << cal.GetMeanErr()
799 << " in Pixel " << cal.GetPixId() << endl;
800 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
801 }
802
803 if (cal.GetSigma() < cal.GetPedRms())
804 {
805 *fLog << warn << GetDescriptor() << ": Sigma of Fitted Charge: " << cal.GetSigma()
806 << " smaller than Pedestal RMS: " << cal.GetPedRms() << " in Pixel " << cal.GetPixId() << endl;
807 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
808 }
809
810 if (!cal.CalcReducedSigma())
811 {
812 *fLog << warn << GetDescriptor()
813 << ": Could not calculate reduced sigmas of pixel: " << cal.GetPixId() << endl;
814 bad.SetUncalibrated(MBadPixelsPix::kChargeIsPedestal);
815 return kFALSE;
816 }
817
818 if (!cal.CalcFFactorMethod())
819 {
820 *fLog << warn << GetDescriptor()
821 << ": Could not calculate F-Factor of pixel: " << cal.GetPixId() << endl;
822 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
823 return kFALSE;
824 }
825
826 return kTRUE;
827}
828
829
830
831// -----------------------------------------------------------------------------------
832//
833// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
834// - MBadPixelsPix::kChargeIsPedestal
835// - MBadPixelsPix::kChargeRelErrNotValid
836// - MBadPixelsPix::kChargeSigmaNotValid
837// - MBadPixelsPix::kMeanTimeInFirstBin
838// - MBadPixelsPix::kMeanTimeInLast2Bins
839// - MBadPixelsPix::kDeviatingNumPhes
840//
841// - Call MCalibrationPix::SetExcluded() for the bad pixels
842//
843void MCalibrationChargeCalc::FinalizeBadPixels()
844{
845
846 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
847 {
848
849 MBadPixelsPix &bad = (*fBadPixels)[i];
850 MCalibrationPix &pix = (*fCam)[i];
851
852 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
853 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
854
855 if (bad.IsUncalibrated( MBadPixelsPix::kChargeErrNotValid ))
856 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
857
858 if (bad.IsUncalibrated( MBadPixelsPix::kChargeRelErrNotValid ))
859 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
860
861 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
862 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
863
864 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
865 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
866
867 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
868 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
869
870 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
871 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
872
873 if (bad.IsUnsuitable( MBadPixelsPix::kUnsuitableRun ))
874 pix.SetExcluded();
875
876 }
877}
878
879// ------------------------------------------------------------------------
880//
881//
882// First loop: Calculate a mean and mean RMS of photo-electrons per area index
883// Include only pixels which are not MBadPixelsPix::kUnsuitableRun or
884// MBadPixelsPix::kUnreliableRun (see FinalizeBadPixels()) and set
885// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
886//
887// Second loop: Get weighted mean number of photo-electrons and its RMS including
888// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
889// and further exclude those deviating by more than fPheErrLimit mean
890// sigmas from the mean (obtained in first loop). Set
891// MBadPixelsPix::kDeviatingNumPhes if excluded.
892//
893// Set weighted mean and variance of photo-electrons per area index in:
894// average area pixels of MCalibrationChargeCam (obtained from:
895// MCalibrationChargeCam::GetAverageArea() )
896//
897// Set weighted mean and variance of photo-electrons per sector in:
898// average sector pixels of MCalibrationChargeCam (obtained from:
899// MCalibrationChargeCam::GetAverageSector() )
900//
901Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
902{
903
904 const UInt_t npixels = fGeom->GetNumPixels();
905 const UInt_t nareas = fGeom->GetNumAreas();
906 const UInt_t nsectors = fGeom->GetNumSectors();
907
908 Float_t lowlim [nareas];
909 Float_t upplim [nareas];
910 Float_t areavars [nareas];
911 Float_t areaweights [nareas], sectorweights [nsectors];
912 Float_t areaphes [nareas], sectorphes [nsectors];
913 Int_t numareavalid[nareas], numsectorvalid[nsectors];
914
915 memset(lowlim ,0, nareas * sizeof(Float_t));
916 memset(upplim ,0, nareas * sizeof(Float_t));
917 memset(areaphes ,0, nareas * sizeof(Float_t));
918 memset(areavars ,0, nareas * sizeof(Float_t));
919 memset(areaweights ,0, nareas * sizeof(Float_t));
920 memset(numareavalid ,0, nareas * sizeof(Int_t ));
921 memset(sectorweights ,0, nsectors * sizeof(Float_t));
922 memset(sectorphes ,0, nsectors * sizeof(Float_t));
923 memset(numsectorvalid,0, nsectors * sizeof(Int_t ));
924
925 //
926 // First loop: Get mean number of photo-electrons and the RMS
927 // The loop is only to recognize later pixels with very deviating numbers
928 //
929 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
930
931 for (UInt_t i=0; i<npixels; i++)
932 {
933
934 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam) [i];
935 MBadPixelsPix &bad = (*fBadPixels)[i];
936
937 if (!pix.IsFFactorMethodValid())
938 continue;
939
940 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
941 {
942 pix.SetFFactorMethodValid(kFALSE);
943 continue;
944 }
945
946 // if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
947 // continue;
948
949 const Float_t nphe = pix.GetPheFFactorMethod();
950 const Float_t nvar = pix.GetPheFFactorMethod()*pix.GetPheFFactorMethod();
951 const Int_t aidx = (*fGeom)[i].GetAidx();
952
953 camphes.Fill(i,nphe);
954 camphes.SetUsed(i);
955
956 areaphes [aidx] += nphe;
957 areavars [aidx] += nvar;
958 numareavalid[aidx] ++;
959 }
960
961 for (UInt_t i=0; i<nareas; i++)
962 {
963 if (numareavalid[i] == 0)
964 {
965 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
966 << "in area index: " << i << endl;
967 continue;
968 }
969
970 areaphes[i] = areaphes[i] / numareavalid[i];
971 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1.);
972
973 if (areavars[i] > 0.)
974 areavars[i] = TMath::Sqrt(areavars[i]);
975 else
976 {
977 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of photo-electrons found "
978 << "in area index: " << i << endl;
979 continue;
980 }
981
982 lowlim [i] = areaphes[i] - fPheErrLimit*TMath::Sqrt(areavars[i]);
983 upplim [i] = areaphes[i] + fPheErrLimit*TMath::Sqrt(areavars[i]);
984
985 TArrayI area(1);
986 area[0] = i;
987
988 TH1D *hist = camphes.ProjectionS(TArrayI(),area,"_py",100);
989 hist->Fit("gaus","Q");
990 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
991 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
992 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
993
994 if (ndf < 2)
995 {
996 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
997 << "in the camera with area index: " << i << endl;
998 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
999 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1000 delete hist;
1001 continue;
1002 }
1003
1004 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1005
1006 if (prob < 0.001)
1007 {
1008 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1009 << "in the camera with area index: " << i << endl;
1010 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1011 << " is smaller than 0.001 " << endl;
1012 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1013 delete hist;
1014 continue;
1015 }
1016
1017 *fLog << inf << endl;
1018 *fLog << inf << GetDescriptor() << ": Mean Number of photo-electrons "
1019 << "in the camera with area index: " << i << ": "
1020 << Form("%4.2f%s%4.2f",mean,"+-",sigma) << endl;
1021 *fLog << inf << endl;
1022
1023 lowlim [i] = mean - fPheErrLimit*sigma;
1024 upplim [i] = mean + fPheErrLimit*sigma;
1025
1026 delete hist;
1027 }
1028
1029 memset(numareavalid,0,nareas*sizeof(Int_t));
1030 memset(areaphes ,0,nareas*sizeof(Int_t));
1031 memset(areavars ,0,nareas*sizeof(Int_t));
1032
1033 //
1034 // Second loop: Get mean number of photo-electrons and its RMS excluding
1035 // pixels deviating by more than fPheErrLimit sigma.
1036 // Set the conversion factor FADC counts to photo-electrons
1037 //
1038 for (UInt_t i=0; i<npixels; i++)
1039 {
1040
1041 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1042
1043 if (!pix.IsFFactorMethodValid())
1044 continue;
1045
1046 const Float_t nvar = pix.GetPheFFactorMethodVar();
1047
1048 if (nvar <= 0.)
1049 {
1050 pix.SetFFactorMethodValid(kFALSE);
1051 continue;
1052 }
1053
1054 MBadPixelsPix &bad = (*fBadPixels)[i];
1055
1056 const Int_t aidx = (*fGeom)[i].GetAidx();
1057 const Int_t sector = (*fGeom)[i].GetSector();
1058 const Float_t nphe = pix.GetPheFFactorMethod();
1059
1060 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1061 {
1062 *fLog << warn << GetDescriptor() << ": Deviating number of photo-electrons: "
1063 << Form("%4.2f",nphe) << " out of accepted limits: ["
1064 << Form("%4.2f%s%4.2f",lowlim[aidx],",",upplim[aidx]) << "] in pixel " << i << endl;
1065 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1066 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1067 pix.SetFFactorMethodValid(kFALSE);
1068 continue;
1069 }
1070
1071 // const Float_t weight = 1./nvar;
1072 // areaweights [aidx] += weight;
1073 // areaphes [aidx] += weight*nphe;
1074 areaweights [aidx] += nphe*nphe;
1075 areaphes [aidx] += nphe;
1076 numareavalid [aidx] ++;
1077
1078 if (aidx == 0)
1079 fNumInnerFFactorMethodUsed++;
1080 // sectorweights [sector] += weight;
1081 // sectorphes [sector] += weight*nphe;
1082 sectorweights [sector] += nphe*nphe;
1083 sectorphes [sector] += nphe;
1084 numsectorvalid[sector] ++;
1085 }
1086
1087 for (UInt_t aidx=0; aidx<nareas; aidx++)
1088 {
1089
1090 MCalibrationChargePix &apix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
1091
1092 areaweights[aidx] -= areaphes[aidx]*areaphes[aidx]/numareavalid[aidx];
1093 areaphes[aidx] /= numareavalid[aidx];
1094 areaweights[aidx] /= numareavalid[aidx]-1.;
1095
1096 if (areaweights[aidx] <= 0. || areaphes[aidx] <= 0.)
1097 {
1098 *fLog << warn << " Mean number of phe's from area index " << aidx << " cannot be calculated: "
1099 << " Mean number of phes: " << areaphes[aidx]
1100 << " Variance: " << areaweights[aidx] << endl;
1101 apix.SetFFactorMethodValid(kFALSE);
1102 continue;
1103 }
1104
1105 *fLog << inf << "Replacing number photo-electrons of average area idx " << aidx << ": "
1106 << Form("%5.3f%s%5.3f",apix.GetPheFFactorMethod()," +- ",apix.GetPheFFactorMethodErr()) << endl;
1107 *fLog << inf << " by average number of photo-electrons from area idx " << aidx << ": "
1108 << Form("%5.3f%s%5.3f",areaphes[aidx]," +- ",TMath::Sqrt(areaweights[aidx])) << endl;
1109
1110 // apix.SetPheFFactorMethod ( areaphes[aidx]/ areaweights[aidx] );
1111 // apix.SetPheFFactorMethodVar( 1. / areaweights[aidx] );
1112 apix.SetPheFFactorMethod ( areaphes[aidx] );
1113 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1114 apix.SetFFactorMethodValid ( kTRUE );
1115
1116 }
1117
1118 for (UInt_t sector=0; sector<nsectors; sector++)
1119 {
1120
1121 sectorweights[sector] -= sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector];
1122 sectorphes[sector] /= numsectorvalid[sector];
1123 sectorweights[sector] /= numsectorvalid[sector]-1.;
1124
1125 MCalibrationChargePix &spix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
1126
1127 if (sectorweights[sector] <= 0. || sectorphes[sector] <= 0.)
1128 {
1129 *fLog << warn << " Mean number of phe's from sector " << sector << " cannot be calculated: "
1130 << " Mean number of phes: " << sectorphes[sector]
1131 << " Variance: " << sectorweights[sector] << endl;
1132 spix.SetFFactorMethodValid(kFALSE);
1133 continue;
1134 }
1135
1136 *fLog << inf << "Replacing number photo-electrons of average sector " << sector << ": "
1137 << Form("%5.3f%s%5.3f",spix.GetPheFFactorMethod()," +- ",spix.GetPheFFactorMethodErr()) << endl;
1138 *fLog << inf << " by average number photo-electrons from sector " << sector << ": "
1139 << Form("%5.3f%s%5.3f",sectorphes[sector]," +- ",TMath::Sqrt(sectorweights[sector])) << endl;
1140
1141 // spix.SetPheFFactorMethod ( sectorphes[sector]/ sectorweights[sector] );
1142 // spix.SetPheFFactorMethodVar( 1. / sectorweights[sector] );
1143 spix.SetPheFFactorMethod ( sectorphes[sector] );
1144 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1145 spix.SetFFactorMethodValid ( kTRUE );
1146
1147 }
1148
1149 return kTRUE;
1150}
1151
1152
1153// ------------------------------------------------------------------------
1154//
1155// Returns kFALSE if pointer to MCalibrationChargeBlindPix is NULL
1156//
1157// The check returns kFALSE if:
1158//
1159// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1160// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1161//
1162// Calls:
1163// - MCalibrationChargeBlindPix::CalcFluxInsidePlexiglass()
1164//
1165Bool_t MCalibrationChargeCalc::FinalizeBlindPixel()
1166{
1167
1168 if (!fBlindPixel)
1169 return kFALSE;
1170
1171 const Float_t lambda = fBlindPixel->GetLambda();
1172 const Float_t lambdaerr = fBlindPixel->GetLambdaErr();
1173 const Float_t lambdacheck = fBlindPixel->GetLambdaCheck();
1174
1175 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1176 {
1177 *fLog << warn << GetDescriptor()
1178 << Form("%s%4.2f%s%4.2f%s%4.2f%s",": Lambda: ",lambda," and Lambda-Check: ",
1179 lambdacheck," differ by more than ",fLambdaCheckLimit," in the Blind Pixel ")
1180 << endl;
1181 return kFALSE;
1182 }
1183
1184 if (lambdaerr > fLambdaErrLimit)
1185 {
1186 *fLog << warn << GetDescriptor()
1187 << Form("%s%4.2f%s%4.2f%s",": Error of Fitted Lambda: ",lambdaerr," is greater than ",
1188 fLambdaErrLimit," in Blind Pixel ") << endl;
1189 return kFALSE;
1190 }
1191
1192 if (!fBlindPixel->CalcFluxInsidePlexiglass())
1193 {
1194 *fLog << warn << "Could not calculate the flux of photons from the Blind Pixel, "
1195 << "will skip Blind Pixel Calibration " << endl;
1196 return kFALSE;
1197 }
1198
1199 return kTRUE;
1200}
1201
1202// ------------------------------------------------------------------------
1203//
1204// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1205//
1206// The check returns kFALSE if:
1207//
1208// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1209// 2) PINDiode has a fit error smaller than fChargeErrLimit
1210// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1211// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1212//
1213// Calls:
1214// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1215//
1216Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1217{
1218
1219 if (!fPINDiode)
1220 return kFALSE;
1221
1222 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1223 {
1224 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1225 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1226 return kFALSE;
1227 }
1228
1229 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1230 {
1231 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1232 << fChargeErrLimit << " in PINDiode " << endl;
1233 return kFALSE;
1234 }
1235
1236 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1237 {
1238 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1239 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1240 return kFALSE;
1241 }
1242
1243 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1244 {
1245 *fLog << warn << GetDescriptor()
1246 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1247 return kFALSE;
1248 }
1249
1250
1251 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1252 {
1253 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1254 << "will skip PIN Diode Calibration " << endl;
1255 return kFALSE;
1256 }
1257
1258 return kTRUE;
1259}
1260
1261// ------------------------------------------------------------------------
1262//
1263// Calculate the average number of photons outside the plexiglass with the
1264// formula:
1265//
1266// av.Num.photons(area index) = av.Num.Phes(area index)
1267// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1268// / MCalibrationQEPix::GetPMTCollectionEff()
1269// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1270// / MCalibrationQECam::GetPlexiglassQE()
1271//
1272// Calculate the variance on the average number of photons assuming that the error on the
1273// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1274// values keeps it ordinary error since it is systematic.
1275//
1276// Loop over pixels:
1277//
1278// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1279// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1280//
1281// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1282// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1283//
1284// - Calculate the quantum efficiency with the formula:
1285//
1286// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1287//
1288// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1289//
1290// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1291// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1292//
1293// - Call MCalibrationQEPix::UpdateFFactorMethod()
1294//
1295void MCalibrationChargeCalc::FinalizeFFactorQECam()
1296{
1297
1298 if (fNumInnerFFactorMethodUsed < 2)
1299 {
1300 *fLog << warn << GetDescriptor()
1301 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1302 return;
1303 }
1304
1305 MCalibrationChargePix &avpix = (MCalibrationChargePix&)fCam->GetAverageArea(0);
1306 MCalibrationQEPix &qepix = (MCalibrationQEPix&) fQECam->GetAverageArea(0);
1307
1308 const Float_t avphotons = avpix.GetPheFFactorMethod()
1309 / qepix.GetDefaultQE(fPulserColor)
1310 / qepix.GetPMTCollectionEff()
1311 / qepix.GetLightGuidesEff(fPulserColor)
1312 / fQECam->GetPlexiglassQE();
1313
1314 const Float_t avphotrelvar = avpix.GetPheFFactorMethodRelVar()
1315 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1316 + qepix.GetPMTCollectionEffRelVar()
1317 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1318 + fQECam->GetPlexiglassQERelVar();
1319
1320 const UInt_t nareas = fGeom->GetNumAreas();
1321
1322 //
1323 // Set the results in the MCalibrationChargeCam
1324 //
1325 fCam->SetNumPhotonsFFactorMethod (avphotons);
1326 if (avphotrelvar > 0.)
1327 fCam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1328
1329 Float_t lowlim [nareas];
1330 Float_t upplim [nareas];
1331 Float_t avffactorphotons [nareas];
1332 Float_t avffactorphotvar [nareas];
1333 Int_t numffactor [nareas];
1334
1335 memset(lowlim ,0, nareas * sizeof(Float_t));
1336 memset(upplim ,0, nareas * sizeof(Float_t));
1337 memset(avffactorphotons,0, nareas * sizeof(Float_t));
1338 memset(avffactorphotvar,0, nareas * sizeof(Float_t));
1339 memset(numffactor ,0, nareas * sizeof(Int_t));
1340
1341 const UInt_t npixels = fGeom->GetNumPixels();
1342
1343 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1344
1345 for (UInt_t i=0; i<npixels; i++)
1346 {
1347
1348 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1349 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1350 MBadPixelsPix &bad = (*fBadPixels)[i];
1351
1352 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1353 continue;
1354
1355 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1356 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1357
1358 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1359 {
1360 (*fBadPixels)[i].SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1361 continue;
1362 }
1363
1364 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1365
1366 qepix.SetQEFFactor ( qe , fPulserColor );
1367 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1368 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1369
1370 if (!qepix.UpdateFFactorMethod())
1371 *fLog << warn << GetDescriptor()
1372 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1373
1374 const Int_t aidx = (*fGeom)[i].GetAidx();
1375 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1376
1377 camffactor.Fill(i,ffactor);
1378 camffactor.SetUsed(i);
1379
1380 avffactorphotons[aidx] += ffactor;
1381 avffactorphotvar[aidx] += pix.GetMeanFFactorFADC2PhotVar();
1382 numffactor[aidx]++;
1383 }
1384
1385 for (UInt_t i=0; i<nareas; i++)
1386 {
1387 avffactorphotons[i] /= numffactor[i];
1388 avffactorphotvar[i] /= numffactor[i];
1389
1390 lowlim [i] = 1.1; // Lowest known F-Factor of a PMT
1391 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1392
1393 TArrayI area(1);
1394 area[0] = i;
1395
1396 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1397 hist->Fit("gaus","Q");
1398 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1399 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1400 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1401
1402 if (ndf < 2)
1403 {
1404 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1405 << "in the camera with area index: " << i << endl;
1406 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1407 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1408 delete hist;
1409 continue;
1410 }
1411
1412 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1413
1414 if (prob < 0.001)
1415 {
1416 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1417 << "in the camera with area index: " << i << endl;
1418 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1419 << " is smaller than 0.001 " << endl;
1420 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1421 delete hist;
1422 continue;
1423 }
1424
1425 *fLog << inf << endl;
1426 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1427 << "in the camera with area index: " << i << ": "
1428 << Form("%4.2f%s%4.2f",mean,"+-",sigma) << endl;
1429 *fLog << inf << endl;
1430
1431 lowlim [i] = mean - fFFactorErrLimit*sigma;
1432 upplim [i] = mean + fFFactorErrLimit*sigma;
1433
1434 delete hist;
1435 }
1436
1437 for (UInt_t i=0; i<npixels; i++)
1438 {
1439
1440 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1441 MBadPixelsPix &bad = (*fBadPixels)[i];
1442
1443 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1444 continue;
1445
1446 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1447 const Int_t aidx = (*fGeom)[i].GetAidx();
1448
1449 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1450 {
1451 *fLog << warn << GetDescriptor() << ": Deviating F-Factor: "
1452 << Form("%4.2f",ffactor) << " out of accepted limits: ["
1453 << Form("%4.2f%s%4.2f",lowlim[aidx],",",upplim[aidx]) << "] in pixel " << i << endl;
1454 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1455 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1456 }
1457 }
1458
1459 for (UInt_t i=0; i<npixels; i++)
1460 {
1461
1462 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1463 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1464 MBadPixelsPix &bad = (*fBadPixels)[i];
1465
1466 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1467 {
1468 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1469 pix.SetFFactorMethodValid(kFALSE);
1470 pix.SetExcluded();
1471 continue;
1472 }
1473 }
1474}
1475
1476
1477// ------------------------------------------------------------------------
1478//
1479// Loop over pixels:
1480//
1481// - Continue, if not MCalibrationChargeBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1482// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1483//
1484// - Calculate the quantum efficiency with the formula:
1485//
1486// QE = Num.Phes / MCalibrationChargeBlindPix::GetFluxInsidePlexiglass()
1487// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1488//
1489// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1490// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1491// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1492//
1493// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1494//
1495void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1496{
1497
1498 const UInt_t npixels = fGeom->GetNumPixels();
1499
1500 //
1501 // Set the results in the MCalibrationChargeCam
1502 //
1503 if (fBlindPixel)
1504 {
1505 if (fBlindPixel->IsFluxInsidePlexiglassAvailable())
1506 {
1507
1508 const Float_t photons = fBlindPixel->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
1509 / fQECam->GetPlexiglassQE();
1510 fCam->SetNumPhotonsBlindPixelMethod(photons);
1511
1512 const Float_t photrelvar = fBlindPixel->GetFluxInsidePlexiglassRelVar()
1513 + fQECam->GetPlexiglassQERelVar();
1514 if (photrelvar > 0.)
1515 fCam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
1516 }
1517 }
1518 //
1519 // With the knowledge of the overall photon flux, calculate the
1520 // quantum efficiencies after the Blind Pixel and PIN Diode method
1521 //
1522 for (UInt_t i=0; i<npixels; i++)
1523 {
1524
1525 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1526
1527 if (!fBlindPixel)
1528 {
1529 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1530 continue;
1531 }
1532
1533 if (!fBlindPixel->IsFluxInsidePlexiglassAvailable())
1534 {
1535 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1536 continue;
1537 }
1538
1539 MBadPixelsPix &bad = (*fBadPixels)[i];
1540 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1541 {
1542 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1543 continue;
1544 }
1545
1546 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1547 MGeomPix &geo = (*fGeom)[i];
1548
1549 const Float_t qe = pix.GetPheFFactorMethod()
1550 / fBlindPixel->GetFluxInsidePlexiglass()
1551 / geo.GetA()
1552 * fQECam->GetPlexiglassQE();
1553
1554 const Float_t qerelvar = fBlindPixel->GetFluxInsidePlexiglassRelVar()
1555 + fQECam->GetPlexiglassQERelVar()
1556 + pix.GetPheFFactorMethodRelVar();
1557
1558 qepix.SetQEBlindPixel ( qe , fPulserColor );
1559 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
1560 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
1561
1562 if (!qepix.UpdateBlindPixelMethod())
1563 *fLog << warn << GetDescriptor()
1564 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
1565 }
1566}
1567
1568// ------------------------------------------------------------------------
1569//
1570// Loop over pixels:
1571//
1572// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
1573// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
1574//
1575// - Calculate the quantum efficiency with the formula:
1576//
1577// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
1578//
1579// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
1580// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
1581// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
1582//
1583// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
1584//
1585void MCalibrationChargeCalc::FinalizePINDiodeQECam()
1586{
1587
1588 const UInt_t npixels = fGeom->GetNumPixels();
1589
1590 //
1591 // With the knowledge of the overall photon flux, calculate the
1592 // quantum efficiencies after the PIN Diode method
1593 //
1594 for (UInt_t i=0; i<npixels; i++)
1595 {
1596
1597 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1598
1599 if (!fPINDiode)
1600 {
1601 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1602 continue;
1603 }
1604
1605 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
1606 {
1607 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1608 continue;
1609 }
1610
1611 MBadPixelsPix &bad = (*fBadPixels)[i];
1612
1613 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1614 {
1615 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1616 continue;
1617 }
1618
1619 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1620 MGeomPix &geo = (*fGeom)[i];
1621
1622 const Float_t qe = pix.GetPheFFactorMethod()
1623 / fPINDiode->GetFluxOutsidePlexiglass()
1624 / geo.GetA();
1625
1626 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
1627
1628 qepix.SetQEPINDiode ( qe , fPulserColor );
1629 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
1630 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
1631
1632 if (!qepix.UpdatePINDiodeMethod())
1633 *fLog << warn << GetDescriptor()
1634 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
1635 }
1636}
1637
1638// -----------------------------------------------------------------------------------------------
1639//
1640// - Print out statistics about BadPixels of type UnsuitableType_t
1641// - store numbers of bad pixels of each type in fCam
1642//
1643void MCalibrationChargeCalc::FinalizeUnsuitablePixels()
1644{
1645
1646 *fLog << inf << endl;
1647 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
1648 *fLog << dec << setfill(' ');
1649
1650 const Int_t nareas = fGeom->GetNumAreas();
1651
1652 Int_t counts[nareas];
1653 memset(counts,0,nareas*sizeof(Int_t));
1654
1655 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1656 {
1657 MBadPixelsPix &bad = (*fBadPixels)[i];
1658 if (!bad.IsBad())
1659 {
1660 const Int_t aidx = (*fGeom)[i].GetAidx();
1661 counts[aidx]++;
1662 }
1663 }
1664
1665 if (fGeom->InheritsFrom("MGeomCamMagic"))
1666 *fLog << " " << setw(7) << "Successfully calibrated Pixels: "
1667 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1668
1669 memset(counts,0,nareas*sizeof(Int_t));
1670
1671 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1672 {
1673 MBadPixelsPix &bad = (*fBadPixels)[i];
1674 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1675 {
1676 const Int_t aidx = (*fGeom)[i].GetAidx();
1677 counts[aidx]++;
1678 }
1679 }
1680
1681 for (Int_t aidx=0; aidx<nareas; aidx++)
1682 fCam->SetNumUnsuitable(counts[aidx], aidx);
1683
1684 if (fGeom->InheritsFrom("MGeomCamMagic"))
1685 *fLog << " " << setw(7) << "Uncalibrated Pixels: "
1686 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1687
1688 memset(counts,0,nareas*sizeof(Int_t));
1689
1690 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1691 {
1692 MBadPixelsPix &bad = (*fBadPixels)[i];
1693 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1694 {
1695 const Int_t aidx = (*fGeom)[i].GetAidx();
1696 counts[aidx]++;
1697 }
1698 }
1699
1700 for (Int_t aidx=0; aidx<nareas; aidx++)
1701 fCam->SetNumUnreliable(counts[aidx], aidx);
1702
1703 *fLog << " " << setw(7) << "Unreliable Pixels: "
1704 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1705
1706}
1707
1708// -----------------------------------------------------------------------------------------------
1709//
1710// Print out statistics about BadPixels of type UncalibratedType_t
1711//
1712void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
1713{
1714
1715 UInt_t countinner = 0;
1716 UInt_t countouter = 0;
1717
1718 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1719 {
1720 MBadPixelsPix &bad = (*fBadPixels)[i];
1721 if (bad.IsUncalibrated(typ))
1722 {
1723 if (fGeom->GetPixRatio(i) == 1.)
1724 countinner++;
1725 else
1726 countouter++;
1727 }
1728 }
1729
1730 *fLog << " " << setw(7) << text
1731 << Form("%s%3i%s%3i","Inner: ",countinner," Outer: ",countouter) << endl;
1732}
1733
1734// --------------------------------------------------------------------------
1735//
1736// Set the path for output file
1737//
1738void MCalibrationChargeCalc::SetOutputPath(TString path)
1739{
1740 fOutputPath = path;
1741 if (fOutputPath.EndsWith("/"))
1742 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
1743}
1744
1745void MCalibrationChargeCalc::SetOutputFile(TString file)
1746{
1747 fOutputFile = file;
1748}
1749
1750// --------------------------------------------------------------------------
1751//
1752// Get the output file
1753//
1754const char* MCalibrationChargeCalc::GetOutputFile()
1755{
1756 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
1757}
Note: See TracBrowser for help on using the repository browser.