source: trunk/MagicSoft/Mars/mcalib/MHCalibrationChargeCam.cc@ 3669

Last change on this file since 3669 was 3669, checked in by gaug, 20 years ago
*** empty log message ***
File size: 25.5 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// MHCalibrationChargeCam
27//
28// Fills the extracted signals of MExtractedSignalCam into the MHGausEvents-classes
29// MHCalibrationChargeHiGainPix and MHCalibrationChargeLoGainPix for every:
30//
31// - Pixel, stored in the TObjArray's MHCalibrationCam::fHiGainArray and
32// MHCalibrationCam::fLoGainArray
33//
34// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
35// stored in the TObjArray's MHCalibrationCam::fAverageHiGainAreas and
36// MHCalibrationCam::fAverageLoGainAreas
37//
38// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
39// stored in the TObjArray's MHCalibrationCam::fAverageHiGainSectors and
40// MHCalibrationCam::fAverageLoGainSectors
41//
42// Every signal is taken from MExtractedSignalCam and filled into a histogram and
43// an array, in order to perform a Fourier analysis (see MHGausEvents).
44// The signals are moreover averaged on an event-by-event basis and written into
45// the corresponding average pixels.
46//
47// Additionally, the (FADC slice) position of the maximum is stored in an Absolute
48// Arrival Time histogram. This histogram serves for a rough cross-check if the
49// signal does not lie at or outside the edges of the extraction window.
50//
51// The Charge histograms are fitted to a Gaussian, mean and sigma with its errors
52// and the fit probability are extracted. If none of these values are NaN's and
53// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
54// the fit is declared valid.
55// Otherwise, the fit is repeated within ranges of the previous mean
56// +- MHGausEvents::fPickupLimit (default: 5) sigma (see MHGausEvents::RepeatFit())
57// In case this does not make the fit valid, the histogram means and RMS's are
58// taken directly (see MHGausEvents::BypassFit()) and the following flags are set:
59// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted ) or
60// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted ) and
61// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
62//
63// Outliers of more than MHGausEvents::fPickupLimit (default: 5) sigmas
64// from the mean are counted as Pickup events (stored in MHGausEvents::fPickup)
65//
66// Unless more than fNumHiGainSaturationLimit (default: 1%) of the overall FADC
67// slices show saturation, the following flag is set:
68// - MCalibrationChargePix::SetHiGainSaturation();
69// In that case, the calibration constants are derived from the low-gain results.
70//
71// If more than fNumLoGainSaturationLimit (default: 1%) of the overall
72// low-gain FADC slices saturate, the following flags are set:
73// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturation ) and
74// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnsuitableRun )
75//
76// The class also fills arrays with the signal vs. event number, creates a fourier
77// spectrum and investigates if the projected fourier components follow an exponential
78// distribution. In case that the probability of the exponential fit is less than
79// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
80// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating ) or
81// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating ) and
82// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
83//
84// This same procedure is performed for the average pixels.
85//
86// The following results are written into MCalibrationChargeCam:
87//
88// - MCalibrationPix::SetHiGainSaturation()
89// - MCalibrationPix::SetHiGainMean()
90// - MCalibrationPix::SetHiGainMeanErr()
91// - MCalibrationPix::SetHiGainSigma()
92// - MCalibrationPix::SetHiGainSigmaErr()
93// - MCalibrationPix::SetHiGainProb()
94// - MCalibrationPix::SetHiGainNumPickup()
95//
96// - MCalibrationPix::SetLoGainMean()
97// - MCalibrationPix::SetLoGainMeanErr()
98// - MCalibrationPix::SetLoGainSigma()
99// - MCalibrationPix::SetLoGainSigmaErr()
100// - MCalibrationPix::SetLoGainProb()
101// - MCalibrationPix::SetLoGainNumPickup()
102//
103// - MCalibrationChargePix::SetAbsTimeMean()
104// - MCalibrationChargePix::SetAbsTimeRms()
105//
106// For all averaged areas, the fitted sigma is multiplied with the square root of
107// the number involved pixels in order to be able to compare it to the average of
108// sigmas in the camera.
109//
110/////////////////////////////////////////////////////////////////////////////
111#include "MHCalibrationChargeCam.h"
112#include "MHCalibrationCam.h"
113
114#include "MLog.h"
115#include "MLogManip.h"
116
117#include "MParList.h"
118
119#include "MHCalibrationChargeHiGainPix.h"
120#include "MHCalibrationChargeLoGainPix.h"
121#include "MHCalibrationChargePix.h"
122
123#include "MCalibrationChargeCam.h"
124#include "MCalibrationChargePix.h"
125
126#include "MGeomCam.h"
127#include "MGeomPix.h"
128
129#include "MHGausEvents.h"
130
131#include "MBadPixelsCam.h"
132#include "MBadPixelsPix.h"
133
134#include "MRawEvtData.h"
135#include "MRawEvtPixelIter.h"
136
137#include "MExtractedSignalCam.h"
138#include "MExtractedSignalPix.h"
139
140ClassImp(MHCalibrationChargeCam);
141
142using namespace std;
143
144const Float_t MHCalibrationChargeCam::fgNumHiGainSaturationLimit = 0.01;
145const Float_t MHCalibrationChargeCam::fgNumLoGainSaturationLimit = 0.005;
146// --------------------------------------------------------------------------
147//
148// Default Constructor.
149//
150// Sets:
151// - all pointers to NULL
152//
153// Initializes:
154// - fNumHiGainSaturationLimit to fgNumHiGainSaturationLimit
155// - fNumLoGainSaturationLimit to fgNumLoGainSaturationLimit
156//
157MHCalibrationChargeCam::MHCalibrationChargeCam(const char *name, const char *title)
158 : fRawEvt(NULL)
159{
160 fName = name ? name : "MHCalibrationChargeCam";
161 fTitle = title ? title : "Class to fill the calibration histograms ";
162
163 SetNumHiGainSaturationLimit();
164 SetNumLoGainSaturationLimit();
165}
166
167// --------------------------------------------------------------------------
168//
169// Gets the pointers to:
170// - MRawEvtData
171//
172Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
173{
174
175 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
176 if (!fRawEvt)
177 {
178 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
179 return kFALSE;
180 }
181
182 return kTRUE;
183}
184
185// --------------------------------------------------------------------------
186//
187// Gets or creates the pointers to:
188// - MCalibrationChargeCam
189//
190// Searches pointer to:
191// - MExtractedSignalCam
192//
193// Initializes, if empty to MGeomCam::GetNumPixels():
194// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
195//
196// Initializes, if empty to MGeomCam::GetNumAreas() for:
197// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
198//
199// Initializes, if empty to MGeomCam::GetNumSectors() for:
200// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
201//
202// Calls MHCalibrationCam::InitHists() for every entry in:
203// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
204// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
205// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
206//
207// Sets Titles and Names for the Charge Histograms:
208// - MHCalibrationCam::fAverageHiGainAreas
209// - MHCalibrationCam::fAverageHiGainSectors
210//
211// Sets number of bins to MHCalibrationCam::fAverageNbins for:
212// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
213// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
214//
215Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
216{
217
218 fCam = (MCalibrationCam*)pList->FindObject("MCalibrationChargeCam");
219 if (!fCam)
220 {
221 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationChargeCam"));
222 if (!fCam)
223 {
224 gLog << err << "Cannot find nor create MCalibrationChargeCam ... abort." << endl;
225 return kFALSE;
226 }
227 else
228 fCam->Init(*fGeom);
229 }
230
231 MExtractedSignalCam *signal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
232 if (!signal)
233 {
234 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
235 return kFALSE;
236 }
237
238 const Int_t npixels = fGeom->GetNumPixels();
239 const Int_t nsectors = fGeom->GetNumSectors();
240 const Int_t nareas = fGeom->GetNumAreas();
241
242 if (fHiGainArray->GetEntries()==0)
243 {
244 fHiGainArray->Expand(npixels);
245 for (Int_t i=0; i<npixels; i++)
246 {
247 (*fHiGainArray)[i] = new MHCalibrationChargeHiGainPix;
248 InitHists((*this)[i],(*fBadPixels)[i],i);
249 }
250 }
251
252 if (fLoGainArray->GetEntries()==0)
253 {
254 fLoGainArray->Expand(npixels);
255
256 for (Int_t i=0; i<npixels; i++)
257 {
258 (*fLoGainArray)[i] = new MHCalibrationChargeLoGainPix;
259 InitHists((*this)(i),(*fBadPixels)[i],i);
260 }
261
262 }
263
264 if (fAverageHiGainAreas->GetEntries()==0)
265 {
266 fAverageHiGainAreas->Expand(nareas);
267
268 for (Int_t j=0; j<nareas; j++)
269 {
270 (*fAverageHiGainAreas)[j] =
271 new MHCalibrationChargeHiGainPix("AverageHiGainArea",
272 "Average HiGain FADC sums area idx ");
273
274 InitHists(GetAverageHiGainArea(j),fCam->GetAverageBadArea(j),j);
275
276 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
277
278 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Area Idx ");
279 hist.SetNbins(fAverageNbins);
280 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Area Idx ");
281 }
282 }
283
284 if (fAverageLoGainAreas->GetEntries()==0)
285 {
286 fAverageLoGainAreas->Expand(nareas);
287
288 for (Int_t j=0; j<nareas; j++)
289 {
290 (*fAverageLoGainAreas)[j] =
291 new MHCalibrationChargeLoGainPix("AverageLoGainArea",
292 "Average LoGain FADC sums of pixel area idx ");
293
294 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
295
296 InitHists(hist,fCam->GetAverageBadArea(j),j);
297
298 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Area Idx ");
299 hist.SetNbins(fAverageNbins);
300 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Area Idx ");
301
302 }
303 }
304
305 if (fAverageHiGainSectors->GetEntries()==0)
306 {
307 fAverageHiGainSectors->Expand(nsectors);
308
309 for (Int_t j=0; j<nsectors; j++)
310 {
311 (*fAverageHiGainSectors)[j] =
312 new MHCalibrationChargeHiGainPix("AverageHiGainSector",
313 "Average HiGain FADC sums of pixel sector ");
314
315 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
316
317 InitHists(hist,fCam->GetAverageBadSector(j),j);
318
319 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Sector ");
320 hist.SetNbins(fAverageNbins);
321 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Sector ");
322
323 }
324 }
325
326 if (fAverageLoGainSectors->GetEntries()==0)
327 {
328 fAverageLoGainSectors->Expand(nsectors);
329
330 for (Int_t j=0; j<nsectors; j++)
331 {
332 (*fAverageLoGainSectors)[j] =
333 new MHCalibrationChargeLoGainPix("AverageLoGainSector",
334 "Average LoGain FADC sums of pixel sector ");
335
336 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
337
338 InitHists(hist,fCam->GetAverageBadSector(j),j);
339
340 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Sector ");
341 hist.SetNbins(fAverageNbins);
342 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Sector ");
343
344 }
345 }
346
347 return kTRUE;
348}
349
350
351// --------------------------------------------------------------------------
352//
353// Retrieves from MExtractedSignalCam:
354// - first used LoGain FADC slice
355//
356// Retrieves from MGeomCam:
357// - number of pixels
358// - number of pixel areas
359// - number of sectors
360//
361// For all TObjArray's (including the averaged ones), the following steps are performed:
362//
363// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
364// - MExtractedSignalPix::GetExtractedSignalHiGain();
365// - MExtractedSignalPix::GetExtractedSignalLoGain();
366//
367// 2) Set number of saturated slices (MHCalibrationChargePix::SetSaturated()) with:
368// - MExtractedSignalPix::GetNumHiGainSaturated();
369// - MExtractedSignalPix::GetNumLoGainSaturated();
370//
371// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
372// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
373// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
374//
375Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
376{
377
378 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
379 if (!signal)
380 {
381 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
382 return kFALSE;
383 }
384
385 const UInt_t npixels = fGeom->GetNumPixels();
386 const UInt_t nareas = fGeom->GetNumAreas();
387 const UInt_t nsectors = fGeom->GetNumSectors();
388 const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
389
390 Float_t sumhiarea [nareas], sumloarea [nareas], timehiarea [nareas], timeloarea [nareas];
391 Float_t sumhisector[nsectors], sumlosector[nsectors], timehisector[nsectors], timelosector[nsectors];
392 Int_t sathiarea [nareas], satloarea [nareas];
393 Int_t sathisector[nsectors], satlosector[nsectors];
394
395 memset(sumhiarea, 0, nareas * sizeof(Float_t));
396 memset(sumloarea, 0, nareas * sizeof(Float_t));
397 memset(timehiarea, 0, nareas * sizeof(Float_t));
398 memset(timeloarea, 0, nareas * sizeof(Float_t));
399 memset(sathiarea, 0, nareas * sizeof(Int_t ));
400 memset(satloarea, 0, nareas * sizeof(Int_t ));
401 memset(sumhisector, 0, nsectors*sizeof(Float_t));
402 memset(sumlosector, 0, nsectors*sizeof(Float_t));
403 memset(timehisector,0, nsectors*sizeof(Float_t));
404 memset(timelosector,0, nsectors*sizeof(Float_t));
405 memset(sathisector, 0, nsectors*sizeof(Int_t ));
406 memset(satlosector, 0, nsectors*sizeof(Int_t ));
407
408 for (UInt_t i=0; i<npixels; i++)
409 {
410
411 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
412 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
413
414 if (histhi.IsExcluded())
415 continue;
416
417 const MExtractedSignalPix &pix = (*signal)[i];
418
419 const Float_t sumhi = pix.GetExtractedSignalHiGain();
420 const Float_t sumlo = pix.GetExtractedSignalLoGain();
421
422 histhi.FillHistAndArray(sumhi);
423 histlo.FillHistAndArray(sumlo);
424
425 const Int_t sathi = (Int_t)pix.GetNumHiGainSaturated();
426 const Int_t satlo = (Int_t)pix.GetNumLoGainSaturated();
427
428 histhi.SetSaturated(sathi);
429 histlo.SetSaturated(satlo);
430
431 const Int_t aidx = (*fGeom)[i].GetAidx();
432 const Int_t sector = (*fGeom)[i].GetSector();
433
434 sumhiarea[aidx] += sumhi;
435 sumloarea[aidx] += sumlo;
436 sathiarea[aidx] += sathi;
437 satloarea[aidx] += satlo;
438
439 sumhisector[sector] += sumhi;
440 sumlosector[sector] += sumlo;
441 sathisector[sector] += sathi;
442 satlosector[sector] += satlo;
443 }
444
445 MRawEvtPixelIter pixel(fRawEvt);
446 while (pixel.Next())
447 {
448
449 const UInt_t pixid = pixel.GetPixelId();
450
451 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
452 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
453
454 if (histhi.IsExcluded())
455 continue;
456
457 const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
458 const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
459
460 histhi.FillAbsTime(timehi);
461 histlo.FillAbsTime(timelo);
462
463 const Int_t aidx = (*fGeom)[pixid].GetAidx();
464 const Int_t sector = (*fGeom)[pixid].GetSector();
465
466 timehiarea[aidx] += timehi;
467 timeloarea[aidx] += timelo;
468
469 timehisector[sector] += timehi;
470 timelosector[sector] += timelo;
471 }
472
473 for (UInt_t j=0; j<nareas; j++)
474 {
475
476 const Int_t npix = fAverageAreaNum[j];
477
478 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
479 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
480
481 hipix.FillHistAndArray(sumhiarea[j]/npix);
482 lopix.FillHistAndArray(sumloarea[j]/npix);
483
484 hipix.SetSaturated((Float_t)sathiarea[j]/npix);
485 lopix.SetSaturated((Float_t)satloarea[j]/npix);
486
487 hipix.FillAbsTime(timehiarea[j]/npix);
488 lopix.FillAbsTime(timeloarea[j]/npix);
489
490 }
491
492 for (UInt_t j=0; j<nsectors; j++)
493 {
494
495 const Int_t npix = fAverageSectorNum[j];
496
497 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
498 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
499
500 hipix.FillHistAndArray(sumhisector[j]/npix);
501 lopix.FillHistAndArray(sumlosector[j]/npix);
502
503 hipix.SetSaturated((Float_t)sathisector[j]/npix);
504 lopix.SetSaturated((Float_t)satlosector[j]/npix);
505
506 hipix.FillAbsTime(timehisector[j]/npix);
507 lopix.FillAbsTime(timelosector[j]/npix);
508
509 }
510
511 return kTRUE;
512}
513
514// --------------------------------------------------------------------------
515//
516// For all TObjArray's (including the averaged ones), the following steps are performed:
517//
518// 1) Returns if the pixel is excluded.
519// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
520// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
521// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
522// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
523// otherwise the Hi-Gain ones.
524// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
525// with the flags:
526// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
527// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
528// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
529// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
530//
531Bool_t MHCalibrationChargeCam::FinalizeHists()
532{
533
534 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
535 {
536
537 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
538 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
539
540 if (histhi.IsExcluded())
541 continue;
542
543 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
544 {
545 pix.SetHiGainSaturation();
546 histhi.CreateFourierSpectrum();
547 continue;
548 }
549
550 pix.SetAbsTimeMean ( histhi.GetAbsTimeMean());
551 pix.SetAbsTimeRms ( histhi.GetAbsTimeRms() );
552 }
553
554 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
555 {
556
557 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
558 MBadPixelsPix &bad = (*fBadPixels)[i];
559
560 if (histlo.IsExcluded())
561 continue;
562
563 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
564 {
565 *fLog << warn << "Saturated Lo Gain histogram in pixel: " << i << endl;
566 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
567 histlo.CreateFourierSpectrum();
568 continue;
569 }
570
571 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
572
573 if (pix.IsHiGainSaturation())
574 {
575 pix.SetAbsTimeMean ( histlo.GetAbsTimeMean());
576 pix.SetAbsTimeRms ( histlo.GetAbsTimeRms() );
577 }
578 }
579
580 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
581 {
582
583 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
584 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
585
586 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
587 {
588 pix.SetHiGainSaturation();
589 histhi.CreateFourierSpectrum();
590 continue;
591 }
592
593 pix.SetAbsTimeMean ( histhi.GetAbsTimeMean());
594 pix.SetAbsTimeRms ( histhi.GetAbsTimeRms() );
595 }
596
597 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
598 {
599
600 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
601 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
602
603 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
604 {
605 *fLog << warn << "Saturated Lo Gain histogram in area idx: " << j << endl;
606 histlo.CreateFourierSpectrum();
607 continue;
608 }
609
610 if (pix.IsHiGainSaturation())
611 {
612 pix.SetAbsTimeMean ( histlo.GetAbsTimeMean());
613 pix.SetAbsTimeRms ( histlo.GetAbsTimeRms() );
614 }
615 }
616
617 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
618 {
619
620 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
621 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
622
623 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
624 {
625 pix.SetHiGainSaturation();
626 histhi.CreateFourierSpectrum();
627 continue;
628 }
629
630 pix.SetAbsTimeMean ( histhi.GetAbsTimeMean());
631 pix.SetAbsTimeRms ( histhi.GetAbsTimeRms() );
632 }
633
634 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
635 {
636
637 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
638 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
639 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
640
641 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
642 {
643 *fLog << warn << "Saturated Lo Gain histogram in sector: " << j << endl;
644 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
645 histlo.CreateFourierSpectrum();
646 continue;
647 }
648
649 if (pix.IsHiGainSaturation())
650 {
651 pix.SetAbsTimeMean ( histlo.GetAbsTimeMean());
652 pix.SetAbsTimeRms ( histlo.GetAbsTimeRms() );
653 }
654 }
655
656 //
657 // Perform the fitting for the High Gain (done in MHCalibrationCam)
658 //
659 FitHiGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
660 MBadPixelsPix::kHiGainNotFitted,
661 MBadPixelsPix::kHiGainOscillating);
662 //
663 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
664 //
665 FitLoGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
666 MBadPixelsPix::kLoGainNotFitted,
667 MBadPixelsPix::kLoGainOscillating);
668
669 return kTRUE;
670}
671
672// --------------------------------------------------------------------------
673//
674// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
675// - MBadPixelsPix::kLoGainSaturation
676//
677// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
678// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
679// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
680// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
681// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
682//
683void MHCalibrationChargeCam::FinalizeBadPixels()
684{
685
686 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
687 {
688
689 MBadPixelsPix &bad = (*fBadPixels)[i];
690 MCalibrationPix &pix = (*fCam)[i];
691
692 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainNotFitted ))
693 if (!pix.IsHiGainSaturation())
694 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
695
696 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOscillating ))
697 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
698
699 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainNotFitted ))
700 if (pix.IsHiGainSaturation())
701 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
702
703 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOscillating ))
704 if (pix.IsHiGainSaturation())
705 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
706
707 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainSaturation ))
708 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
709 }
710}
711
712// --------------------------------------------------------------------------
713//
714// Dummy, needed by MCamEvent
715//
716Bool_t MHCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
717{
718 return kTRUE;
719}
720
721// --------------------------------------------------------------------------
722//
723// Calls MHGausEvents::DrawClone() for pixel idx
724//
725void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
726{
727 (*this)[idx].DrawClone();
728}
729
Note: See TracBrowser for help on using the repository browser.