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

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