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

Last change on this file since 4537 was 4537, checked in by gaug, 20 years ago
*** empty log message ***
File size: 36.4 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
140#include <TPad.h>
141#include <TVirtualPad.h>
142#include <TCanvas.h>
143#include <TStyle.h>
144#include <TF1.h>
145#include <TH2D.h>
146#include <TLine.h>
147#include <TLatex.h>
148#include <TLegend.h>
149
150ClassImp(MHCalibrationChargeCam);
151
152using namespace std;
153
154const Float_t MHCalibrationChargeCam::fgNumHiGainSaturationLimit = 0.01;
155const Float_t MHCalibrationChargeCam::fgNumLoGainSaturationLimit = 0.005;
156const Float_t MHCalibrationChargeCam::fgTimeLowerLimit = 1.;
157const Float_t MHCalibrationChargeCam::fgTimeUpperLimit = 2.;
158// 1Led Green, 1 LED blue, 5 LEDs blue, 10 LEDs blue, 10 LEDs UV, CT1, 5Leds Green
159const Float_t MHCalibrationChargeCam::gkHiGainInnerRefLines[7] = { 245., 323. , 1065., 1467., 180., 211. , 533.5};
160const Float_t MHCalibrationChargeCam::gkHiGainOuterRefLines[7] = { 217., 307.5, 932. , 1405., 167., 183.5, 405.5};
161const Float_t MHCalibrationChargeCam::gkLoGainInnerRefLines[7] = { 20.8, 28.0 , 121. , 200.2, 16.5, 13.5 , 41.7 };
162const Float_t MHCalibrationChargeCam::gkLoGainOuterRefLines[7] = { 18.9, 26.0 , 108.3, 198. , 14.0, 11. , 42. };
163// --------------------------------------------------------------------------
164//
165// Default Constructor.
166//
167// Sets:
168// - all pointers to NULL
169//
170// Initializes:
171// - fNumHiGainSaturationLimit to fgNumHiGainSaturationLimit
172// - fNumLoGainSaturationLimit to fgNumLoGainSaturationLimit
173// - fTimeLowerLimit to fgTimeLowerLimit
174// - fTimeUpperLimit to fgTimeUpperLimit
175//
176MHCalibrationChargeCam::MHCalibrationChargeCam(const char *name, const char *title)
177 : fRawEvt(NULL)
178{
179 fName = name ? name : "MHCalibrationChargeCam";
180 fTitle = title ? title : "Class to fill the calibration histograms ";
181
182 SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
183 SetNumLoGainSaturationLimit(fgNumLoGainSaturationLimit);
184 SetTimeLowerLimit();
185 SetTimeUpperLimit();
186
187
188}
189
190// --------------------------------------------------------------------------
191//
192// Gets the pointers to:
193// - MRawEvtData
194//
195Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
196{
197
198 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
199 if (!fRawEvt)
200 {
201 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
202 return kFALSE;
203 }
204
205 return kTRUE;
206}
207
208// --------------------------------------------------------------------------
209//
210// Gets or creates the pointers to:
211// - MExtractedSignalCam
212// - MCalibrationChargeCam
213// - MBadPixelsCam
214//
215// Initializes the number of used FADC slices from MExtractedSignalCam
216// into MCalibrationChargeCam and test for changes in that variable
217//
218// Initializes, if empty to MGeomCam::GetNumPixels():
219// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
220//
221// Initializes, if empty to MGeomCam::GetNumAreas() for:
222// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
223//
224// Initializes, if empty to MGeomCam::GetNumSectors() for:
225// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
226//
227// Calls MHCalibrationCam::InitHists() for every entry in:
228// - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
229// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
230// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
231//
232// Sets Titles and Names for the Charge Histograms:
233// - MHCalibrationCam::fAverageHiGainAreas
234// - MHCalibrationCam::fAverageHiGainSectors
235//
236// Sets number of bins to MHCalibrationCam::fAverageNbins for:
237// - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
238// - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
239//
240Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
241{
242
243 MExtractedSignalCam *signal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
244 if (!signal)
245 {
246 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
247 return kFALSE;
248 }
249
250 fCam = (MCalibrationCam*)pList->FindObject("MCalibrationChargeCam");
251 if (!fCam)
252 {
253 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationChargeCam"));
254 if (!fCam)
255 {
256 gLog << err << "Cannot find nor create MCalibrationChargeCam ... abort." << endl;
257 return kFALSE;
258 }
259 else
260 fCam->Init(*fGeom);
261 }
262
263 fFirstHiGain = signal->GetFirstUsedSliceHiGain();
264 fLastHiGain = signal->GetLastUsedSliceHiGain();
265 fFirstLoGain = signal->GetFirstUsedSliceLoGain();
266 fLastLoGain = signal->GetLastUsedSliceLoGain();
267
268 const Float_t numhigain = signal->GetNumUsedHiGainFADCSlices();
269 const Float_t numlogain = signal->GetNumUsedLoGainFADCSlices();
270
271 if (fCam->GetNumHiGainFADCSlices() == 0.)
272 fCam->SetNumHiGainFADCSlices ( numhigain );
273 else if (fCam->GetNumHiGainFADCSlices() != numhigain)
274 {
275 *fLog << err << GetDescriptor()
276 << ": Number of High Gain FADC extraction slices has changed, abort..." << endl;
277 return kFALSE;
278 }
279
280 if (fCam->GetNumLoGainFADCSlices() == 0.)
281 fCam->SetNumLoGainFADCSlices ( numlogain );
282 else if (fCam->GetNumLoGainFADCSlices() != numlogain)
283 {
284 *fLog << err << GetDescriptor()
285 << ": Number of Low Gain FADC extraction slices has changes, abort..." << endl;
286 return kFALSE;
287 }
288
289 const Int_t npixels = fGeom->GetNumPixels();
290 const Int_t nsectors = fGeom->GetNumSectors();
291 const Int_t nareas = fGeom->GetNumAreas();
292
293 if (fHiGainArray->GetEntries()==0)
294 {
295 fHiGainArray->Expand(npixels);
296 for (Int_t i=0; i<npixels; i++)
297 {
298 (*fHiGainArray)[i] = new MHCalibrationChargeHiGainPix;
299 InitHists((*this)[i],(*fBadPixels)[i],i);
300 }
301 }
302
303 if (fLoGainArray->GetEntries()==0)
304 {
305 fLoGainArray->Expand(npixels);
306
307 for (Int_t i=0; i<npixels; i++)
308 {
309 (*fLoGainArray)[i] = new MHCalibrationChargeLoGainPix;
310 InitHists((*this)(i),(*fBadPixels)[i],i);
311 }
312
313 }
314
315 if (fAverageHiGainAreas->GetEntries()==0)
316 {
317 fAverageHiGainAreas->Expand(nareas);
318
319 for (Int_t j=0; j<nareas; j++)
320 {
321 (*fAverageHiGainAreas)[j] =
322 new MHCalibrationChargeHiGainPix("AverageHiGainArea",
323 "Average HiGain FADC sums area idx ");
324
325 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
326
327 hist.SetNbins(fAverageNbins);
328 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Area Idx ");
329
330 if (fGeom->InheritsFrom("MGeomCamMagic"))
331 {
332 hist.GetHGausHist()->SetTitle(Form("%s%s%s","Signal averaged on event-by-event basis ",
333 j==0 ? "Inner Pixels " : "Outer Pixels ","High Gain Runs: "));
334 hist.InitBins();
335 hist.SetEventFrequency(fPulserFrequency);
336 }
337 else
338 {
339 hist.GetHGausHist()->SetTitle("Signal averaged on event-by-event basis High Gain Area Idx ");
340 InitHists(hist,fCam->GetAverageBadArea(j),j);
341 }
342 }
343 }
344
345
346 if (fAverageLoGainAreas->GetEntries()==0)
347 {
348 fAverageLoGainAreas->Expand(nareas);
349
350 for (Int_t j=0; j<nareas; j++)
351 {
352 (*fAverageLoGainAreas)[j] =
353 new MHCalibrationChargeLoGainPix("AverageLoGainArea",
354 "Average LoGain FADC sums of pixel area idx ");
355
356 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
357
358 hist.SetNbins(fAverageNbins);
359 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Area Idx ");
360
361 if (fGeom->InheritsFrom("MGeomCamMagic"))
362 {
363 hist.GetHGausHist()->SetTitle(Form("%s%s%s","Signal averaged on event-by-event basis ",
364 j==0 ? "Inner Pixels " : "Outer Pixels ","High Gain Runs: "));
365 hist.InitBins();
366 hist.SetEventFrequency(fPulserFrequency);
367 }
368 else
369 {
370 hist.GetHGausHist()->SetTitle("Signal averaged on event-by-event basis High Gain Area Idx ");
371 InitHists(hist,fCam->GetAverageBadArea(j),j);
372 }
373 }
374 }
375
376 if (fAverageHiGainSectors->GetEntries()==0)
377 {
378 fAverageHiGainSectors->Expand(nsectors);
379
380 for (Int_t j=0; j<nsectors; j++)
381 {
382 (*fAverageHiGainSectors)[j] =
383 new MHCalibrationChargeHiGainPix("AverageHiGainSector",
384 "Average HiGain FADC sums of pixel sector ");
385
386 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
387
388 hist.GetHGausHist()->SetTitle("Summed FADC slices average HiGain Sector ");
389 hist.SetNbins(fAverageNbins);
390 hist.SetLast (2.*hist.GetLast());
391 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average HiGain Sector ");
392
393 InitHists(hist,fCam->GetAverageBadSector(j),j);
394
395 }
396 }
397
398 if (fAverageLoGainSectors->GetEntries()==0)
399 {
400 fAverageLoGainSectors->Expand(nsectors);
401
402 for (Int_t j=0; j<nsectors; j++)
403 {
404 (*fAverageLoGainSectors)[j] =
405 new MHCalibrationChargeLoGainPix("AverageLoGainSector",
406 "Average LoGain FADC sums of pixel sector ");
407
408 MHCalibrationChargePix &hist = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
409
410 hist.GetHGausHist()->SetTitle("Summed FADC slices average LoGain Sector ");
411 hist.SetNbins(fAverageNbins);
412 hist.GetHAbsTime()->SetTitle("Absolute Arrival Time average LoGain Sector ");
413
414 InitHists(hist,fCam->GetAverageBadSector(j),j);
415
416 }
417 }
418
419 return kTRUE;
420}
421
422
423// --------------------------------------------------------------------------
424//
425// Retrieves from MExtractedSignalCam:
426// - first used LoGain FADC slice
427//
428// Retrieves from MGeomCam:
429// - number of pixels
430// - number of pixel areas
431// - number of sectors
432//
433// For all TObjArray's (including the averaged ones), the following steps are performed:
434//
435// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
436// - MExtractedSignalPix::GetExtractedSignalHiGain();
437// - MExtractedSignalPix::GetExtractedSignalLoGain();
438//
439// 2) Set number of saturated slices (MHCalibrationChargePix::SetSaturated()) with:
440// - MExtractedSignalPix::GetNumHiGainSaturated();
441// - MExtractedSignalPix::GetNumLoGainSaturated();
442//
443// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
444// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
445// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
446//
447Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
448{
449
450 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
451 if (!signal)
452 {
453 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
454 return kFALSE;
455 }
456
457 const UInt_t npixels = fGeom->GetNumPixels();
458 const UInt_t nareas = fGeom->GetNumAreas();
459 const UInt_t nsectors = fGeom->GetNumSectors();
460 const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
461
462 Float_t sumhiarea [nareas], sumloarea [nareas], timehiarea [nareas], timeloarea [nareas];
463 Float_t sumhisector[nsectors], sumlosector[nsectors], timehisector[nsectors], timelosector[nsectors];
464 Int_t sathiarea [nareas], satloarea [nareas];
465 Int_t sathisector[nsectors], satlosector[nsectors];
466
467 memset(sumhiarea, 0, nareas * sizeof(Float_t));
468 memset(sumloarea, 0, nareas * sizeof(Float_t));
469 memset(timehiarea, 0, nareas * sizeof(Float_t));
470 memset(timeloarea, 0, nareas * sizeof(Float_t));
471 memset(sathiarea, 0, nareas * sizeof(Int_t ));
472 memset(satloarea, 0, nareas * sizeof(Int_t ));
473 memset(sumhisector, 0, nsectors*sizeof(Float_t));
474 memset(sumlosector, 0, nsectors*sizeof(Float_t));
475 memset(timehisector,0, nsectors*sizeof(Float_t));
476 memset(timelosector,0, nsectors*sizeof(Float_t));
477 memset(sathisector, 0, nsectors*sizeof(Int_t ));
478 memset(satlosector, 0, nsectors*sizeof(Int_t ));
479
480 for (UInt_t i=0; i<npixels; i++)
481 {
482
483 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
484 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
485
486 if (histhi.IsExcluded())
487 continue;
488
489 const MExtractedSignalPix &pix = (*signal)[i];
490
491 const Float_t sumhi = pix.GetExtractedSignalHiGain();
492 const Float_t sumlo = pix.GetExtractedSignalLoGain();
493
494 if (!histhi.FillHistAndArray(sumhi))
495 fHiGainOverFlow++;
496 if (!histlo.FillHistAndArray(sumlo))
497 fLoGainOverFlow++;
498
499 const Int_t sathi = (Int_t)pix.GetNumHiGainSaturated();
500 const Int_t satlo = (Int_t)pix.GetNumLoGainSaturated();
501
502 histhi.SetSaturated(sathi);
503 histlo.SetSaturated(satlo);
504
505 const Int_t aidx = (*fGeom)[i].GetAidx();
506 const Int_t sector = (*fGeom)[i].GetSector();
507
508 sumhiarea[aidx] += sumhi;
509 sumloarea[aidx] += sumlo;
510 sathiarea[aidx] += sathi;
511 satloarea[aidx] += satlo;
512
513 sumhisector[sector] += sumhi;
514 sumlosector[sector] += sumlo;
515 sathisector[sector] += sathi;
516 satlosector[sector] += satlo;
517 }
518
519 MRawEvtPixelIter pixel(fRawEvt);
520 while (pixel.Next())
521 {
522
523 const UInt_t pixid = pixel.GetPixelId();
524
525 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
526 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
527
528 if (histhi.IsExcluded())
529 continue;
530
531 const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
532 const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
533
534 histhi.FillAbsTime(timehi);
535 histlo.FillAbsTime(timelo);
536
537 const Int_t aidx = (*fGeom)[pixid].GetAidx();
538 const Int_t sector = (*fGeom)[pixid].GetSector();
539
540 timehiarea[aidx] += timehi;
541 timeloarea[aidx] += timelo;
542
543 timehisector[sector] += timehi;
544 timelosector[sector] += timelo;
545 }
546
547 for (UInt_t j=0; j<nareas; j++)
548 {
549
550 const Int_t npix = fAverageAreaNum[j];
551
552 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
553 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
554
555 hipix.FillHistAndArray(sumhiarea[j]/npix);
556 lopix.FillHistAndArray(sumloarea[j]/npix);
557
558 hipix.SetSaturated((Float_t)sathiarea[j]/npix);
559 lopix.SetSaturated((Float_t)satloarea[j]/npix);
560
561 hipix.FillAbsTime(timehiarea[j]/npix);
562 lopix.FillAbsTime(timeloarea[j]/npix);
563
564 }
565
566 for (UInt_t j=0; j<nsectors; j++)
567 {
568
569 const Int_t npix = fAverageSectorNum[j];
570
571 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
572 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
573
574 hipix.FillHistAndArray(sumhisector[j]/npix);
575 lopix.FillHistAndArray(sumlosector[j]/npix);
576
577 hipix.SetSaturated((Float_t)sathisector[j]/npix);
578 lopix.SetSaturated((Float_t)satlosector[j]/npix);
579
580 hipix.FillAbsTime(timehisector[j]/npix);
581 lopix.FillAbsTime(timelosector[j]/npix);
582
583 }
584
585 return kTRUE;
586}
587
588// --------------------------------------------------------------------------
589//
590// For all TObjArray's (including the averaged ones), the following steps are performed:
591//
592// 1) Returns if the pixel is excluded.
593// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
594// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
595// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
596// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
597// otherwise the Hi-Gain ones.
598// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
599// with the flags:
600// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
601// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
602// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
603// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
604//
605Bool_t MHCalibrationChargeCam::FinalizeHists()
606{
607
608 if (fHiGainOverFlow)
609 *fLog << warn << GetDescriptor()
610 << ": WARNING: Histogram Overflow has occurred " << fHiGainOverFlow << " in the High-Gain! " << endl;
611 if (fLoGainOverFlow)
612 *fLog << warn << GetDescriptor()
613 << ": WARNING: Histogram Overflow has occurred " << fLoGainOverFlow << " in the Low-Gain! " << endl;
614
615 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
616 {
617
618 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
619 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
620 MBadPixelsPix &bad = (*fBadPixels)[i];
621
622 if (histhi.IsExcluded())
623 continue;
624
625 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
626 {
627 pix.SetHiGainSaturation();
628 histhi.CreateFourierSpectrum();
629 continue;
630 }
631
632 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
633 }
634
635 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
636 {
637
638 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
639 MBadPixelsPix &bad = (*fBadPixels)[i];
640
641 if (histlo.IsExcluded())
642 continue;
643
644 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
645 {
646 *fLog << warn << "Saturated Lo Gain histogram in pixel: " << i << endl;
647 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
648 histlo.CreateFourierSpectrum();
649 continue;
650 }
651
652 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
653
654 if (pix.IsHiGainSaturation())
655 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
656 }
657
658 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
659 {
660
661 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
662 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
663 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
664
665 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
666 {
667 pix.SetHiGainSaturation();
668 histhi.CreateFourierSpectrum();
669 continue;
670 }
671
672 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
673 }
674
675 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
676 {
677
678 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
679 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
680 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
681
682 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
683 {
684 *fLog << warn << "Saturated Lo Gain histogram in area idx: " << j << endl;
685 histlo.CreateFourierSpectrum();
686 continue;
687 }
688
689 if (pix.IsHiGainSaturation())
690 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
691 }
692
693 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
694 {
695
696 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
697 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
698 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
699
700 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
701 {
702 pix.SetHiGainSaturation();
703 histhi.CreateFourierSpectrum();
704 continue;
705 }
706
707 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
708 }
709
710 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
711 {
712
713 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
714 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
715 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
716
717 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
718 {
719 *fLog << warn << "Saturated Lo Gain histogram in sector: " << j << endl;
720 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
721 histlo.CreateFourierSpectrum();
722 continue;
723 }
724
725 if (pix.IsHiGainSaturation())
726 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
727 }
728
729 //
730 // Perform the fitting for the High Gain (done in MHCalibrationCam)
731 //
732 FitHiGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
733 MBadPixelsPix::kHiGainNotFitted,
734 MBadPixelsPix::kHiGainOscillating);
735 //
736 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
737 //
738 FitLoGainArrays((MCalibrationCam&)(*fCam),(*fBadPixels),
739 MBadPixelsPix::kLoGainNotFitted,
740 MBadPixelsPix::kLoGainOscillating);
741
742 return kTRUE;
743}
744
745// --------------------------------------------------------------------------------
746//
747// Fill the absolute time results into MCalibrationChargePix
748//
749// Check absolute time validity:
750// - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
751// - Mean arrival time is at least fUpperLimit slices from the upper edge
752//
753void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
754 Byte_t first, Byte_t last)
755{
756
757 const Float_t mean = hist.GetAbsTimeMean();
758 const Float_t rms = hist.GetAbsTimeRms();
759
760 pix.SetAbsTimeMean ( mean );
761 pix.SetAbsTimeRms ( rms );
762
763 const Float_t lowerlimit = (Float_t)first + fTimeLowerLimit;
764 const Float_t upperlimit = (Float_t)last + fTimeUpperLimit;
765
766 if ( mean < lowerlimit)
767 {
768 *fLog << warn << GetDescriptor()
769 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," smaller than ",fTimeLowerLimit,
770 " FADC slices from lower edge in pixel ",hist.GetPixId()) << endl;
771 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin );
772 }
773
774 if ( mean > upperlimit )
775 {
776 *fLog << warn << GetDescriptor()
777 << Form("%s%3.1f%s%2.1f%s%4i",": Mean ArrivalTime: ",mean," greater than ",fTimeUpperLimit,
778 " FADC slices from upper edge in pixel ",hist.GetPixId()) << endl;
779 bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins );
780 }
781}
782
783// --------------------------------------------------------------------------
784//
785// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
786// - MBadPixelsPix::kLoGainSaturation
787//
788// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
789// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
790// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
791// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
792// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
793//
794void MHCalibrationChargeCam::FinalizeBadPixels()
795{
796
797 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
798 {
799
800 MBadPixelsPix &bad = (*fBadPixels)[i];
801 MCalibrationPix &pix = (*fCam)[i];
802
803 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainNotFitted ))
804 if (!pix.IsHiGainSaturation())
805 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
806
807 if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOscillating ))
808 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
809
810 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainNotFitted ))
811 if (pix.IsHiGainSaturation())
812 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
813
814 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOscillating ))
815 if (pix.IsHiGainSaturation())
816 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
817
818 if (bad.IsUncalibrated( MBadPixelsPix::kLoGainSaturation ))
819 bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
820 }
821}
822
823// --------------------------------------------------------------------------
824//
825// Dummy, needed by MCamEvent
826//
827Bool_t MHCalibrationChargeCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
828{
829 return kTRUE;
830}
831
832// --------------------------------------------------------------------------
833//
834// Calls MHGausEvents::DrawClone() for pixel idx
835//
836void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
837{
838 (*this)[idx].DrawClone();
839}
840
841
842// -----------------------------------------------------------------------------
843//
844// Default draw:
845//
846// Displays the averaged areas, both High Gain and Low Gain
847//
848// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
849//
850void MHCalibrationChargeCam::Draw(const Option_t *opt)
851{
852
853 const Int_t nareas = fAverageHiGainAreas->GetEntries();
854 if (nareas == 0)
855 return;
856
857 TString option(opt);
858 option.ToLower();
859
860 if (!option.Contains("datacheck"))
861 {
862 MHCalibrationCam::Draw(opt);
863 return;
864 }
865
866 //
867 // From here on , the datacheck - Draw
868 //
869 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
870 pad->SetBorderMode(0);
871
872 pad->Divide(2,nareas);
873
874 for (Int_t i=0; i<nareas;i++)
875 {
876 pad->cd(2*(i+1)-1);
877 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(i);
878
879 if (i==0)
880 DrawDataCheckPixel(hipix,gkHiGainInnerRefLines);
881 else
882 DrawDataCheckPixel(hipix,gkHiGainOuterRefLines);
883
884 pad->cd(2*(i+1));
885
886 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(i);
887
888 if (i==0)
889 DrawDataCheckPixel(lopix,gkLoGainInnerRefLines);
890 else
891 DrawDataCheckPixel(lopix,gkLoGainOuterRefLines);
892
893 }
894}
895
896
897// --------------------------------------------------------------------------
898//
899// Our own clone function is necessary since root 3.01/06 or Mars 0.4
900// I don't know the reason.
901//
902// Creates new MHCalibrationCam
903//
904TObject *MHCalibrationChargeCam::Clone(const char *) const
905{
906
907 const Int_t navhi = fAverageHiGainAreas->GetEntries();
908 const Int_t navlo = fAverageLoGainAreas->GetEntries();
909 const Int_t nsehi = fAverageHiGainSectors->GetEntries();
910 const Int_t nselo = fAverageLoGainSectors->GetEntries();
911
912 //
913 // FIXME, this might be done faster and more elegant, by direct copy.
914 //
915 MHCalibrationChargeCam *cam = new MHCalibrationChargeCam();
916
917 cam->fAverageHiGainAreas->Expand(navhi);
918 cam->fAverageLoGainAreas->Expand(navlo);
919 cam->fAverageHiGainSectors->Expand(nsehi);
920 cam->fAverageLoGainSectors->Expand(nselo);
921
922 cam->fAverageHiGainAreas->Expand(navhi);
923 cam->fAverageLoGainAreas->Expand(navlo);
924 cam->fAverageHiGainSectors->Expand(nsehi);
925 cam->fAverageLoGainSectors->Expand(nselo);
926
927 for (int i=0; i<navhi; i++)
928 (*cam->fAverageHiGainAreas)[i] = (*fAverageHiGainAreas)[i]->Clone();
929 for (int i=0; i<navlo; i++)
930 (*cam->fAverageLoGainAreas)[i] = (*fAverageLoGainAreas)[i]->Clone();
931 for (int i=0; i<nsehi; i++)
932 (*cam->fAverageHiGainSectors)[i] = (*fAverageHiGainSectors)[i]->Clone();
933 for (int i=0; i<nselo; i++)
934 (*cam->fAverageLoGainSectors)[i] = (*fAverageLoGainSectors)[i]->Clone();
935
936 cam->fAverageAreaNum = fAverageAreaNum;
937 cam->fAverageAreaSat = fAverageAreaSat;
938 cam->fAverageAreaSigma = fAverageAreaSigma;
939 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
940 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
941 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
942 cam->fAverageSectorNum = fAverageSectorNum;
943 cam->fRunNumbers = fRunNumbers;
944
945 cam->fPulserFrequency = fPulserFrequency;
946 cam->fAverageNbins = fAverageNbins;
947
948 return cam;
949
950}
951
952void MHCalibrationChargeCam::DrawDataCheckPixel(MHCalibrationChargePix &pix, const Float_t refline[])
953{
954
955 TVirtualPad *newpad = gPad;
956 newpad->Divide(1,2);
957 newpad->cd(1);
958
959 gPad->SetTicks();
960 if (!pix.IsEmpty())
961 gPad->SetLogy();
962
963 gStyle->SetOptStat(0);
964
965 TH1F *hist = pix.GetHGausHist();
966
967
968 TH2D *null = new TH2D("Null",hist->GetTitle(),100,pix.GetFirst(),pix.GetLast(),
969 100,0.,hist->GetEntries()/10.);
970
971 null->SetDirectory(NULL);
972 null->SetBit(kCanDelete);
973 null->GetXaxis()->SetTitle(hist->GetXaxis()->GetTitle());
974 null->GetYaxis()->SetTitle(hist->GetYaxis()->GetTitle());
975 null->GetXaxis()->CenterTitle();
976 null->GetYaxis()->CenterTitle();
977 null->Draw();
978 hist->Draw("same");
979
980 gStyle->SetOptFit();
981
982 if (pix.GetFGausFit())
983 {
984 pix.GetFGausFit()->SetLineColor(pix.IsGausFitOK() ? kGreen : kRed);
985 pix.GetFGausFit()->Draw("same");
986 }
987
988 DisplayRefLines(null,refline);
989
990
991 newpad->cd(2);
992 gPad->SetTicks();
993
994 pix.DrawEvents();
995 return;
996
997}
998
999
1000void MHCalibrationChargeCam::DisplayRefLines(const TH2D *hist, const Float_t refline[]) const
1001{
1002
1003 TLine *green1 = new TLine(refline[0],0.,refline[0],hist->GetYaxis()->GetXmax());
1004 green1->SetBit(kCanDelete);
1005 green1->SetLineColor(kGreen);
1006 green1->SetLineStyle(2);
1007 green1->SetLineWidth(3);
1008 green1->Draw();
1009
1010 TLine *green5 = new TLine(refline[6],0.,refline[6],hist->GetYaxis()->GetXmax());
1011 green5->SetBit(kCanDelete);
1012 green5->SetLineColor(8);
1013 green5->SetLineStyle(2);
1014 green5->SetLineWidth(3);
1015 green5->Draw();
1016
1017 TLine *blue1 = new TLine(refline[1],0.,refline[1],hist->GetYaxis()->GetXmax());
1018 blue1->SetBit(kCanDelete);
1019 blue1->SetLineColor(007);
1020 blue1->SetLineStyle(2);
1021 blue1->SetLineWidth(3);
1022 blue1->Draw();
1023
1024 TLine *blue5 = new TLine(refline[2],0.,refline[2],hist->GetYaxis()->GetXmax());
1025 blue5->SetBit(kCanDelete);
1026 blue5->SetLineColor(062);
1027 blue5->SetLineStyle(2);
1028 blue5->SetLineWidth(3);
1029 blue5->Draw();
1030
1031 TLine *blue10 = new TLine(refline[3],0.,refline[3],hist->GetYaxis()->GetXmax());
1032 blue10->SetBit(kCanDelete);
1033 blue10->SetLineColor(004);
1034 blue10->SetLineStyle(2);
1035 blue10->SetLineWidth(3);
1036 blue10->Draw();
1037
1038 TLine *uv10 = new TLine(refline[4],0.,refline[4],hist->GetYaxis()->GetXmax());
1039 uv10->SetBit(kCanDelete);
1040 uv10->SetLineColor(106);
1041 uv10->SetLineStyle(2);
1042 uv10->SetLineWidth(3);
1043 uv10->Draw();
1044
1045 TLine *ct1 = new TLine(refline[5],0.,refline[5],hist->GetYaxis()->GetXmax());
1046 ct1->SetBit(kCanDelete);
1047 ct1->SetLineColor(006);
1048 ct1->SetLineStyle(2);
1049 ct1->SetLineWidth(3);
1050 ct1->Draw();
1051
1052 TLegend *leg = new TLegend(0.4,0.75,0.7,0.99);
1053 leg->SetBit(kCanDelete);
1054 leg->AddEntry(green1,"1 Led GREEN","l");
1055 leg->AddEntry(green5,"5 Leds GREEN","l");
1056 leg->AddEntry(blue1,"1 Led BLUE","l");
1057 leg->AddEntry(blue5,"5 Leds BLUE","l");
1058 leg->AddEntry(blue10,"10 Leds BLUE","l");
1059 leg->AddEntry(uv10,"10 Leds UV","l");
1060 leg->AddEntry(ct1,"CT1-Pulser","l");
1061
1062 leg->Draw();
1063}
Note: See TracBrowser for help on using the repository browser.