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

Last change on this file since 4379 was 4360, checked in by gaug, 21 years ago
*** empty log message ***
File size: 60.9 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MCalibrationChargeCalc
28//
29// Task to calculate the calibration conversion factors and quantum efficiencies
30// from the fit results to the summed FADC slice distributions delivered by
31// MCalibrationChargeCam, MCalibrationChargePix, MCalibrationChargeBlindPix and
32// MCalibrationChargePINDiode, calculated and filled by MHCalibrationChargeCam,
33// MHCalibrationChargePix, MHCalibrationChargeBlindPix and MHCalibrationChargePINDiode.
34//
35// PreProcess(): Initialize pointers to MCalibrationChargeCam, MCalibrationChargeBlindPix
36// MCalibrationChargePINDiode and MCalibrationQECam
37//
38// Initialize pulser light wavelength
39//
40// ReInit(): MCalibrationCam::InitSize(NumPixels) is called from MGeomApply (which allocates
41// memory in a TClonesArray of type MCalibrationChargePix)
42// Initializes pointer to MBadPixelsCam
43//
44// Process(): Nothing to be done, histograms getting filled by MHCalibrationChargeCam
45//
46// PostProcess(): - FinalizePedestals()
47// - FinalizeCharges()
48// - FinalizeFFactorMethod()
49// - FinalizeBadPixels()
50// - FinalizeBlindPixel()
51// - 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 SetDebug ( kFALSE );
285
286 Clear();
287
288}
289
290// --------------------------------------------------------------------------
291//
292// Sets:
293// - all variables to 0.,
294// - all flags to kFALSE
295//
296void MCalibrationChargeCalc::Clear(const Option_t *o)
297{
298
299 fNumHiGainSamples = 0.;
300 fNumLoGainSamples = 0.;
301 fSqrtHiGainSamples = 0.;
302 fSqrtLoGainSamples = 0.;
303 fNumInnerFFactorMethodUsed = 0;
304 SkipHiLoGainCalibration( kFALSE );
305}
306
307
308// -----------------------------------------------------------------------------------
309//
310// The following container are searched for and execution aborted if not in MParList:
311// - MPedestalCam
312//
313// The following containers are searched and created if they were not found:
314//
315// - MCalibrationQECam
316// - MBadPixelsCam
317//
318// The following output containers are only searched, but not created. If they
319// cannot be found, the corresponding calibration part is only skipped.
320//
321// - MTime
322//
323Int_t MCalibrationChargeCalc::PreProcess(MParList *pList)
324{
325
326 //
327 // Containers that have to be there.
328 //
329 fPedestals = (MPedestalCam*)pList->FindObject("MPedestalCam");
330 if (!fPedestals)
331 {
332 *fLog << err << "MPedestalCam not found... aborting" << endl;
333 return kFALSE;
334 }
335
336 //
337 // Containers that are created in case that they are not there.
338 //
339 fQECam = (MCalibrationQECam*)pList->FindCreateObj("MCalibrationQECam");
340 if (!fQECam)
341 {
342 *fLog << err << "Cannot find nor create MCalibrationQECam... aborting" << endl;
343 return kFALSE;
344 }
345
346 fBadPixels = (MBadPixelsCam*)pList->FindCreateObj("MBadPixelsCam");
347 if (!fBadPixels)
348 {
349 *fLog << err << "Could not find or create MBadPixelsCam ... aborting." << endl;
350 return kFALSE;
351 }
352
353
354 fEvtTime = (MTime*)pList->FindObject("MTime");
355
356 //
357 // Check the pulser colour --> FIXME: this solution is only valid until the arrival of the DM's
358 //
359 if (fPulserColor == MCalibrationCam::kNONE)
360 {
361 *fLog << endl;
362 *fLog << err << GetDescriptor()
363 << ": No Pulser colour has been chosen. Since the installation of the IFAE pulser box,"
364 << " you HAVE to provide the LEDs colour, otherwise there is no calibration. " << endl;
365 *fLog << "See e.g. the macro calibration.C " << endl;
366 return kFALSE;
367 }
368
369 return kTRUE;
370}
371
372
373// --------------------------------------------------------------------------
374//
375// Search for the following input containers and abort if not existing:
376// - MGeomCam
377// - MCalibrationChargeCam
378//
379// Search for the following input containers and give a warning if not existing:
380// - MCalibrationChargeBlindPix
381// - MCalibrationChargePINDiode
382//
383// It retrieves the following variables from MCalibrationChargeCam:
384//
385// - fNumHiGainSamples
386// - fNumLoGainSamples
387//
388// It defines the PixId of every pixel in:
389//
390// - MCalibrationChargeCam
391// - MCalibrationQECam
392//
393// It sets all pixels in excluded which have the flag fBadBixelsPix::IsBad() set in:
394//
395// - MCalibrationChargePix
396// - MCalibrationQEPix
397//
398// Sets the pulser colour and tests if it has not changed w.r.t. fPulserColor in:
399//
400// - MCalibrationChargeCam
401// - MCalibrationChargeBlindPix (if existing)
402// - MCalibrationChargePINDiode (if existing)
403//
404Bool_t MCalibrationChargeCalc::ReInit(MParList *pList )
405{
406
407 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
408 if (!fGeom)
409 {
410 *fLog << err << "No MGeomCam found... aborting." << endl;
411 return kFALSE;
412 }
413
414 fCam = (MCalibrationChargeCam*)pList->FindObject("MCalibrationChargeCam");
415 if (!fCam)
416 {
417 *fLog << err << "Cannot find MCalibrationChargeCam... aborting" << endl;
418 *fLog << err << "Maybe you forget to call an MFillH for the MHCalibrationChargeCam before..." << endl;
419 return kFALSE;
420 }
421
422 //
423 // Optional Containers
424 //
425 fBlindPixel = (MCalibrationChargeBlindPix*)pList->FindObject("MCalibrationChargeBlindPix");
426 if (!fBlindPixel)
427 {
428 *fLog << endl;
429 *fLog << warn << GetDescriptor()
430 << ": MCalibrationChargeBlindPix not found... no Blind Pixel method! " << endl;
431 }
432
433 fPINDiode = (MCalibrationChargePINDiode*)pList->FindObject("MCalibrationChargePINDiode");
434 if (!fPINDiode)
435 {
436 *fLog << endl;
437 *fLog << warn << GetDescriptor()
438 << ": MCalibrationChargePINDiode not found... no PIN Diode method! " << endl;
439 }
440
441
442 //
443 // Initialize the pulser colours
444 //
445 if (fCam->GetPulserColor() == MCalibrationCam::kNONE)
446 {
447 fCam->SetPulserColor( fPulserColor );
448
449 if (fBlindPixel)
450 fBlindPixel->SetColor( fPulserColor );
451
452 if (fPINDiode)
453 fPINDiode->SetColor( fPulserColor );
454 }
455
456 if (fPulserColor != fCam->GetPulserColor())
457 {
458 *fLog << err << GetDescriptor()
459 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargeCam" << endl;
460 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
461 return kFALSE;
462 }
463
464 if (fBlindPixel)
465 if (fPulserColor != fBlindPixel->GetColor())
466 {
467 *fLog << err << GetDescriptor()
468 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargeBlindPix." << endl;
469 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
470 return kFALSE;
471 }
472
473 if (fPINDiode)
474 if (fPulserColor != fPINDiode->GetColor())
475 {
476 *fLog << err << GetDescriptor()
477 << ": Pulser colour has changed w.r.t. last file in MCalibrationChargePINDiode." << endl;
478 *fLog << err << "This feature is not yet implemented, sorry ... aborting " << endl;
479 return kFALSE;
480 }
481
482
483 fNumHiGainSamples = fCam->GetNumHiGainFADCSlices();
484 fNumLoGainSamples = fCam->GetNumLoGainFADCSlices();
485
486 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
487 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
488
489 UInt_t npixels = fGeom->GetNumPixels();
490
491 for (UInt_t i=0; i<npixels; i++)
492 {
493
494 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam) [i];
495 MCalibrationQEPix &pqe = (MCalibrationQEPix&) (*fQECam)[i];
496 MBadPixelsPix &bad = (*fBadPixels)[i];
497
498 pix.SetPixId(i);
499 pqe.SetPixId(i);
500
501 if (bad.IsBad())
502 {
503 pix.SetExcluded();
504 pqe.SetExcluded();
505 continue;
506 }
507
508 if (IsDebug())
509 pix.SetDebug();
510 }
511
512 return kTRUE;
513}
514
515// ----------------------------------------------------------------------------------
516//
517// Nothing to be done in Process, but have a look at MHCalibrationChargeCam, instead
518//
519Int_t MCalibrationChargeCalc::Process()
520{
521 return kTRUE;
522}
523
524// -----------------------------------------------------------------------
525//
526// Return if number of executions is null.
527//
528// First loop over pixels, average areas and sectors, call:
529// - FinalizePedestals()
530// - FinalizeCharges()
531// for every entry. Count number of valid pixels in loop and return kFALSE
532// if there are none (the "Michele check").
533//
534// Call FinalizeBadPixels()
535//
536// Call FinalizeFFactorMethod() (second and third loop over pixels and areas)
537//
538// Call FinalizeBlindPixel()
539// Call FinalizePINDiode()
540//
541// Call FinalizeFFactorQECam() (fourth loop over pixels and areas)
542// Call FinalizeBlindPixelQECam() (fifth loop over pixels and areas)
543// Call FinalizePINDiodeQECam() (sixth loop over pixels and areas)
544//
545// Call FinalizeUnsuitablePixels()
546//
547// Call MParContainer::SetReadyToSave() for fCam, fQECam, fBadPixels and
548// fBlindPixel and fPINDiode if they exist
549//
550// Print out some statistics
551//
552Int_t MCalibrationChargeCalc::PostProcess()
553{
554
555 if (GetNumExecutions()==0)
556 return kFALSE;
557
558 if (fPINDiode)
559 if (!fPINDiode->IsValid())
560 {
561 *fLog << warn << GetDescriptor()
562 << ": MCalibrationChargePINDiode is declared not valid... no PIN Diode method! " << endl;
563 fPINDiode = NULL;
564 }
565
566 if (fBlindPixel)
567 if (!fBlindPixel->IsValid())
568 {
569 *fLog << warn << GetDescriptor()
570 << ": MCalibrationChargeBlindPix is declared not valid... no Blind Pixel method! " << endl;
571 fBlindPixel = NULL;
572 }
573
574 //
575 // First loop over pixels, call FinalizePedestals and FinalizeCharges
576 //
577 Int_t nvalid = 0;
578
579 for (Int_t pixid=0; pixid<fPedestals->GetSize(); pixid++)
580 {
581
582 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[pixid];
583 //
584 // Check if the pixel has been excluded from the fits
585 //
586 if (pix.IsExcluded())
587 continue;
588
589 MPedestalPix &ped = (*fPedestals)[pixid];
590 MBadPixelsPix &bad = (*fBadPixels)[pixid];
591 const Int_t aidx = (*fGeom)[pixid].GetAidx();
592
593 FinalizePedestals(ped,pix,aidx);
594
595 if (FinalizeCharges(pix,bad))
596 nvalid++;
597 }
598
599 //
600 // The Michele check ...
601 //
602 if (nvalid == 0)
603 {
604 *fLog << err << GetDescriptor() << ": All pixels have non-valid calibration. "
605 << "Did you forget to fill the histograms "
606 << "(filling MHCalibrationChargeCam from MExtractedSignalCam using MFillH) ? " << endl;
607 *fLog << err << GetDescriptor() << ": Or, maybe, you have used a pedestal run "
608 << "instead of a calibration run " << endl;
609 return kFALSE;
610 }
611
612 for (UInt_t aidx=0; aidx<fGeom->GetNumAreas(); aidx++)
613 {
614
615 const MPedestalPix &ped = fPedestals->GetAverageArea(aidx);
616 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
617
618 FinalizePedestals(ped,pix,aidx);
619 FinalizeCharges(pix, fCam->GetAverageBadArea(aidx));
620 }
621
622 for (UInt_t sector=0; sector<fGeom->GetNumSectors(); sector++)
623 {
624
625 const MPedestalPix &ped = fPedestals->GetAverageSector(sector);
626
627 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
628 FinalizePedestals(ped,pix, 0);
629 FinalizeCharges(pix, fCam->GetAverageBadSector(sector));
630 }
631
632 //
633 // Finalize Bad Pixels
634 //
635 FinalizeBadPixels();
636
637 //
638 // Finalize F-Factor method
639 //
640 if (!FinalizeFFactorMethod())
641 {
642 *fLog << warn << "Could not calculate the photons flux from the F-Factor method " << endl;
643 fCam->SetFFactorMethodValid(kFALSE);
644 return kFALSE;
645 }
646 else
647 fCam->SetFFactorMethodValid(kTRUE);
648
649 //
650 // Finalize Blind Pixel
651 //
652 if (FinalizeBlindPixel())
653 fQECam->SetBlindPixelMethodValid(kTRUE);
654 else
655 fQECam->SetBlindPixelMethodValid(kFALSE);
656
657 //
658 // Finalize PIN Diode
659 //
660 if (FinalizePINDiode())
661 fQECam->SetPINDiodeMethodValid(kTRUE);
662 else
663 fQECam->SetPINDiodeMethodValid(kFALSE);
664
665 //
666 // Finalize QE Cam
667 //
668 FinalizeFFactorQECam();
669 FinalizeBlindPixelQECam();
670 FinalizePINDiodeQECam();
671
672 //
673 // Re-direct the output to an ascii-file from now on:
674 //
675 MLog asciilog;
676 asciilog.SetOutputFile(GetOutputFile(),kTRUE);
677 SetLogStream(&asciilog);
678 //
679 // Finalize calibration statistics
680 //
681 FinalizeUnsuitablePixels();
682
683 fCam ->SetReadyToSave();
684 fQECam ->SetReadyToSave();
685 fBadPixels->SetReadyToSave();
686
687 if (fBlindPixel)
688 fBlindPixel->SetReadyToSave();
689 if (fPINDiode)
690 fPINDiode->SetReadyToSave();
691
692 *fLog << inf << endl;
693 *fLog << GetDescriptor() << ": Fatal errors statistics:" << endl;
694
695 PrintUncalibrated(MBadPixelsPix::kChargeIsPedestal,
696 Form("%s%2.1f%s","Signal less than ",fChargeLimit," Pedestal RMS: "));
697 PrintUncalibrated(MBadPixelsPix::kChargeRelErrNotValid,
698 Form("%s%2.1f%s","Signal Error bigger than ",fChargeRelErrLimit," times Mean Signal: "));
699 PrintUncalibrated(MBadPixelsPix::kChargeSigmaNotValid,
700 "Signal Sigma smaller than Pedestal RMS: ");
701 PrintUncalibrated(MBadPixelsPix::kLoGainSaturation,
702 "Pixels with Low Gain Saturation: ");
703 PrintUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin,
704 Form("%s%2.1f%s","Mean Abs. Arr. Time in First ",1.," Bin(s): "));
705 PrintUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins,
706 Form("%s%2.1f%s","Mean Abs. Arr. Time in Last ",2.," Bin(s): "));
707 PrintUncalibrated(MBadPixelsPix::kDeviatingNumPhes,
708 "Pixels with deviating number of phes: ");
709 PrintUncalibrated(MBadPixelsPix::kDeviatingFFactor,
710 "Pixels with deviating F-Factor: ");
711
712 *fLog << inf << endl;
713 *fLog << GetDescriptor() << ": Unreliable errors statistics:" << endl;
714
715 PrintUncalibrated(MBadPixelsPix::kHiGainOscillating,
716 "Pixels with changing Hi Gain signal over time: ");
717 PrintUncalibrated(MBadPixelsPix::kLoGainOscillating,
718 "Pixels with changing Lo Gain signal over time: ");
719 PrintUncalibrated(MBadPixelsPix::kHiGainNotFitted,
720 "Pixels with unsuccesful Gauss fit to the Hi Gain: ");
721 PrintUncalibrated(MBadPixelsPix::kLoGainNotFitted,
722 "Pixels with unsuccesful Gauss fit to the Lo Gain: ");
723
724 SetLogStream(&gLog);
725
726 return kTRUE;
727}
728
729// ----------------------------------------------------------------------------------
730//
731// Retrieves pedestal and pedestal RMS from MPedestalPix
732// Retrieves total entries from MPedestalCam
733// Sets pedestal*fNumHiGainSamples and pedestal*fNumLoGainSamples in MCalibrationChargePix
734// Sets pedRMS *fSqrtHiGainSamples and pedRMS *fSqrtLoGainSamples in MCalibrationChargePix
735//
736// If the flag MCalibrationPix::IsHiGainSaturation() is set, call also:
737// - MCalibrationChargePix::CalcLoGainPedestal()
738//
739void MCalibrationChargeCalc::FinalizePedestals(const MPedestalPix &ped, MCalibrationChargePix &cal, const Int_t aidx)
740{
741
742 //
743 // get the pedestals
744 //
745 const Float_t pedes = ped.GetPedestal();
746 const Float_t prms = ped.GetPedestalRms();
747 const Float_t num = TMath::Sqrt((Float_t)fPedestals->GetTotalEntries());
748
749
750 //
751 // set them in the calibration camera
752 //
753 if (cal.IsHiGainSaturation())
754 {
755 cal.SetPedestal(pedes* fNumLoGainSamples,
756 prms * fSqrtLoGainSamples,
757 prms * fNumLoGainSamples / num);
758 cal.CalcLoGainPedestal((Float_t)fNumLoGainSamples, aidx);
759 }
760 else
761 {
762 cal.SetPedestal(pedes* fNumHiGainSamples,
763 prms * fSqrtHiGainSamples,
764 prms * fNumHiGainSamples / num);
765 }
766
767}
768
769// ----------------------------------------------------------------------------------------------------
770//
771// Check fit results validity. Bad Pixels flags are set if:
772//
773// 1) Pixel has a mean smaller than fChargeLimit*PedRMS ( Flag: MBadPixelsPix::kChargeIsPedestal)
774// 2) Pixel has a mean error smaller than fChargeErrLimit ( Flag: MBadPixelsPix::kChargeErrNotValid)
775// 3) Pixel has mean smaller than fChargeRelVarLimit times its mean error
776// ( Flag: MBadPixelsPix::kChargeRelErrNotValid)
777// 4) Pixel has a sigma bigger than its Pedestal RMS ( Flag: MBadPixelsPix::kChargeSigmaNotValid )
778//
779// Further returns if flags: MBadPixelsPix::kUnsuitableRun is set
780//
781// Calls MCalibrationChargePix::CalcReducedSigma() and sets flag: MBadPixelsPix::kChargeIsPedestal
782// if not succesful.
783//
784// Calls MCalibrationChargePix::CalcFFactorMethod() and sets flag: MBadPixelsPix::kDeviatingNumPhes)
785// if not succesful.
786//
787Bool_t MCalibrationChargeCalc::FinalizeCharges(MCalibrationChargePix &cal, MBadPixelsPix &bad)
788{
789
790 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
791 return kFALSE;
792
793 if (cal.GetMean() < fChargeLimit*cal.GetPedRms())
794 {
795 *fLog << warn << GetDescriptor() << ": Fitted Charge: " << cal.GetMean() << " is smaller than "
796 << fChargeLimit << " Pedestal RMS: " << cal.GetPedRms() << " in Pixel " << cal.GetPixId() << endl;
797 bad.SetUncalibrated( MBadPixelsPix::kChargeIsPedestal);
798 }
799
800 if (cal.GetMean() < fChargeRelErrLimit*cal.GetMeanErr())
801 {
802 *fLog << warn << GetDescriptor() << ": Fitted Charge: " << cal.GetMean() << " is smaller than "
803 << fChargeRelErrLimit << "* its error: " << cal.GetMeanErr()
804 << " in Pixel " << cal.GetPixId() << endl;
805 bad.SetUncalibrated( MBadPixelsPix::kChargeRelErrNotValid );
806 }
807
808 if (cal.GetSigma() < cal.GetPedRms())
809 {
810 *fLog << warn << GetDescriptor() << ": Sigma of Fitted Charge: " << cal.GetSigma()
811 << " smaller than Pedestal RMS: " << cal.GetPedRms() << " in Pixel " << cal.GetPixId() << endl;
812 bad.SetUncalibrated( MBadPixelsPix::kChargeSigmaNotValid );
813 }
814
815 if (!cal.CalcReducedSigma())
816 {
817 *fLog << warn << GetDescriptor()
818 << ": Could not calculate reduced sigmas of pixel: " << cal.GetPixId() << endl;
819 bad.SetUncalibrated(MBadPixelsPix::kChargeIsPedestal);
820 return kFALSE;
821 }
822
823 if (!cal.CalcFFactorMethod())
824 {
825 *fLog << warn << GetDescriptor()
826 << ": Could not calculate F-Factor of pixel: " << cal.GetPixId() << endl;
827 bad.SetUncalibrated(MBadPixelsPix::kDeviatingNumPhes);
828 return kFALSE;
829 }
830
831 return kTRUE;
832}
833
834
835
836// -----------------------------------------------------------------------------------
837//
838// Sets pixel to MBadPixelsPix::kUnsuitableRun, if one of the following flags is set:
839// - MBadPixelsPix::kChargeIsPedestal
840// - MBadPixelsPix::kChargeRelErrNotValid
841// - MBadPixelsPix::kChargeSigmaNotValid
842// - MBadPixelsPix::kMeanTimeInFirstBin
843// - MBadPixelsPix::kMeanTimeInLast2Bins
844// - MBadPixelsPix::kDeviatingNumPhes
845//
846// - Call MCalibrationPix::SetExcluded() for the bad pixels
847//
848void MCalibrationChargeCalc::FinalizeBadPixels()
849{
850
851 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
852 {
853
854 MBadPixelsPix &bad = (*fBadPixels)[i];
855 MCalibrationPix &pix = (*fCam)[i];
856
857 if (bad.IsUncalibrated( MBadPixelsPix::kChargeIsPedestal))
858 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
859
860 if (bad.IsUncalibrated( MBadPixelsPix::kChargeErrNotValid ))
861 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
862
863 if (bad.IsUncalibrated( MBadPixelsPix::kChargeRelErrNotValid ))
864 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
865
866 if (bad.IsUncalibrated( MBadPixelsPix::kChargeSigmaNotValid ))
867 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
868
869 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin ))
870 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
871
872 if (bad.IsUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins ))
873 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
874
875 if (bad.IsUncalibrated( MBadPixelsPix::kDeviatingNumPhes ))
876 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
877
878 if (bad.IsUnsuitable( MBadPixelsPix::kUnsuitableRun ))
879 pix.SetExcluded();
880
881 }
882}
883
884// ------------------------------------------------------------------------
885//
886//
887// First loop: Calculate a mean and mean RMS of photo-electrons per area index
888// Include only pixels which are not MBadPixelsPix::kUnsuitableRun or
889// MBadPixelsPix::kUnreliableRun (see FinalizeBadPixels()) and set
890// MCalibrationChargePix::SetFFactorMethodValid(kFALSE) in that case.
891//
892// Second loop: Get weighted mean number of photo-electrons and its RMS including
893// only pixels with flag MCalibrationChargePix::IsFFactorMethodValid()
894// and further exclude those deviating by more than fPheErrLimit mean
895// sigmas from the mean (obtained in first loop). Set
896// MBadPixelsPix::kDeviatingNumPhes if excluded.
897//
898// Set weighted mean and variance of photo-electrons per area index in:
899// average area pixels of MCalibrationChargeCam (obtained from:
900// MCalibrationChargeCam::GetAverageArea() )
901//
902// Set weighted mean and variance of photo-electrons per sector in:
903// average sector pixels of MCalibrationChargeCam (obtained from:
904// MCalibrationChargeCam::GetAverageSector() )
905//
906Bool_t MCalibrationChargeCalc::FinalizeFFactorMethod()
907{
908
909 const UInt_t npixels = fGeom->GetNumPixels();
910 const UInt_t nareas = fGeom->GetNumAreas();
911 const UInt_t nsectors = fGeom->GetNumSectors();
912
913 Float_t lowlim [nareas];
914 Float_t upplim [nareas];
915 Float_t areavars [nareas];
916 Float_t areaweights [nareas], sectorweights [nsectors];
917 Float_t areaphes [nareas], sectorphes [nsectors];
918 Int_t numareavalid[nareas], numsectorvalid[nsectors];
919
920 memset(lowlim ,0, nareas * sizeof(Float_t));
921 memset(upplim ,0, nareas * sizeof(Float_t));
922 memset(areaphes ,0, nareas * sizeof(Float_t));
923 memset(areavars ,0, nareas * sizeof(Float_t));
924 memset(areaweights ,0, nareas * sizeof(Float_t));
925 memset(numareavalid ,0, nareas * sizeof(Int_t ));
926 memset(sectorweights ,0, nsectors * sizeof(Float_t));
927 memset(sectorphes ,0, nsectors * sizeof(Float_t));
928 memset(numsectorvalid,0, nsectors * sizeof(Int_t ));
929
930 //
931 // First loop: Get mean number of photo-electrons and the RMS
932 // The loop is only to recognize later pixels with very deviating numbers
933 //
934 MHCamera camphes(*fGeom,"Camphes","Phes in Camera");
935
936 for (UInt_t i=0; i<npixels; i++)
937 {
938
939 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam) [i];
940 MBadPixelsPix &bad = (*fBadPixels)[i];
941
942 if (!pix.IsFFactorMethodValid())
943 continue;
944
945 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
946 {
947 pix.SetFFactorMethodValid(kFALSE);
948 continue;
949 }
950
951 // if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
952 // continue;
953
954 const Float_t nphe = pix.GetPheFFactorMethod();
955 const Float_t nvar = pix.GetPheFFactorMethod()*pix.GetPheFFactorMethod();
956 const Int_t aidx = (*fGeom)[i].GetAidx();
957
958 camphes.Fill(i,nphe);
959 camphes.SetUsed(i);
960
961 areaphes [aidx] += nphe;
962 areavars [aidx] += nvar;
963 numareavalid[aidx] ++;
964 }
965
966 for (UInt_t i=0; i<nareas; i++)
967 {
968 if (numareavalid[i] == 0)
969 {
970 *fLog << warn << GetDescriptor() << ": No pixels with valid number of photo-electrons found "
971 << "in area index: " << i << endl;
972 continue;
973 }
974
975 areaphes[i] = areaphes[i] / numareavalid[i];
976 areavars[i] = (areavars[i] - areaphes[i]*areaphes[i]/numareavalid[i]) / (numareavalid[i]-1.);
977
978 if (areavars[i] > 0.)
979 areavars[i] = TMath::Sqrt(areavars[i]);
980 else
981 {
982 *fLog << warn << GetDescriptor() << ": No pixels with valid variance of photo-electrons found "
983 << "in area index: " << i << endl;
984 continue;
985 }
986
987 lowlim [i] = areaphes[i] - fPheErrLimit*TMath::Sqrt(areavars[i]);
988 upplim [i] = areaphes[i] + fPheErrLimit*TMath::Sqrt(areavars[i]);
989
990 TArrayI area(1);
991 area[0] = i;
992
993 TH1D *hist = camphes.ProjectionS(TArrayI(),area,"_py",100);
994 hist->Fit("gaus","Q");
995 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
996 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
997 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
998
999 if (ndf < 2)
1000 {
1001 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1002 << "in the camera with area index: " << i << endl;
1003 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1004 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1005 delete hist;
1006 continue;
1007 }
1008
1009 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1010
1011 if (prob < 0.001)
1012 {
1013 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the number of photo-electrons "
1014 << "in the camera with area index: " << i << endl;
1015 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1016 << " is smaller than 0.001 " << endl;
1017 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1018 delete hist;
1019 continue;
1020 }
1021
1022 *fLog << inf << endl;
1023 *fLog << inf << GetDescriptor() << ": Mean Number of photo-electrons "
1024 << "in the camera with area index: " << i << ": "
1025 << Form("%4.2f%s%4.2f",mean,"+-",sigma) << endl;
1026 *fLog << inf << endl;
1027
1028 lowlim [i] = mean - fPheErrLimit*sigma;
1029 upplim [i] = mean + fPheErrLimit*sigma;
1030
1031 delete hist;
1032 }
1033
1034 memset(numareavalid,0,nareas*sizeof(Int_t));
1035 memset(areaphes ,0,nareas*sizeof(Int_t));
1036 memset(areavars ,0,nareas*sizeof(Int_t));
1037
1038 //
1039 // Second loop: Get mean number of photo-electrons and its RMS excluding
1040 // pixels deviating by more than fPheErrLimit sigma.
1041 // Set the conversion factor FADC counts to photo-electrons
1042 //
1043 for (UInt_t i=0; i<npixels; i++)
1044 {
1045
1046 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1047
1048 if (!pix.IsFFactorMethodValid())
1049 continue;
1050
1051 const Float_t nvar = pix.GetPheFFactorMethodVar();
1052
1053 if (nvar <= 0.)
1054 {
1055 pix.SetFFactorMethodValid(kFALSE);
1056 continue;
1057 }
1058
1059 MBadPixelsPix &bad = (*fBadPixels)[i];
1060
1061 const Int_t aidx = (*fGeom)[i].GetAidx();
1062 const Int_t sector = (*fGeom)[i].GetSector();
1063 const Float_t nphe = pix.GetPheFFactorMethod();
1064
1065 if ( nphe < lowlim[aidx] || nphe > upplim[aidx] )
1066 {
1067 *fLog << warn << GetDescriptor() << ": Deviating number of photo-electrons: "
1068 << Form("%4.2f",nphe) << " out of accepted limits: ["
1069 << Form("%4.2f%s%4.2f",lowlim[aidx],",",upplim[aidx]) << "] in pixel " << i << endl;
1070 bad.SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1071 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1072 pix.SetFFactorMethodValid(kFALSE);
1073 continue;
1074 }
1075
1076 // const Float_t weight = 1./nvar;
1077 // areaweights [aidx] += weight;
1078 // areaphes [aidx] += weight*nphe;
1079 areaweights [aidx] += nphe*nphe;
1080 areaphes [aidx] += nphe;
1081 numareavalid [aidx] ++;
1082
1083 if (aidx == 0)
1084 fNumInnerFFactorMethodUsed++;
1085 // sectorweights [sector] += weight;
1086 // sectorphes [sector] += weight*nphe;
1087 sectorweights [sector] += nphe*nphe;
1088 sectorphes [sector] += nphe;
1089 numsectorvalid[sector] ++;
1090 }
1091
1092 for (UInt_t aidx=0; aidx<nareas; aidx++)
1093 {
1094
1095 MCalibrationChargePix &apix = (MCalibrationChargePix&)fCam->GetAverageArea(aidx);
1096
1097 areaweights[aidx] -= areaphes[aidx]*areaphes[aidx]/numareavalid[aidx];
1098 areaphes[aidx] /= numareavalid[aidx];
1099 areaweights[aidx] /= numareavalid[aidx]-1.;
1100
1101 if (areaweights[aidx] <= 0. || areaphes[aidx] <= 0.)
1102 {
1103 *fLog << warn << " Mean number of phe's from area index " << aidx << " cannot be calculated: "
1104 << " Mean number of phes: " << areaphes[aidx]
1105 << " Variance: " << areaweights[aidx] << endl;
1106 apix.SetFFactorMethodValid(kFALSE);
1107 continue;
1108 }
1109
1110 *fLog << inf << "Replacing number photo-electrons of average area idx " << aidx << ": "
1111 << Form("%5.3f%s%5.3f",apix.GetPheFFactorMethod()," +- ",apix.GetPheFFactorMethodErr()) << endl;
1112 *fLog << inf << " by average number of photo-electrons from area idx " << aidx << ": "
1113 << Form("%5.3f%s%5.3f",areaphes[aidx]," +- ",TMath::Sqrt(areaweights[aidx])) << endl;
1114
1115 // apix.SetPheFFactorMethod ( areaphes[aidx]/ areaweights[aidx] );
1116 // apix.SetPheFFactorMethodVar( 1. / areaweights[aidx] );
1117 apix.SetPheFFactorMethod ( areaphes[aidx] );
1118 apix.SetPheFFactorMethodVar( areaweights[aidx] / numareavalid[aidx] );
1119 apix.SetFFactorMethodValid ( kTRUE );
1120
1121 }
1122
1123 for (UInt_t sector=0; sector<nsectors; sector++)
1124 {
1125
1126 sectorweights[sector] -= sectorphes[sector]*sectorphes[sector]/numsectorvalid[sector];
1127 sectorphes[sector] /= numsectorvalid[sector];
1128 sectorweights[sector] /= numsectorvalid[sector]-1.;
1129
1130 MCalibrationChargePix &spix = (MCalibrationChargePix&)fCam->GetAverageSector(sector);
1131
1132 if (sectorweights[sector] <= 0. || sectorphes[sector] <= 0.)
1133 {
1134 *fLog << warn << " Mean number of phe's from sector " << sector << " cannot be calculated: "
1135 << " Mean number of phes: " << sectorphes[sector]
1136 << " Variance: " << sectorweights[sector] << endl;
1137 spix.SetFFactorMethodValid(kFALSE);
1138 continue;
1139 }
1140
1141 *fLog << inf << "Replacing number photo-electrons of average sector " << sector << ": "
1142 << Form("%5.3f%s%5.3f",spix.GetPheFFactorMethod()," +- ",spix.GetPheFFactorMethodErr()) << endl;
1143 *fLog << inf << " by average number photo-electrons from sector " << sector << ": "
1144 << Form("%5.3f%s%5.3f",sectorphes[sector]," +- ",TMath::Sqrt(sectorweights[sector])) << endl;
1145
1146 // spix.SetPheFFactorMethod ( sectorphes[sector]/ sectorweights[sector] );
1147 // spix.SetPheFFactorMethodVar( 1. / sectorweights[sector] );
1148 spix.SetPheFFactorMethod ( sectorphes[sector] );
1149 spix.SetPheFFactorMethodVar( sectorweights[sector] / numsectorvalid[sector]);
1150 spix.SetFFactorMethodValid ( kTRUE );
1151
1152 }
1153
1154 return kTRUE;
1155}
1156
1157
1158// ------------------------------------------------------------------------
1159//
1160// Returns kFALSE if pointer to MCalibrationChargeBlindPix is NULL
1161//
1162// The check returns kFALSE if:
1163//
1164// 1) fLambda and fLambdaCheck are separated relatively to each other by more than fLambdaCheckLimit
1165// 2) BlindPixel has an fLambdaErr greater than fLambdaErrLimit
1166//
1167// Calls:
1168// - MCalibrationChargeBlindPix::CalcFluxInsidePlexiglass()
1169//
1170Bool_t MCalibrationChargeCalc::FinalizeBlindPixel()
1171{
1172
1173 if (!fBlindPixel)
1174 return kFALSE;
1175
1176 const Float_t lambda = fBlindPixel->GetLambda();
1177 const Float_t lambdaerr = fBlindPixel->GetLambdaErr();
1178 const Float_t lambdacheck = fBlindPixel->GetLambdaCheck();
1179
1180 if (2.*(lambdacheck-lambda)/(lambdacheck+lambda) > fLambdaCheckLimit)
1181 {
1182 *fLog << warn << GetDescriptor()
1183 << Form("%s%4.2f%s%4.2f%s%4.2f%s",": Lambda: ",lambda," and Lambda-Check: ",
1184 lambdacheck," differ by more than ",fLambdaCheckLimit," in the Blind Pixel ")
1185 << endl;
1186 return kFALSE;
1187 }
1188
1189 if (lambdaerr > fLambdaErrLimit)
1190 {
1191 *fLog << warn << GetDescriptor()
1192 << Form("%s%4.2f%s%4.2f%s",": Error of Fitted Lambda: ",lambdaerr," is greater than ",
1193 fLambdaErrLimit," in Blind Pixel ") << endl;
1194 return kFALSE;
1195 }
1196
1197 if (!fBlindPixel->CalcFluxInsidePlexiglass())
1198 {
1199 *fLog << warn << "Could not calculate the flux of photons from the Blind Pixel, "
1200 << "will skip Blind Pixel Calibration " << endl;
1201 return kFALSE;
1202 }
1203
1204 return kTRUE;
1205}
1206
1207// ------------------------------------------------------------------------
1208//
1209// Returns kFALSE if pointer to MCalibrationChargePINDiode is NULL
1210//
1211// The check returns kFALSE if:
1212//
1213// 1) PINDiode has a fitted charge smaller than fChargeLimit*PedRMS
1214// 2) PINDiode has a fit error smaller than fChargeErrLimit
1215// 3) PINDiode has a fitted charge smaller its fChargeRelErrLimit times its charge error
1216// 4) PINDiode has a charge sigma smaller than its Pedestal RMS
1217//
1218// Calls:
1219// - MCalibrationChargePINDiode::CalcFluxOutsidePlexiglass()
1220//
1221Bool_t MCalibrationChargeCalc::FinalizePINDiode()
1222{
1223
1224 if (!fPINDiode)
1225 return kFALSE;
1226
1227 if (fPINDiode->GetMean() < fChargeLimit*fPINDiode->GetPedRms())
1228 {
1229 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1230 << fChargeLimit << " Pedestal RMS in PINDiode " << endl;
1231 return kFALSE;
1232 }
1233
1234 if (fPINDiode->GetMeanErr() < fChargeErrLimit)
1235 {
1236 *fLog << warn << GetDescriptor() << ": Error of Fitted Charge is smaller than "
1237 << fChargeErrLimit << " in PINDiode " << endl;
1238 return kFALSE;
1239 }
1240
1241 if (fPINDiode->GetMean() < fChargeRelErrLimit*fPINDiode->GetMeanErr())
1242 {
1243 *fLog << warn << GetDescriptor() << ": Fitted Charge is smaller than "
1244 << fChargeRelErrLimit << "* its error in PINDiode " << endl;
1245 return kFALSE;
1246 }
1247
1248 if (fPINDiode->GetSigma() < fPINDiode->GetPedRms())
1249 {
1250 *fLog << warn << GetDescriptor()
1251 << ": Sigma of Fitted Charge smaller than Pedestal RMS in PINDiode " << endl;
1252 return kFALSE;
1253 }
1254
1255
1256 if (!fPINDiode->CalcFluxOutsidePlexiglass())
1257 {
1258 *fLog << warn << "Could not calculate the flux of photons from the PIN Diode, "
1259 << "will skip PIN Diode Calibration " << endl;
1260 return kFALSE;
1261 }
1262
1263 return kTRUE;
1264}
1265
1266// ------------------------------------------------------------------------
1267//
1268// Calculate the average number of photons outside the plexiglass with the
1269// formula:
1270//
1271// av.Num.photons(area index) = av.Num.Phes(area index)
1272// / MCalibrationQEPix::GetDefaultQE(fPulserColor)
1273// / MCalibrationQEPix::GetPMTCollectionEff()
1274// / MCalibrationQEPix::GetLightGuidesEff(fPulserColor)
1275// / MCalibrationQECam::GetPlexiglassQE()
1276//
1277// Calculate the variance on the average number of photons assuming that the error on the
1278// Quantum efficiency is reduced by the number of used inner pixels, but the rest of the
1279// values keeps it ordinary error since it is systematic.
1280//
1281// Loop over pixels:
1282//
1283// - Continue, if not MCalibrationChargePix::IsFFactorMethodValid() and set:
1284// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor)
1285//
1286// - Call MCalibrationChargePix::CalcMeanFFactor(av.Num.photons) and set:
1287// MCalibrationQEPix::SetFFactorMethodValid(kFALSE,fPulserColor) if not succesful
1288//
1289// - Calculate the quantum efficiency with the formula:
1290//
1291// QE = ( Num.Phes / av.Num.photons ) * MGeomCam::GetPixRatio()
1292//
1293// - Set QE in MCalibrationQEPix::SetQEFFactor ( QE, fPulserColor );
1294//
1295// - Set Variance of QE in MCalibrationQEPix::SetQEFFactorVar ( Variance, fPulserColor );
1296// - Set bit MCalibrationQEPix::SetFFactorMethodValid(kTRUE,fPulserColor)
1297//
1298// - Call MCalibrationQEPix::UpdateFFactorMethod()
1299//
1300void MCalibrationChargeCalc::FinalizeFFactorQECam()
1301{
1302
1303 if (fNumInnerFFactorMethodUsed < 2)
1304 {
1305 *fLog << warn << GetDescriptor()
1306 << ": Could not calculate F-Factor Method: Less than 2 inner pixels valid! " << endl;
1307 return;
1308 }
1309
1310 MCalibrationChargePix &avpix = (MCalibrationChargePix&)fCam->GetAverageArea(0);
1311 MCalibrationQEPix &qepix = (MCalibrationQEPix&) fQECam->GetAverageArea(0);
1312
1313 const Float_t avphotons = avpix.GetPheFFactorMethod()
1314 / qepix.GetDefaultQE(fPulserColor)
1315 / qepix.GetPMTCollectionEff()
1316 / qepix.GetLightGuidesEff(fPulserColor)
1317 / fQECam->GetPlexiglassQE();
1318
1319 const Float_t avphotrelvar = avpix.GetPheFFactorMethodRelVar()
1320 + qepix.GetDefaultQERelVar(fPulserColor) / fNumInnerFFactorMethodUsed
1321 + qepix.GetPMTCollectionEffRelVar()
1322 + qepix.GetLightGuidesEffRelVar(fPulserColor)
1323 + fQECam->GetPlexiglassQERelVar();
1324
1325 const UInt_t nareas = fGeom->GetNumAreas();
1326
1327 //
1328 // Set the results in the MCalibrationChargeCam
1329 //
1330 fCam->SetNumPhotonsFFactorMethod (avphotons);
1331 if (avphotrelvar > 0.)
1332 fCam->SetNumPhotonsFFactorMethodErr(TMath::Sqrt( avphotrelvar * avphotons * avphotons));
1333
1334 Float_t lowlim [nareas];
1335 Float_t upplim [nareas];
1336 Float_t avffactorphotons [nareas];
1337 Float_t avffactorphotvar [nareas];
1338 Int_t numffactor [nareas];
1339
1340 memset(lowlim ,0, nareas * sizeof(Float_t));
1341 memset(upplim ,0, nareas * sizeof(Float_t));
1342 memset(avffactorphotons,0, nareas * sizeof(Float_t));
1343 memset(avffactorphotvar,0, nareas * sizeof(Float_t));
1344 memset(numffactor ,0, nareas * sizeof(Int_t));
1345
1346 const UInt_t npixels = fGeom->GetNumPixels();
1347
1348 MHCamera camffactor(*fGeom,"Camffactor","F-Factor in Camera");
1349
1350 for (UInt_t i=0; i<npixels; i++)
1351 {
1352
1353 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1354 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1355 MBadPixelsPix &bad = (*fBadPixels)[i];
1356
1357 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1358 continue;
1359
1360 const Float_t photons = avphotons / fGeom->GetPixRatio(i);
1361 const Float_t qe = pix.GetPheFFactorMethod() / photons ;
1362
1363 if (!pix.CalcMeanFFactor( photons , avphotrelvar ))
1364 {
1365 (*fBadPixels)[i].SetUncalibrated( MBadPixelsPix::kDeviatingNumPhes );
1366 continue;
1367 }
1368
1369 const Float_t qerelvar = avphotrelvar + pix.GetPheFFactorMethodRelVar();
1370
1371 qepix.SetQEFFactor ( qe , fPulserColor );
1372 qepix.SetQEFFactorVar ( qerelvar*qe*qe, fPulserColor );
1373 qepix.SetFFactorMethodValid( kTRUE , fPulserColor );
1374
1375 if (!qepix.UpdateFFactorMethod())
1376 *fLog << warn << GetDescriptor()
1377 << ": Cannot update Quantum efficiencies with the F-Factor Method" << endl;
1378
1379 const Int_t aidx = (*fGeom)[i].GetAidx();
1380 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1381
1382 camffactor.Fill(i,ffactor);
1383 camffactor.SetUsed(i);
1384
1385 avffactorphotons[aidx] += ffactor;
1386 avffactorphotvar[aidx] += pix.GetMeanFFactorFADC2PhotVar();
1387 numffactor[aidx]++;
1388 }
1389
1390 for (UInt_t i=0; i<nareas; i++)
1391 {
1392 avffactorphotons[i] /= numffactor[i];
1393 avffactorphotvar[i] /= numffactor[i];
1394
1395 lowlim [i] = 1.1; // Lowest known F-Factor of a PMT
1396 upplim [i] = avffactorphotons[i] + fFFactorErrLimit*TMath::Sqrt(avffactorphotvar[i]);
1397
1398 TArrayI area(1);
1399 area[0] = i;
1400
1401 TH1D *hist = camffactor.ProjectionS(TArrayI(),area,"_py",100);
1402 hist->Fit("gaus","Q");
1403 const Float_t mean = hist->GetFunction("gaus")->GetParameter(1);
1404 const Float_t sigma = hist->GetFunction("gaus")->GetParameter(2);
1405 const Int_t ndf = hist->GetFunction("gaus")->GetNDF();
1406
1407 if (ndf < 2)
1408 {
1409 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1410 << "in the camera with area index: " << i << endl;
1411 *fLog << warn << GetDescriptor() << ": Number of dof.: " << ndf << " is smaller than 2 " << endl;
1412 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1413 delete hist;
1414 continue;
1415 }
1416
1417 const Double_t prob = hist->GetFunction("gaus")->GetProb();
1418
1419 if (prob < 0.001)
1420 {
1421 *fLog << warn << GetDescriptor() << ": Cannot use a Gauss fit to the F-Factor "
1422 << "in the camera with area index: " << i << endl;
1423 *fLog << warn << GetDescriptor() << ": Fit probability " << prob
1424 << " is smaller than 0.001 " << endl;
1425 *fLog << warn << GetDescriptor() << ": Will use the simple mean and rms " << endl;
1426 delete hist;
1427 continue;
1428 }
1429
1430 *fLog << inf << endl;
1431 *fLog << inf << GetDescriptor() << ": Mean F-Factor "
1432 << "in the camera with area index: " << i << ": "
1433 << Form("%4.2f%s%4.2f",mean,"+-",sigma) << endl;
1434 *fLog << inf << endl;
1435
1436 lowlim [i] = mean - fFFactorErrLimit*sigma;
1437 upplim [i] = mean + fFFactorErrLimit*sigma;
1438
1439 delete hist;
1440 }
1441
1442 for (UInt_t i=0; i<npixels; i++)
1443 {
1444
1445 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1446 MBadPixelsPix &bad = (*fBadPixels)[i];
1447
1448 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1449 continue;
1450
1451 const Float_t ffactor = pix.GetMeanFFactorFADC2Phot();
1452 const Int_t aidx = (*fGeom)[i].GetAidx();
1453
1454 if ( ffactor < lowlim[aidx] || ffactor > upplim[aidx] )
1455 {
1456 *fLog << warn << GetDescriptor() << ": Deviating F-Factor: "
1457 << Form("%4.2f",ffactor) << " out of accepted limits: ["
1458 << Form("%4.2f%s%4.2f",lowlim[aidx],",",upplim[aidx]) << "] in pixel " << i << endl;
1459 bad.SetUncalibrated( MBadPixelsPix::kDeviatingFFactor );
1460 bad.SetUnsuitable ( MBadPixelsPix::kUnsuitableRun );
1461 }
1462 }
1463
1464 for (UInt_t i=0; i<npixels; i++)
1465 {
1466
1467 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1468 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1469 MBadPixelsPix &bad = (*fBadPixels)[i];
1470
1471 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1472 {
1473 qepix.SetFFactorMethodValid(kFALSE,fPulserColor);
1474 pix.SetFFactorMethodValid(kFALSE);
1475 pix.SetExcluded();
1476 continue;
1477 }
1478 }
1479}
1480
1481
1482// ------------------------------------------------------------------------
1483//
1484// Loop over pixels:
1485//
1486// - Continue, if not MCalibrationChargeBlindPix::IsFluxInsidePlexiglassAvailable() and set:
1487// MCalibrationQEPix::SetBlindPixelMethodValid(kFALSE,fPulserColor)
1488//
1489// - Calculate the quantum efficiency with the formula:
1490//
1491// QE = Num.Phes / MCalibrationChargeBlindPix::GetFluxInsidePlexiglass()
1492// / MGeomPix::GetA() * MCalibrationQECam::GetPlexiglassQE()
1493//
1494// - Set QE in MCalibrationQEPix::SetQEBlindPixel ( QE, fPulserColor );
1495// - Set Variance of QE in MCalibrationQEPix::SetQEBlindPixelVar ( Variance, fPulserColor );
1496// - Set bit MCalibrationQEPix::SetBlindPixelMethodValid(kTRUE,fPulserColor)
1497//
1498// - Call MCalibrationQEPix::UpdateBlindPixelMethod()
1499//
1500void MCalibrationChargeCalc::FinalizeBlindPixelQECam()
1501{
1502
1503 const UInt_t npixels = fGeom->GetNumPixels();
1504
1505 //
1506 // Set the results in the MCalibrationChargeCam
1507 //
1508 if (fBlindPixel)
1509 {
1510 if (fBlindPixel->IsFluxInsidePlexiglassAvailable())
1511 {
1512
1513 const Float_t photons = fBlindPixel->GetFluxInsidePlexiglass() * (*fGeom)[0].GetA()
1514 / fQECam->GetPlexiglassQE();
1515 fCam->SetNumPhotonsBlindPixelMethod(photons);
1516
1517 const Float_t photrelvar = fBlindPixel->GetFluxInsidePlexiglassRelVar()
1518 + fQECam->GetPlexiglassQERelVar();
1519 if (photrelvar > 0.)
1520 fCam->SetNumPhotonsBlindPixelMethodErr(TMath::Sqrt( photrelvar * photons * photons));
1521 }
1522 }
1523 //
1524 // With the knowledge of the overall photon flux, calculate the
1525 // quantum efficiencies after the Blind Pixel and PIN Diode method
1526 //
1527 for (UInt_t i=0; i<npixels; i++)
1528 {
1529
1530 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1531
1532 if (!fBlindPixel)
1533 {
1534 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1535 continue;
1536 }
1537
1538 if (!fBlindPixel->IsFluxInsidePlexiglassAvailable())
1539 {
1540 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1541 continue;
1542 }
1543
1544 MBadPixelsPix &bad = (*fBadPixels)[i];
1545 if (bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1546 {
1547 qepix.SetBlindPixelMethodValid(kFALSE, fPulserColor);
1548 continue;
1549 }
1550
1551 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1552 MGeomPix &geo = (*fGeom)[i];
1553
1554 const Float_t qe = pix.GetPheFFactorMethod()
1555 / fBlindPixel->GetFluxInsidePlexiglass()
1556 / geo.GetA()
1557 * fQECam->GetPlexiglassQE();
1558
1559 const Float_t qerelvar = fBlindPixel->GetFluxInsidePlexiglassRelVar()
1560 + fQECam->GetPlexiglassQERelVar()
1561 + pix.GetPheFFactorMethodRelVar();
1562
1563 qepix.SetQEBlindPixel ( qe , fPulserColor );
1564 qepix.SetQEBlindPixelVar ( qerelvar*qe*qe, fPulserColor );
1565 qepix.SetBlindPixelMethodValid( kTRUE , fPulserColor );
1566
1567 if (!qepix.UpdateBlindPixelMethod())
1568 *fLog << warn << GetDescriptor()
1569 << ": Cannot update Quantum efficiencies with the Blind Pixel Method" << endl;
1570 }
1571}
1572
1573// ------------------------------------------------------------------------
1574//
1575// Loop over pixels:
1576//
1577// - Continue, if not MCalibrationChargePINDiode::IsFluxOutsidePlexiglassAvailable() and set:
1578// MCalibrationQEPix::SetPINDiodeMethodValid(kFALSE,fPulserColor)
1579//
1580// - Calculate the quantum efficiency with the formula:
1581//
1582// QE = Num.Phes / MCalibrationChargePINDiode::GetFluxOutsidePlexiglass() / MGeomPix::GetA()
1583//
1584// - Set QE in MCalibrationQEPix::SetQEPINDiode ( QE, fPulserColor );
1585// - Set Variance of QE in MCalibrationQEPix::SetQEPINDiodeVar ( Variance, fPulserColor );
1586// - Set bit MCalibrationQEPix::SetPINDiodeMethodValid(kTRUE,fPulserColor)
1587//
1588// - Call MCalibrationQEPix::UpdatePINDiodeMethod()
1589//
1590void MCalibrationChargeCalc::FinalizePINDiodeQECam()
1591{
1592
1593 const UInt_t npixels = fGeom->GetNumPixels();
1594
1595 //
1596 // With the knowledge of the overall photon flux, calculate the
1597 // quantum efficiencies after the PIN Diode method
1598 //
1599 for (UInt_t i=0; i<npixels; i++)
1600 {
1601
1602 MCalibrationQEPix &qepix = (MCalibrationQEPix&) (*fQECam)[i];
1603
1604 if (!fPINDiode)
1605 {
1606 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1607 continue;
1608 }
1609
1610 if (!fPINDiode->IsFluxOutsidePlexiglassAvailable())
1611 {
1612 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1613 continue;
1614 }
1615
1616 MBadPixelsPix &bad = (*fBadPixels)[i];
1617
1618 if (!bad.IsUnsuitable (MBadPixelsPix::kUnsuitableRun))
1619 {
1620 qepix.SetPINDiodeMethodValid(kFALSE, fPulserColor);
1621 continue;
1622 }
1623
1624 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1625 MGeomPix &geo = (*fGeom)[i];
1626
1627 const Float_t qe = pix.GetPheFFactorMethod()
1628 / fPINDiode->GetFluxOutsidePlexiglass()
1629 / geo.GetA();
1630
1631 const Float_t qerelvar = fPINDiode->GetFluxOutsidePlexiglassRelVar() + pix.GetPheFFactorMethodRelVar();
1632
1633 qepix.SetQEPINDiode ( qe , fPulserColor );
1634 qepix.SetQEPINDiodeVar ( qerelvar*qe*qe, fPulserColor );
1635 qepix.SetPINDiodeMethodValid( kTRUE , fPulserColor );
1636
1637 if (!qepix.UpdatePINDiodeMethod())
1638 *fLog << warn << GetDescriptor()
1639 << ": Cannot update Quantum efficiencies with the PIN Diode Method" << endl;
1640 }
1641}
1642
1643// -----------------------------------------------------------------------------------------------
1644//
1645// - Print out statistics about BadPixels of type UnsuitableType_t
1646// - store numbers of bad pixels of each type in fCam
1647//
1648void MCalibrationChargeCalc::FinalizeUnsuitablePixels()
1649{
1650
1651 *fLog << inf << endl;
1652 *fLog << GetDescriptor() << ": Charge Calibration status:" << endl;
1653 *fLog << dec << setfill(' ');
1654
1655 const Int_t nareas = fGeom->GetNumAreas();
1656
1657 Int_t counts[nareas];
1658 memset(counts,0,nareas*sizeof(Int_t));
1659
1660 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1661 {
1662 MBadPixelsPix &bad = (*fBadPixels)[i];
1663 if (!bad.IsBad())
1664 {
1665 const Int_t aidx = (*fGeom)[i].GetAidx();
1666 counts[aidx]++;
1667 }
1668 }
1669
1670 if (fGeom->InheritsFrom("MGeomCamMagic"))
1671 *fLog << " " << setw(7) << "Successfully calibrated Pixels: "
1672 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1673
1674 memset(counts,0,nareas*sizeof(Int_t));
1675
1676 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1677 {
1678 MBadPixelsPix &bad = (*fBadPixels)[i];
1679 if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
1680 {
1681 const Int_t aidx = (*fGeom)[i].GetAidx();
1682 counts[aidx]++;
1683 }
1684 }
1685
1686 for (Int_t aidx=0; aidx<nareas; aidx++)
1687 fCam->SetNumUnsuitable(counts[aidx], aidx);
1688
1689 if (fGeom->InheritsFrom("MGeomCamMagic"))
1690 *fLog << " " << setw(7) << "Uncalibrated Pixels: "
1691 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1692
1693 memset(counts,0,nareas*sizeof(Int_t));
1694
1695 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1696 {
1697 MBadPixelsPix &bad = (*fBadPixels)[i];
1698 if (bad.IsUnsuitable(MBadPixelsPix::kUnreliableRun))
1699 {
1700 const Int_t aidx = (*fGeom)[i].GetAidx();
1701 counts[aidx]++;
1702 }
1703 }
1704
1705 for (Int_t aidx=0; aidx<nareas; aidx++)
1706 fCam->SetNumUnreliable(counts[aidx], aidx);
1707
1708 *fLog << " " << setw(7) << "Unreliable Pixels: "
1709 << Form("%s%3i%s%3i","Inner: ",counts[0]," Outer: ",counts[1]) << endl;
1710
1711}
1712
1713// -----------------------------------------------------------------------------------------------
1714//
1715// Print out statistics about BadPixels of type UncalibratedType_t
1716//
1717void MCalibrationChargeCalc::PrintUncalibrated(MBadPixelsPix::UncalibratedType_t typ, const char *text) const
1718{
1719
1720 UInt_t countinner = 0;
1721 UInt_t countouter = 0;
1722
1723 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1724 {
1725 MBadPixelsPix &bad = (*fBadPixels)[i];
1726 if (bad.IsUncalibrated(typ))
1727 {
1728 if (fGeom->GetPixRatio(i) == 1.)
1729 countinner++;
1730 else
1731 countouter++;
1732 }
1733 }
1734
1735 *fLog << " " << setw(7) << text
1736 << Form("%s%3i%s%3i","Inner: ",countinner," Outer: ",countouter) << endl;
1737}
1738
1739// --------------------------------------------------------------------------
1740//
1741// Set the path for output file
1742//
1743void MCalibrationChargeCalc::SetOutputPath(TString path)
1744{
1745 fOutputPath = path;
1746 if (fOutputPath.EndsWith("/"))
1747 fOutputPath = fOutputPath(0, fOutputPath.Length()-1);
1748}
1749
1750void MCalibrationChargeCalc::SetOutputFile(TString file)
1751{
1752 fOutputFile = file;
1753}
1754
1755// --------------------------------------------------------------------------
1756//
1757// Get the output file
1758//
1759const char* MCalibrationChargeCalc::GetOutputFile()
1760{
1761 return Form("%s/%s", (const char*)fOutputPath, (const char*)fOutputFile);
1762}
Note: See TracBrowser for help on using the repository browser.