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

Last change on this file since 3839 was 3722, checked in by gaug, 21 years ago
*** empty log message ***
File size: 28.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// 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;
146const Float_t MHCalibrationChargeCam::fgTimeLowerLimit = 1.;
147const Float_t MHCalibrationChargeCam::fgTimeUpperLimit = 2.;
148// --------------------------------------------------------------------------
149//
150// Default Constructor.
151//
152// Sets:
153// - all pointers to NULL
154//
155// Initializes:
156// - fNumHiGainSaturationLimit to fgNumHiGainSaturationLimit
157// - fNumLoGainSaturationLimit to fgNumLoGainSaturationLimit
158// - fTimeLowerLimit to fgTimeLowerLimit
159// - fTimeUpperLimit to fgTimeUpperLimit
160//
161MHCalibrationChargeCam::MHCalibrationChargeCam(const char *name, const char *title)
162 : fRawEvt(NULL)
163{
164 fName = name ? name : "MHCalibrationChargeCam";
165 fTitle = title ? title : "Class to fill the calibration histograms ";
166
167 SetNumHiGainSaturationLimit();
168 SetNumLoGainSaturationLimit();
169 SetTimeLowerLimit();
170 SetTimeUpperLimit();
171}
172
173// --------------------------------------------------------------------------
174//
175// Gets the pointers to:
176// - MRawEvtData
177//
178Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
179{
180
181 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
182 if (!fRawEvt)
183 {
184 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
185 return kFALSE;
186 }
187
188 return kTRUE;
189}
190
191// --------------------------------------------------------------------------
192//
193// Gets or creates the pointers to:
194// - MExtractedSignalCam
195// - MCalibrationChargeCam
196// - MBadPixelsCam
197//
198// Initializes the number of used FADC slices from MExtractedSignalCam
199// into MCalibrationChargeCam and test for changes in that variable
200//
201// Initializes, if empty to MGeomCam::GetNumPixels():
202// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
203//
204// Initializes, if empty to MGeomCam::GetNumAreas() for:
205// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
206//
207// Initializes, if empty to MGeomCam::GetNumSectors() for:
208// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
209//
210// Calls MHCalibrationCam::InitHists() for every entry in:
211// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
212// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
213// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
214//
215// Sets Titles and Names for the Charge Histograms:
216// - MHCalibrationCam::fAverageHiGainAreas
217// - MHCalibrationCam::fAverageHiGainSectors
218//
219// Sets number of bins to MHCalibrationCam::fAverageNbins for:
220// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
221// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
222//
223Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
224{
225
226 MExtractedSignalCam *signal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
227 if (!signal)
228 {
229 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
230 return kFALSE;
231 }
232
233 fCam = (MCalibrationCam*)pList->FindObject("MCalibrationChargeCam");
234 if (!fCam)
235 {
236 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationChargeCam"));
237 if (!fCam)
238 {
239 gLog << err << "Cannot find nor create MCalibrationChargeCam ... abort." << endl;
240 return kFALSE;
241 }
242 else
243 fCam->Init(*fGeom);
244 }
245
246 fFirstHiGain = signal->GetFirstUsedSliceHiGain();
247 fLastHiGain = signal->GetLastUsedSliceHiGain();
248 fFirstLoGain = signal->GetFirstUsedSliceLoGain();
249 fLastLoGain = signal->GetLastUsedSliceLoGain();
250
251 const Float_t numhigain = signal->GetNumUsedHiGainFADCSlices();
252 const Float_t numlogain = signal->GetNumUsedLoGainFADCSlices();
253
254 if (fCam->GetNumHiGainFADCSlices() == 0.)
255 fCam->SetNumHiGainFADCSlices ( numhigain );
256 else if (fCam->GetNumHiGainFADCSlices() != numhigain)
257 {
258 *fLog << err << GetDescriptor()
259 << ": Number of High Gain FADC extraction slices has changed, abort..." << endl;
260 return kFALSE;
261 }
262
263 if (fCam->GetNumLoGainFADCSlices() == 0.)
264 fCam->SetNumLoGainFADCSlices ( numlogain );
265 else if (fCam->GetNumLoGainFADCSlices() != numlogain)
266 {
267 *fLog << err << GetDescriptor()
268 << ": Number of Low Gain FADC extraction slices has changes, abort..." << endl;
269 return kFALSE;
270 }
271
272 const Int_t npixels = fGeom->GetNumPixels();
273 const Int_t nsectors = fGeom->GetNumSectors();
274 const Int_t nareas = fGeom->GetNumAreas();
275
276 if (fHiGainArray->GetEntries()==0)
277 {
278 fHiGainArray->Expand(npixels);
279 for (Int_t i=0; i<npixels; i++)
280 {
281 (*fHiGainArray)[i] = new MHCalibrationChargeHiGainPix;
282 InitHists((*this)[i],(*fBadPixels)[i],i);
283 }
284 }
285
286 if (fLoGainArray->GetEntries()==0)
287 {
288 fLoGainArray->Expand(npixels);
289
290 for (Int_t i=0; i<npixels; i++)
291 {
292 (*fLoGainArray)[i] = new MHCalibrationChargeLoGainPix;
293 InitHists((*this)(i),(*fBadPixels)[i],i);
294 }
295
296 }
297
298 if (fAverageHiGainAreas->GetEntries()==0)
299 {
300 fAverageHiGainAreas->Expand(nareas);
301
302 for (Int_t j=0; j<nareas; j++)
303 {
304 (*fAverageHiGainAreas)[j] =
305 new MHCalibrationChargeHiGainPix("AverageHiGainArea",
306 "Average HiGain FADC sums area idx ");
307
308 InitHists(GetAverageHiGainArea(j),fCam->GetAverageBadArea(j),j);
309
310 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
311
312 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Area Idx ");
313 hist.SetNbins(fAverageNbins);
314 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Area Idx ");
315 }
316 }
317
318
319 if (fAverageLoGainAreas->GetEntries()==0)
320 {
321 fAverageLoGainAreas->Expand(nareas);
322
323 for (Int_t j=0; j<nareas; j++)
324 {
325 (*fAverageLoGainAreas)[j] =
326 new MHCalibrationChargeLoGainPix("AverageLoGainArea",
327 "Average LoGain FADC sums of pixel area idx ");
328
329 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
330
331 InitHists(hist,fCam->GetAverageBadArea(j),j);
332
333 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Area Idx ");
334 hist.SetNbins(fAverageNbins);
335 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Area Idx ");
336
337 }
338 }
339
340 if (fAverageHiGainSectors->GetEntries()==0)
341 {
342 fAverageHiGainSectors->Expand(nsectors);
343
344 for (Int_t j=0; j<nsectors; j++)
345 {
346 (*fAverageHiGainSectors)[j] =
347 new MHCalibrationChargeHiGainPix("AverageHiGainSector",
348 "Average HiGain FADC sums of pixel sector ");
349
350 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
351
352 InitHists(hist,fCam->GetAverageBadSector(j),j);
353
354 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Sector ");
355 hist.SetNbins(fAverageNbins);
356 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Sector ");
357
358 }
359 }
360
361 if (fAverageLoGainSectors->GetEntries()==0)
362 {
363 fAverageLoGainSectors->Expand(nsectors);
364
365 for (Int_t j=0; j<nsectors; j++)
366 {
367 (*fAverageLoGainSectors)[j] =
368 new MHCalibrationChargeLoGainPix("AverageLoGainSector",
369 "Average LoGain FADC sums of pixel sector ");
370
371 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
372
373 InitHists(hist,fCam->GetAverageBadSector(j),j);
374
375 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Sector ");
376 hist.SetNbins(fAverageNbins);
377 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Sector ");
378
379 }
380 }
381
382 return kTRUE;
383}
384
385
386// --------------------------------------------------------------------------
387//
388// Retrieves from MExtractedSignalCam:
389// - first used LoGain FADC slice
390//
391// Retrieves from MGeomCam:
392// - number of pixels
393// - number of pixel areas
394// - number of sectors
395//
396// For all TObjArray's (including the averaged ones), the following steps are performed:
397//
398// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
399// - MExtractedSignalPix::GetExtractedSignalHiGain();
400// - MExtractedSignalPix::GetExtractedSignalLoGain();
401//
402// 2) Set number of saturated slices (MHCalibrationChargePix::SetSaturated()) with:
403// - MExtractedSignalPix::GetNumHiGainSaturated();
404// - MExtractedSignalPix::GetNumLoGainSaturated();
405//
406// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
407// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
408// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
409//
410Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
411{
412
413 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
414 if (!signal)
415 {
416 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
417 return kFALSE;
418 }
419
420 const UInt_t npixels = fGeom->GetNumPixels();
421 const UInt_t nareas = fGeom->GetNumAreas();
422 const UInt_t nsectors = fGeom->GetNumSectors();
423 const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
424
425 Float_t sumhiarea [nareas], sumloarea [nareas], timehiarea [nareas], timeloarea [nareas];
426 Float_t sumhisector[nsectors], sumlosector[nsectors], timehisector[nsectors], timelosector[nsectors];
427 Int_t sathiarea [nareas], satloarea [nareas];
428 Int_t sathisector[nsectors], satlosector[nsectors];
429
430 memset(sumhiarea, 0, nareas * sizeof(Float_t));
431 memset(sumloarea, 0, nareas * sizeof(Float_t));
432 memset(timehiarea, 0, nareas * sizeof(Float_t));
433 memset(timeloarea, 0, nareas * sizeof(Float_t));
434 memset(sathiarea, 0, nareas * sizeof(Int_t ));
435 memset(satloarea, 0, nareas * sizeof(Int_t ));
436 memset(sumhisector, 0, nsectors*sizeof(Float_t));
437 memset(sumlosector, 0, nsectors*sizeof(Float_t));
438 memset(timehisector,0, nsectors*sizeof(Float_t));
439 memset(timelosector,0, nsectors*sizeof(Float_t));
440 memset(sathisector, 0, nsectors*sizeof(Int_t ));
441 memset(satlosector, 0, nsectors*sizeof(Int_t ));
442
443 for (UInt_t i=0; i<npixels; i++)
444 {
445
446 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
447 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
448
449 if (histhi.IsExcluded())
450 continue;
451
452 const MExtractedSignalPix &pix = (*signal)[i];
453
454 const Float_t sumhi = pix.GetExtractedSignalHiGain();
455 const Float_t sumlo = pix.GetExtractedSignalLoGain();
456
457 histhi.FillHistAndArray(sumhi);
458 histlo.FillHistAndArray(sumlo);
459
460 const Int_t sathi = (Int_t)pix.GetNumHiGainSaturated();
461 const Int_t satlo = (Int_t)pix.GetNumLoGainSaturated();
462
463 histhi.SetSaturated(sathi);
464 histlo.SetSaturated(satlo);
465
466 const Int_t aidx = (*fGeom)[i].GetAidx();
467 const Int_t sector = (*fGeom)[i].GetSector();
468
469 sumhiarea[aidx] += sumhi;
470 sumloarea[aidx] += sumlo;
471 sathiarea[aidx] += sathi;
472 satloarea[aidx] += satlo;
473
474 sumhisector[sector] += sumhi;
475 sumlosector[sector] += sumlo;
476 sathisector[sector] += sathi;
477 satlosector[sector] += satlo;
478 }
479
480 MRawEvtPixelIter pixel(fRawEvt);
481 while (pixel.Next())
482 {
483
484 const UInt_t pixid = pixel.GetPixelId();
485
486 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
487 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
488
489 if (histhi.IsExcluded())
490 continue;
491
492 const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
493 const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
494
495 histhi.FillAbsTime(timehi);
496 histlo.FillAbsTime(timelo);
497
498 const Int_t aidx = (*fGeom)[pixid].GetAidx();
499 const Int_t sector = (*fGeom)[pixid].GetSector();
500
501 timehiarea[aidx] += timehi;
502 timeloarea[aidx] += timelo;
503
504 timehisector[sector] += timehi;
505 timelosector[sector] += timelo;
506 }
507
508 for (UInt_t j=0; j<nareas; j++)
509 {
510
511 const Int_t npix = fAverageAreaNum[j];
512
513 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
514 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
515
516 hipix.FillHistAndArray(sumhiarea[j]/npix);
517 lopix.FillHistAndArray(sumloarea[j]/npix);
518
519 hipix.SetSaturated((Float_t)sathiarea[j]/npix);
520 lopix.SetSaturated((Float_t)satloarea[j]/npix);
521
522 hipix.FillAbsTime(timehiarea[j]/npix);
523 lopix.FillAbsTime(timeloarea[j]/npix);
524
525 }
526
527 for (UInt_t j=0; j<nsectors; j++)
528 {
529
530 const Int_t npix = fAverageSectorNum[j];
531
532 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
533 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
534
535 hipix.FillHistAndArray(sumhisector[j]/npix);
536 lopix.FillHistAndArray(sumlosector[j]/npix);
537
538 hipix.SetSaturated((Float_t)sathisector[j]/npix);
539 lopix.SetSaturated((Float_t)satlosector[j]/npix);
540
541 hipix.FillAbsTime(timehisector[j]/npix);
542 lopix.FillAbsTime(timelosector[j]/npix);
543
544 }
545
546 return kTRUE;
547}
548
549// --------------------------------------------------------------------------
550//
551// For all TObjArray's (including the averaged ones), the following steps are performed:
552//
553// 1) Returns if the pixel is excluded.
554// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
555// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
556// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
557// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
558// otherwise the Hi-Gain ones.
559// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
560// with the flags:
561// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
562// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
563// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
564// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
565//
566Bool_t MHCalibrationChargeCam::FinalizeHists()
567{
568
569 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
570 {
571
572 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
573 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
574 MBadPixelsPix &bad = (*fBadPixels)[i];
575
576 if (histhi.IsExcluded())
577 continue;
578
579 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
580 {
581 pix.SetHiGainSaturation();
582 histhi.CreateFourierSpectrum();
583 continue;
584 }
585
586 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
587 }
588
589 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
590 {
591
592 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
593 MBadPixelsPix &bad = (*fBadPixels)[i];
594
595 if (histlo.IsExcluded())
596 continue;
597
598 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
599 {
600 *fLog << warn << "Saturated Lo Gain histogram in pixel: " << i << endl;
601 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
602 histlo.CreateFourierSpectrum();
603 continue;
604 }
605
606 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
607
608 if (pix.IsHiGainSaturation())
609 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
610 }
611
612 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
613 {
614
615 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
616 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
617 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
618
619 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
620 {
621 pix.SetHiGainSaturation();
622 histhi.CreateFourierSpectrum();
623 continue;
624 }
625
626 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
627 }
628
629 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
630 {
631
632 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
633 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
634 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
635
636 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
637 {
638 *fLog << warn << "Saturated Lo Gain histogram in area idx: " << j << endl;
639 histlo.CreateFourierSpectrum();
640 continue;
641 }
642
643 if (pix.IsHiGainSaturation())
644 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
645 }
646
647 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
648 {
649
650 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
651 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
652 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
653
654 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
655 {
656 pix.SetHiGainSaturation();
657 histhi.CreateFourierSpectrum();
658 continue;
659 }
660
661 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
662 }
663
664 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
665 {
666
667 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
668 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
669 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
670
671 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
672 {
673 *fLog << warn << "Saturated Lo Gain histogram in sector: " << j << endl;
674 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
675 histlo.CreateFourierSpectrum();
676 continue;
677 }
678
679 if (pix.IsHiGainSaturation())
680 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
681 }
682
683 //
684 // Perform the fitting for the High Gain (done in MHCalibrationCam)
685 //
686 FitHiGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
687 MBadPixelsPix::kHiGainNotFitted,
688 MBadPixelsPix::kHiGainOscillating);
689 //
690 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
691 //
692 FitLoGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
693 MBadPixelsPix::kLoGainNotFitted,
694 MBadPixelsPix::kLoGainOscillating);
695
696 return kTRUE;
697}
698
699// --------------------------------------------------------------------------------
700//
701// Fill the absolute time results into MCalibrationChargePix
702//
703// Check absolute time validity:
704// - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
705// - Mean arrival time is at least fUpperLimit slices from the upper edge
706//
707void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
708 Byte_t first, Byte_t last)
709{
710
711 const Float_t mean = hist.GetAbsTimeMean();
712 const Float_t rms = hist.GetAbsTimeRms();
713
714 pix.SetAbsTimeMean ( mean );
715 pix.SetAbsTimeRms ( rms );
716
717 const Float_t lowerlimit = (Float_t)first + fTimeLowerLimit;
718 const Float_t upperlimit = (Float_t)last + fTimeUpperLimit;
719
720 if ( mean < lowerlimit)
721 {
722 *fLog << warn << GetDescriptor()
723 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," smaller than ",fTimeLowerLimit,
724 " FADC slices from lower edge in pixel ",hist.GetPixId()) << endl;
725 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin );
726 }
727
728 if ( mean > upperlimit )
729 {
730 *fLog << warn << GetDescriptor()
731 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," greater than ",fTimeUpperLimit,
732 " FADC slices from upper edge in pixel ",hist.GetPixId()) << endl;
733 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins );
734 }
735}
736
737// --------------------------------------------------------------------------
738//
739// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
740// - MBadPixelsPix::kLoGainSaturation
741//
742// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
743// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
744// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
745// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
746// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
747//
748void MHCalibrationChargeCam::FinalizeBadPixels()
749{
750
751 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
752 {
753
754 MBadPixelsPix &bad = (*fBadPixels)[i];
755 MCalibrationPix &pix = (*fCam)[i];
756
757 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainNotFitted ))
758 if (!pix.IsHiGainSaturation())
759 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
760
761 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOscillating ))
762 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
763
764 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainNotFitted ))
765 if (pix.IsHiGainSaturation())
766 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
767
768 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOscillating ))
769 if (pix.IsHiGainSaturation())
770 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
771
772 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainSaturation ))
773 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
774 }
775}
776
777// --------------------------------------------------------------------------
778//
779// Dummy, needed by MCamEvent
780//
781Bool_t MHCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
782{
783 return kTRUE;
784}
785
786// --------------------------------------------------------------------------
787//
788// Calls MHGausEvents::DrawClone() for pixel idx
789//
790void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
791{
792 (*this)[idx].DrawClone();
793}
794
Note: See TracBrowser for help on using the repository browser.