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

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