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

Last change on this file since 4121 was 3939, checked in by gaug, 21 years ago
*** empty log message ***
File size: 28.2 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 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
309
310 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Area Idx ");
311 hist.SetNbins(fAverageNbins);
312 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Area Idx ");
313
314 InitHists(hist,fCam->GetAverageBadArea(j),j);
315
316 }
317 }
318
319
320 if (fAverageLoGainAreas->GetEntries()==0)
321 {
322 fAverageLoGainAreas->Expand(nareas);
323
324 for (Int_t j=0; j<nareas; j++)
325 {
326 (*fAverageLoGainAreas)[j] =
327 new MHCalibrationChargeLoGainPix("AverageLoGainArea",
328 "Average LoGain FADC sums of pixel area idx ");
329
330 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
331
332 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Area Idx ");
333 hist.SetNbins(fAverageNbins);
334 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Area Idx ");
335
336 InitHists(hist,fCam->GetAverageBadArea(j),j);
337
338 }
339 }
340
341 if (fAverageHiGainSectors->GetEntries()==0)
342 {
343 fAverageHiGainSectors->Expand(nsectors);
344
345 for (Int_t j=0; j<nsectors; j++)
346 {
347 (*fAverageHiGainSectors)[j] =
348 new MHCalibrationChargeHiGainPix("AverageHiGainSector",
349 "Average HiGain FADC sums of pixel sector ");
350
351 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
352
353 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Sector ");
354 hist.SetNbins(fAverageNbins);
355 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Sector ");
356
357 InitHists(hist,fCam->GetAverageBadSector(j),j);
358
359 }
360 }
361
362 if (fAverageLoGainSectors->GetEntries()==0)
363 {
364 fAverageLoGainSectors->Expand(nsectors);
365
366 for (Int_t j=0; j<nsectors; j++)
367 {
368 (*fAverageLoGainSectors)[j] =
369 new MHCalibrationChargeLoGainPix("AverageLoGainSector",
370 "Average LoGain FADC sums of pixel sector ");
371
372 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
373
374 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Sector ");
375 hist.SetNbins(fAverageNbins);
376 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Sector ");
377
378 InitHists(hist,fCam->GetAverageBadSector(j),j);
379
380 }
381 }
382
383 return kTRUE;
384}
385
386
387// --------------------------------------------------------------------------
388//
389// Retrieves from MExtractedSignalCam:
390// - first used LoGain FADC slice
391//
392// Retrieves from MGeomCam:
393// - number of pixels
394// - number of pixel areas
395// - number of sectors
396//
397// For all TObjArray's (including the averaged ones), the following steps are performed:
398//
399// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
400// - MExtractedSignalPix::GetExtractedSignalHiGain();
401// - MExtractedSignalPix::GetExtractedSignalLoGain();
402//
403// 2) Set number of saturated slices (MHCalibrationChargePix::SetSaturated()) with:
404// - MExtractedSignalPix::GetNumHiGainSaturated();
405// - MExtractedSignalPix::GetNumLoGainSaturated();
406//
407// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
408// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
409// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
410//
411Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
412{
413
414 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
415 if (!signal)
416 {
417 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
418 return kFALSE;
419 }
420
421 const UInt_t npixels = fGeom->GetNumPixels();
422 const UInt_t nareas = fGeom->GetNumAreas();
423 const UInt_t nsectors = fGeom->GetNumSectors();
424 const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
425
426 Float_t sumhiarea [nareas], sumloarea [nareas], timehiarea [nareas], timeloarea [nareas];
427 Float_t sumhisector[nsectors], sumlosector[nsectors], timehisector[nsectors], timelosector[nsectors];
428 Int_t sathiarea [nareas], satloarea [nareas];
429 Int_t sathisector[nsectors], satlosector[nsectors];
430
431 memset(sumhiarea, 0, nareas * sizeof(Float_t));
432 memset(sumloarea, 0, nareas * sizeof(Float_t));
433 memset(timehiarea, 0, nareas * sizeof(Float_t));
434 memset(timeloarea, 0, nareas * sizeof(Float_t));
435 memset(sathiarea, 0, nareas * sizeof(Int_t ));
436 memset(satloarea, 0, nareas * sizeof(Int_t ));
437 memset(sumhisector, 0, nsectors*sizeof(Float_t));
438 memset(sumlosector, 0, nsectors*sizeof(Float_t));
439 memset(timehisector,0, nsectors*sizeof(Float_t));
440 memset(timelosector,0, nsectors*sizeof(Float_t));
441 memset(sathisector, 0, nsectors*sizeof(Int_t ));
442 memset(satlosector, 0, nsectors*sizeof(Int_t ));
443
444 for (UInt_t i=0; i<npixels; i++)
445 {
446
447 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
448 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
449
450 if (histhi.IsExcluded())
451 continue;
452
453 const MExtractedSignalPix &pix = (*signal)[i];
454
455 const Float_t sumhi = pix.GetExtractedSignalHiGain();
456 const Float_t sumlo = pix.GetExtractedSignalLoGain();
457
458 histhi.FillHistAndArray(sumhi);
459 histlo.FillHistAndArray(sumlo);
460
461 const Int_t sathi = (Int_t)pix.GetNumHiGainSaturated();
462 const Int_t satlo = (Int_t)pix.GetNumLoGainSaturated();
463
464 histhi.SetSaturated(sathi);
465 histlo.SetSaturated(satlo);
466
467 const Int_t aidx = (*fGeom)[i].GetAidx();
468 const Int_t sector = (*fGeom)[i].GetSector();
469
470 sumhiarea[aidx] += sumhi;
471 sumloarea[aidx] += sumlo;
472 sathiarea[aidx] += sathi;
473 satloarea[aidx] += satlo;
474
475 sumhisector[sector] += sumhi;
476 sumlosector[sector] += sumlo;
477 sathisector[sector] += sathi;
478 satlosector[sector] += satlo;
479 }
480
481 MRawEvtPixelIter pixel(fRawEvt);
482 while (pixel.Next())
483 {
484
485 const UInt_t pixid = pixel.GetPixelId();
486
487 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
488 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
489
490 if (histhi.IsExcluded())
491 continue;
492
493 const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
494 const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
495
496 histhi.FillAbsTime(timehi);
497 histlo.FillAbsTime(timelo);
498
499 const Int_t aidx = (*fGeom)[pixid].GetAidx();
500 const Int_t sector = (*fGeom)[pixid].GetSector();
501
502 timehiarea[aidx] += timehi;
503 timeloarea[aidx] += timelo;
504
505 timehisector[sector] += timehi;
506 timelosector[sector] += timelo;
507 }
508
509 for (UInt_t j=0; j<nareas; j++)
510 {
511
512 const Int_t npix = fAverageAreaNum[j];
513
514 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
515 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
516
517 hipix.FillHistAndArray(sumhiarea[j]/npix);
518 lopix.FillHistAndArray(sumloarea[j]/npix);
519
520 hipix.SetSaturated((Float_t)sathiarea[j]/npix);
521 lopix.SetSaturated((Float_t)satloarea[j]/npix);
522
523 hipix.FillAbsTime(timehiarea[j]/npix);
524 lopix.FillAbsTime(timeloarea[j]/npix);
525
526 }
527
528 for (UInt_t j=0; j<nsectors; j++)
529 {
530
531 const Int_t npix = fAverageSectorNum[j];
532
533 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
534 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
535
536 hipix.FillHistAndArray(sumhisector[j]/npix);
537 lopix.FillHistAndArray(sumlosector[j]/npix);
538
539 hipix.SetSaturated((Float_t)sathisector[j]/npix);
540 lopix.SetSaturated((Float_t)satlosector[j]/npix);
541
542 hipix.FillAbsTime(timehisector[j]/npix);
543 lopix.FillAbsTime(timelosector[j]/npix);
544
545 }
546
547 return kTRUE;
548}
549
550// --------------------------------------------------------------------------
551//
552// For all TObjArray's (including the averaged ones), the following steps are performed:
553//
554// 1) Returns if the pixel is excluded.
555// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
556// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
557// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
558// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
559// otherwise the Hi-Gain ones.
560// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
561// with the flags:
562// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
563// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
564// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
565// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
566//
567Bool_t MHCalibrationChargeCam::FinalizeHists()
568{
569
570 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
571 {
572
573 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
574 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
575 MBadPixelsPix &bad = (*fBadPixels)[i];
576
577 if (histhi.IsExcluded())
578 continue;
579
580 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
581 {
582 pix.SetHiGainSaturation();
583 histhi.CreateFourierSpectrum();
584 continue;
585 }
586
587 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
588 }
589
590 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
591 {
592
593 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
594 MBadPixelsPix &bad = (*fBadPixels)[i];
595
596 if (histlo.IsExcluded())
597 continue;
598
599 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
600 {
601 *fLog << warn << "Saturated Lo Gain histogram in pixel: " << i << endl;
602 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
603 histlo.CreateFourierSpectrum();
604 continue;
605 }
606
607 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
608
609 if (pix.IsHiGainSaturation())
610 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
611 }
612
613 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
614 {
615
616 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
617 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
618 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
619
620 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
621 {
622 pix.SetHiGainSaturation();
623 histhi.CreateFourierSpectrum();
624 continue;
625 }
626
627 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
628 }
629
630 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
631 {
632
633 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
634 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
635 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
636
637 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
638 {
639 *fLog << warn << "Saturated Lo Gain histogram in area idx: " << j << endl;
640 histlo.CreateFourierSpectrum();
641 continue;
642 }
643
644 if (pix.IsHiGainSaturation())
645 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
646 }
647
648 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
649 {
650
651 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
652 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
653 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
654
655 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
656 {
657 pix.SetHiGainSaturation();
658 histhi.CreateFourierSpectrum();
659 continue;
660 }
661
662 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
663 }
664
665 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
666 {
667
668 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
669 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
670 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
671
672 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
673 {
674 *fLog << warn << "Saturated Lo Gain histogram in sector: " << j << endl;
675 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
676 histlo.CreateFourierSpectrum();
677 continue;
678 }
679
680 if (pix.IsHiGainSaturation())
681 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
682 }
683
684 //
685 // Perform the fitting for the High Gain (done in MHCalibrationCam)
686 //
687 FitHiGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
688 MBadPixelsPix::kHiGainNotFitted,
689 MBadPixelsPix::kHiGainOscillating);
690 //
691 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
692 //
693 FitLoGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
694 MBadPixelsPix::kLoGainNotFitted,
695 MBadPixelsPix::kLoGainOscillating);
696
697 return kTRUE;
698}
699
700// --------------------------------------------------------------------------------
701//
702// Fill the absolute time results into MCalibrationChargePix
703//
704// Check absolute time validity:
705// - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
706// - Mean arrival time is at least fUpperLimit slices from the upper edge
707//
708void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
709 Byte_t first, Byte_t last)
710{
711
712 const Float_t mean = hist.GetAbsTimeMean();
713 const Float_t rms = hist.GetAbsTimeRms();
714
715 pix.SetAbsTimeMean ( mean );
716 pix.SetAbsTimeRms ( rms );
717
718 const Float_t lowerlimit = (Float_t)first + fTimeLowerLimit;
719 const Float_t upperlimit = (Float_t)last + fTimeUpperLimit;
720
721 if ( mean < lowerlimit)
722 {
723 *fLog << warn << GetDescriptor()
724 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," smaller than ",fTimeLowerLimit,
725 " FADC slices from lower edge in pixel ",hist.GetPixId()) << endl;
726 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin );
727 }
728
729 if ( mean > upperlimit )
730 {
731 *fLog << warn << GetDescriptor()
732 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," greater than ",fTimeUpperLimit,
733 " FADC slices from upper edge in pixel ",hist.GetPixId()) << endl;
734 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins );
735 }
736}
737
738// --------------------------------------------------------------------------
739//
740// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
741// - MBadPixelsPix::kLoGainSaturation
742//
743// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
744// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
745// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
746// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
747// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
748//
749void MHCalibrationChargeCam::FinalizeBadPixels()
750{
751
752 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
753 {
754
755 MBadPixelsPix &bad = (*fBadPixels)[i];
756 MCalibrationPix &pix = (*fCam)[i];
757
758 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainNotFitted ))
759 if (!pix.IsHiGainSaturation())
760 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
761
762 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOscillating ))
763 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
764
765 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainNotFitted ))
766 if (pix.IsHiGainSaturation())
767 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
768
769 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOscillating ))
770 if (pix.IsHiGainSaturation())
771 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
772
773 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainSaturation ))
774 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
775 }
776}
777
778// --------------------------------------------------------------------------
779//
780// Dummy, needed by MCamEvent
781//
782Bool_t MHCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
783{
784 return kTRUE;
785}
786
787// --------------------------------------------------------------------------
788//
789// Calls MHGausEvents::DrawClone() for pixel idx
790//
791void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
792{
793 (*this)[idx].DrawClone();
794}
795
Note: See TracBrowser for help on using the repository browser.