source: releases/Mars.2014.05.26/mhcalib/MHCalibrationChargeCam.cc@ 19858

Last change on this file since 19858 was 10166, checked in by tbretz, 14 years ago
Removed the old obsolete cvs header line.
File size: 47.7 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! Author(s): Thomas Bretz <mailto:tbretz@astro.uni-wuerzburg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2007
22!
23!
24\* ======================================================================== */
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHCalibrationChargeCam
28//
29// Fills the extracted signals of MExtractedSignalCam into the MHCalibrationPix-classes
30// MHCalibrationChargeHiGainPix and MHCalibrationChargeLoGainPix for every:
31//
32// - Pixel, stored in the TOrdCollection's MHCalibrationCam::fHiGainArray and
33// MHCalibrationCam::fLoGainArray
34//
35// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
36// stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainAreas and
37// MHCalibrationCam::fAverageLoGainAreas
38//
39// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
40// stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainSectors and
41// MHCalibrationCam::fAverageLoGainSectors
42//
43// Every signal is taken from MExtractedSignalCam and filled into a histogram and
44// an array, in order to perform a Fourier analysis (see MHGausEvents).
45// The signals are moreover averaged on an event-by-event basis and written into
46// the corresponding average pixels.
47//
48// Additionally, the (FADC slice) position of the maximum is stored in an Absolute
49// Arrival Time histogram. This histogram serves for a rough cross-check if the
50// signal does not lie at or outside the edges of the extraction window.
51//
52// The Charge histograms are fitted to a Gaussian, mean and sigma with its errors
53// and the fit probability are extracted. If none of these values are NaN's and
54// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
55// the fit is declared valid.
56// Otherwise, the fit is repeated within ranges of the previous mean
57// +- MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
58// In case this does not make the fit valid, the histogram means and RMS's are
59// taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
60// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted ) or
61// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted ) and
62// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
63//
64// Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
65// from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
66//
67// If more than fNumHiGainSaturationLimit (default: 15%) of the overall FADC
68// slices show saturation, the following flag is set:
69// - MCalibrationChargePix::SetHiGainSaturation();
70// In that case, the calibration constants are derived from the low-gain results.
71//
72// If more than fNumLoGainSaturationLimit (default: 1%) of the overall
73// low-gain FADC slices saturate, the following flags are set:
74// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturation ) and
75// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnsuitableRun )
76//
77// The class also fills arrays with the signal vs. event number, creates a fourier
78// spectrum and investigates if the projected fourier components follow an exponential
79// distribution. In case that the probability of the exponential fit is less than
80// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
81// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating ) or
82// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating ) and
83// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
84//
85// This same procedure is performed for the average pixels.
86//
87// The following results are written into MCalibrationChargeCam:
88//
89// - MCalibrationPix::SetHiGainSaturation()
90// - MCalibrationPix::SetHiGainMean()
91// - MCalibrationPix::SetHiGainMeanErr()
92// - MCalibrationPix::SetHiGainSigma()
93// - MCalibrationPix::SetHiGainSigmaErr()
94// - MCalibrationPix::SetHiGainProb()
95// - MCalibrationPix::SetHiGainNumPickup()
96//
97// - MCalibrationPix::SetLoGainMean()
98// - MCalibrationPix::SetLoGainMeanErr()
99// - MCalibrationPix::SetLoGainSigma()
100// - MCalibrationPix::SetLoGainSigmaErr()
101// - MCalibrationPix::SetLoGainProb()
102// - MCalibrationPix::SetLoGainNumPickup()
103//
104// - MCalibrationChargePix::SetAbsTimeMean()
105// - MCalibrationChargePix::SetAbsTimeRms()
106//
107// For all averaged areas, the fitted sigma is multiplied with the square root of
108// the number involved pixels in order to be able to compare it to the average of
109// sigmas in the camera.
110//
111// Class Version 2:
112// + Float_t fNumLoGainBlackoutLimit; // Rel. amount blackout logain events until pixel is declared unsuitable
113//
114// Class Version 3:
115// - Byte_t fFirstHiGain; // First used slice High Gain
116// - Byte_t fLastHiGain; // Last used slice High Gain
117// - Byte_t fFirstLoGain; // First used slice Low Gain
118// - Byte_t fLastLoGain; // Last used slice Low Gain
119//
120/////////////////////////////////////////////////////////////////////////////
121#include "MHCalibrationChargeCam.h"
122
123#include <TOrdCollection.h>
124#include <TPad.h>
125#include <TVirtualPad.h>
126#include <TCanvas.h>
127#include <TStyle.h>
128#include <TF1.h>
129#include <TLine.h>
130#include <TLatex.h>
131#include <TLegend.h>
132#include <TGraph.h>
133#include <TEnv.h>
134
135#include "MLog.h"
136#include "MLogManip.h"
137
138#include "MParList.h"
139
140#include "MHCalibrationChargePix.h"
141#include "MHCalibrationPix.h"
142
143#include "MCalibrationChargeCam.h"
144#include "MCalibrationChargePix.h"
145
146#include "MGeomCam.h"
147#include "MGeom.h"
148
149#include "MBadPixelsCam.h"
150#include "MBadPixelsPix.h"
151
152#include "MRawRunHeader.h"
153#include "MRawEvtPixelIter.h"
154
155#include "MPedestalSubtractedEvt.h"
156
157#include "MExtractedSignalCam.h"
158#include "MExtractedSignalPix.h"
159
160ClassImp(MHCalibrationChargeCam);
161
162using namespace std;
163
164const Int_t MHCalibrationChargeCam::fgChargeHiGainNbins = 500;
165const Axis_t MHCalibrationChargeCam::fgChargeHiGainFirst = -98.;
166const Axis_t MHCalibrationChargeCam::fgChargeHiGainLast = 1902.;
167const Int_t MHCalibrationChargeCam::fgChargeLoGainNbins = 500;
168const Axis_t MHCalibrationChargeCam::fgChargeLoGainFirst = -99.;
169const Axis_t MHCalibrationChargeCam::fgChargeLoGainLast = 901.;
170const Float_t MHCalibrationChargeCam::fgProbLimit = 0.00000001;
171const TString MHCalibrationChargeCam::gsHistName = "Charge";
172const TString MHCalibrationChargeCam::gsHistTitle = "Signals";
173const TString MHCalibrationChargeCam::gsHistXTitle = "Signal [FADC counts]";
174const TString MHCalibrationChargeCam::gsHistYTitle = "Nr. events";
175const TString MHCalibrationChargeCam::gsAbsHistName = "AbsTime";
176const TString MHCalibrationChargeCam::gsAbsHistTitle = "Abs. Arr. Times";
177const TString MHCalibrationChargeCam::gsAbsHistXTitle = "Time [FADC slices]";
178const TString MHCalibrationChargeCam::gsAbsHistYTitle = "Nr. events";
179const Float_t MHCalibrationChargeCam::fgNumHiGainSaturationLimit = 0.15;
180const Float_t MHCalibrationChargeCam::fgNumLoGainSaturationLimit = 0.005;
181const Float_t MHCalibrationChargeCam::fgNumLoGainBlackoutLimit = 0.05;
182const Float_t MHCalibrationChargeCam::fgLoGainBlackoutLimit = 3.5;
183const Float_t MHCalibrationChargeCam::fgLoGainPickupLimit = 3.5;
184const Float_t MHCalibrationChargeCam::fgTimeLowerLimit = 1.;
185const Float_t MHCalibrationChargeCam::fgTimeUpperLimit = 3.;
186const TString MHCalibrationChargeCam::fgReferenceFile = "mjobs/calibrationref.rc";
187
188// --------------------------------------------------------------------------
189//
190// Default Constructor.
191//
192// Sets:
193// - all pointers to NULL
194//
195// Initializes:
196// - fNumHiGainSaturationLimit to fgNumHiGainSaturationLimit
197// - fNumLoGainSaturationLimit to fgNumLoGainSaturationLimit
198// - fTimeLowerLimit to fgTimeLowerLimit
199// - fTimeUpperLimit to fgTimeUpperLimit
200// - fNumLoGainBlackoutLimit to fgNumLoGainBlackoutLimit
201//
202// - fNbins to fgChargeHiGainNbins
203// - fFirst to fgChargeHiGainFirst
204// - fLast to fgChargeHiGainLast
205//
206// - fLoGainNbins to fgChargeLoGainNbins
207// - fLoGainFirst to fgChargeLoGainFirst
208// - fLoGainLast to fgChargeLoGainLast
209//
210// - fHistName to gsHistName
211// - fHistTitle to gsHistTitle
212// - fHistXTitle to gsHistXTitle
213// - fHistYTitle to gsHistYTitle
214//
215// - fAbsHistName to gsAbsHistName
216// - fAbsHistTitle to gsAbsHistTitle
217// - fAbsHistXTitle to gsAbsHistXTitle
218// - fAbsHistYTitle to gsAbsHistYTitle
219//
220MHCalibrationChargeCam::MHCalibrationChargeCam(const char *name, const char *title)
221 : fRawEvt(NULL), fSignal(NULL)
222{
223
224 fName = name ? name : "MHCalibrationChargeCam";
225 fTitle = title ? title : "Class to fill the calibration histograms ";
226
227 SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
228 SetNumLoGainSaturationLimit(fgNumLoGainSaturationLimit);
229
230 SetNumLoGainBlackoutLimit (fgNumLoGainBlackoutLimit);
231
232 SetTimeLowerLimit();
233 SetTimeUpperLimit();
234
235 SetBinning(fgChargeHiGainNbins, fgChargeHiGainFirst, fgChargeHiGainLast );
236 SetBinningLoGain(fgChargeLoGainNbins, fgChargeLoGainFirst, fgChargeLoGainLast );
237
238 SetProbLimit(fgProbLimit);
239
240 SetHistName (gsHistName .Data());
241 SetHistTitle (gsHistTitle .Data());
242 SetHistXTitle(gsHistXTitle.Data());
243 SetHistYTitle(gsHistYTitle.Data());
244
245 SetAbsHistName (gsAbsHistName .Data());
246 SetAbsHistTitle (gsAbsHistTitle .Data());
247 SetAbsHistXTitle(gsAbsHistXTitle.Data());
248 SetAbsHistYTitle(gsAbsHistYTitle.Data());
249
250 SetReferenceFile();
251
252 fInnerRefCharge = 278.;
253 fOuterRefCharge = 282.;
254}
255
256// --------------------------------------------------------------------------
257//
258// Creates new MHCalibrationChargeCam only with the averaged areas:
259// the rest has to be retrieved directly, e.g. via:
260// MHCalibrationChargeCam *cam = MParList::FindObject("MHCalibrationChargeCam");
261// - cam->GetAverageSector(5).DrawClone();
262// - (*cam)[100].DrawClone()
263//
264TObject *MHCalibrationChargeCam::Clone(const char *) const
265{
266
267 MHCalibrationChargeCam *cam = new MHCalibrationChargeCam();
268
269 //
270 // Copy the data members
271 //
272 cam->fColor = fColor;
273 cam->fRunNumbers = fRunNumbers;
274 cam->fPulserFrequency = fPulserFrequency;
275 cam->fFlags = fFlags;
276 cam->fNbins = fNbins;
277 cam->fFirst = fFirst;
278 cam->fLast = fLast;
279 cam->fLoGainNbins = fLoGainNbins;
280 cam->fLoGainFirst = fLoGainFirst;
281 cam->fLoGainLast = fLoGainLast;
282 cam->fReferenceFile = fReferenceFile;
283 cam->fInnerRefCharge = fInnerRefCharge;
284 cam->fOuterRefCharge = fOuterRefCharge;
285
286 //
287 // Copy the MArrays
288 //
289 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
290 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
291 cam->fAverageAreaSat = fAverageAreaSat;
292 cam->fAverageAreaSigma = fAverageAreaSigma;
293 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
294 cam->fAverageAreaNum = fAverageAreaNum;
295 cam->fAverageSectorNum = fAverageSectorNum;
296
297 if (!IsAverageing())
298 return cam;
299
300 const Int_t navhi = fAverageHiGainAreas->GetSize();
301
302 for (int i=0; i<navhi; i++)
303 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
304
305 if (IsLoGain())
306 {
307
308 const Int_t navlo = fAverageLoGainAreas->GetSize();
309 for (int i=0; i<navlo; i++)
310 cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
311
312 }
313
314 return cam;
315}
316
317// --------------------------------------------------------------------------
318//
319// Gets the pointers to:
320// - MRawEvtData
321//
322Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
323{
324 fRawEvt = (MPedestalSubtractedEvt*)pList->FindObject("MPedestalSubtractedEvt");
325 if (!fRawEvt)
326 {
327 *fLog << err << dbginf << "MPedestalSubtractedEvt not found... aborting." << endl;
328 return kFALSE;
329 }
330
331 return kTRUE;
332}
333
334// --------------------------------------------------------------------------
335//
336// Gets or creates the pointers to:
337// - MExtractedSignalCam
338// - MCalibrationChargeCam
339// - MBadPixelsCam
340//
341// Initializes the number of used FADC slices from MExtractedSignalCam
342// into MCalibrationChargeCam and test for changes in that variable
343//
344// Calls:
345// - InitHiGainArrays()
346// - InitLoGainArrays()
347//
348// Sets:
349// - fSumhiarea to nareas
350// - fSumloarea to nareas
351// - fTimehiarea to nareas
352// - fTimeloarea to nareas
353// - fSumhisector to nsectors
354// - fSumlosector to nsectors
355// - fTimehisector to nsectors
356// - fTimelosector to nsectors
357// - fSathiarea to nareas
358// - fSatloarea to nareas
359// - fSathisector to nsectors
360// - fSatlosector to nsectors
361//
362Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
363{
364
365 fSignal = (MExtractedSignalCam*)pList->FindObject(AddSerialNumber("MExtractedSignalCam"));
366 if (!fSignal)
367 {
368 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
369 return kFALSE;
370 }
371
372 if (!InitCams(pList,"Charge"))
373 return kFALSE;
374
375 const Int_t npixels = fGeom->GetNumPixels();
376 const Int_t nsectors = fGeom->GetNumSectors();
377 const Int_t nareas = fGeom->GetNumAreas();
378
379 //
380 // In case of the intense blue, double the range
381 //
382 if (fGeom->InheritsFrom("MGeomCamMagic"))
383 if ( fColor == MCalibrationCam::kBLUE)
384 SetBinning(fLoGainNbins, fLoGainFirst, 2.*fLoGainLast - fLoGainFirst);
385
386 InitHiGainArrays(npixels,nareas,nsectors);
387 InitLoGainArrays(npixels,nareas,nsectors);
388
389 fSumhiarea .Set(nareas);
390 fSumloarea .Set(nareas);
391 fTimehiarea .Set(nareas);
392 fTimeloarea .Set(nareas);
393 fSumhisector.Set(nsectors);
394 fSumlosector.Set(nsectors);
395 fTimehisector.Set(nsectors);
396 fTimelosector.Set(nsectors);
397
398 fSathiarea .Set(nareas);
399 fSatloarea .Set(nareas);
400 fSathisector.Set(nsectors);
401 fSatlosector.Set(nsectors);
402
403 return kTRUE;
404}
405
406// --------------------------------------------------------------------------
407//
408// Retrieve:
409// - fRunHeader->GetNumSamplesHiGain();
410//
411// Initializes the High Gain Arrays:
412//
413// - For every entry in the expanded arrays:
414// * Initialize an MHCalibrationChargePix
415// * Set Binning from fNbins, fFirst and fLast
416// * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
417// * Set Histgram names and titles from fHistName and fHistTitle
418// * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
419// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
420// * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
421// * Call InitHists
422//
423//
424void MHCalibrationChargeCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
425{
426
427 TH1F *h;
428
429 const Int_t higainsamples = fRunHeader->GetNumSamplesHiGain();
430
431 if (fHiGainArray->GetSize()==0)
432 {
433 for (Int_t i=0; i<npixels; i++)
434 {
435 fHiGainArray->AddAt(new MHCalibrationChargePix(Form("%sHiGainPix%04d",fHistName.Data(),i),
436 Form("%s High Gain Pixel%04d",fHistTitle.Data(),i)),i);
437
438 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)[i];
439
440 pix.SetBinning(fNbins, fFirst, fLast);
441 pix.SetProbLimit(fProbLimit);
442
443 pix.SetBinningAbsTime(higainsamples, -0.5, higainsamples-0.5);
444
445 InitHists(pix,(*fBadPixels)[i], i);
446
447 h = pix.GetHAbsTime();
448
449 h->SetName (Form("H%sHiGainPix%04d",fAbsHistName.Data(),i));
450 h->SetTitle(Form("%s High Gain Pixel %04d",fAbsHistTitle.Data(),i));
451 h->SetXTitle(fAbsHistXTitle.Data());
452 h->SetYTitle(fAbsHistYTitle.Data());
453 }
454 }
455
456
457 if (fAverageHiGainAreas->GetSize()==0)
458 {
459 for (Int_t j=0; j<nareas; j++)
460 {
461 fAverageHiGainAreas->AddAt(new MHCalibrationChargePix(Form("%sHiGainArea%d",fHistName.Data(),j),
462 Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
463
464 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
465
466 pix.SetBinning(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas),
467 fFirst, fLast);
468
469 pix.SetBinningAbsTime(higainsamples, -0.5, higainsamples-0.5);
470
471 InitHists(pix, fCam->GetAverageBadArea(j),j);
472
473 h = pix.GetHAbsTime();
474
475 h->SetName (Form("H%sHiGainArea%d",fAbsHistName.Data(),j));
476 h->SetTitle(Form("%s%s%d",fAbsHistTitle.Data(),
477 " averaged on event-by-event basis High Gain Area Idx ",j));
478 h->SetXTitle(fAbsHistXTitle.Data());
479 h->SetYTitle(fAbsHistYTitle.Data());
480 }
481 }
482
483 if (fAverageHiGainSectors->GetSize()==0)
484 {
485 for (Int_t j=0; j<nsectors; j++)
486 {
487 fAverageHiGainSectors->AddAt(new MHCalibrationChargePix(Form("%sHiGainSector%02d",fHistName.Data(),j),
488 Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
489
490 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
491
492 pix.SetBinning(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas),
493 fFirst, fLast);
494
495 pix.SetBinningAbsTime(higainsamples, -0.5, higainsamples-0.5);
496
497 InitHists(pix, fCam->GetAverageBadSector(j),j);
498
499 h = pix.GetHAbsTime();
500
501 h->SetName (Form("H%sHiGainSector%02d",fAbsHistName.Data(),j));
502 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
503 " averaged on event-by-event basis High Gain Area Sector ",j));
504 h->SetXTitle(fAbsHistXTitle.Data());
505 h->SetYTitle(fAbsHistYTitle.Data());
506 }
507 }
508}
509
510//--------------------------------------------------------------------------------------
511//
512// Return, if IsLoGain() is kFALSE
513//
514// Retrieve:
515// - fRunHeader->GetNumSamplesHiGain();
516//
517// Initializes the Low Gain Arrays:
518//
519// - For every entry in the expanded arrays:
520// * Initialize an MHCalibrationChargePix
521// * Set Binning from fNbins, fFirst and fLast
522// * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
523// * Set Histgram names and titles from fHistName and fHistTitle
524// * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
525// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
526// * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
527// * Call InitHists
528//
529void MHCalibrationChargeCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
530{
531
532 if (!IsLoGain())
533 return;
534
535
536 const Int_t logainsamples = fRunHeader->GetNumSamplesLoGain();
537
538 TH1F *h;
539
540 if (fLoGainArray->GetSize()==0 )
541 {
542 for (Int_t i=0; i<npixels; i++)
543 {
544 fLoGainArray->AddAt(new MHCalibrationChargePix(Form("%sLoGainPix%04d",fHistName.Data(),i),
545 Form("%s Low Gain Pixel %04d",fHistTitle.Data(),i)),i);
546
547 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)(i);
548
549 pix.SetBinning(fLoGainNbins, fLoGainFirst, fLoGainLast);
550 pix.SetProbLimit(fProbLimit);
551
552 pix.SetBinningAbsTime(logainsamples, -0.5, logainsamples-0.5);
553 pix.SetPickupLimit(fgLoGainPickupLimit);
554 pix.SetBlackoutLimit(fgLoGainBlackoutLimit);
555
556 InitHists(pix,(*fBadPixels)[i], i);
557
558 h = pix.GetHAbsTime();
559
560 h->SetName (Form("H%sLoGainPix%04d",fAbsHistName.Data(),i));
561 h->SetTitle(Form("%s Low Gain Pixel %04d",fAbsHistTitle.Data(),i));
562 h->SetXTitle(fAbsHistXTitle.Data());
563 h->SetYTitle(fAbsHistYTitle.Data());
564 }
565 }
566
567 if (fAverageLoGainAreas->GetSize()==0)
568 {
569 for (Int_t j=0; j<nareas; j++)
570 {
571 fAverageLoGainAreas->AddAt(new MHCalibrationChargePix(Form("%sLoGainArea%d",fHistName.Data(),j),
572 Form("%s Low Gain Area Idx %d",fHistTitle.Data(),j)),j);
573
574 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
575
576 pix.SetBinning(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas),
577 fLoGainFirst, fLoGainLast);
578
579 pix.SetBinningAbsTime(logainsamples, -0.5, logainsamples-0.5);
580
581 InitHists(pix, fCam->GetAverageBadArea(j),j);
582
583 h = pix.GetHAbsTime();
584
585 h->SetName (Form("H%sLoGainArea%02d",fAbsHistName.Data(),j));
586 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
587 " averaged on event-by-event basis Low Gain Area Idx ",j));
588 h->SetXTitle(fAbsHistXTitle.Data());
589 h->SetYTitle(fAbsHistYTitle.Data());
590 }
591 }
592
593
594 if (fAverageLoGainSectors->GetSize()==0 && IsLoGain())
595 {
596 for (Int_t j=0; j<nsectors; j++)
597 {
598 fAverageLoGainSectors->AddAt(new MHCalibrationChargePix(Form("%sLoGainSector%02d",fHistName.Data(),j),
599 Form("%s Low Gain Sector %02d",fHistTitle.Data(),j)),j);
600
601 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
602
603 pix.SetBinning(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas),
604 fLoGainFirst, fLoGainLast);
605
606 pix.SetBinningAbsTime(logainsamples, -0.5, logainsamples-0.5);
607
608 InitHists(pix, fCam->GetAverageBadSector(j),j);
609
610 h = pix.GetHAbsTime();
611
612 h->SetName (Form("H%sLoGainSector%02d",fAbsHistName.Data(),j));
613 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
614 " averaged on event-by-event basis Low Gain Area Sector ",j));
615 h->SetXTitle(fAbsHistXTitle.Data());
616 h->SetYTitle(fAbsHistYTitle.Data());
617 }
618 }
619}
620
621
622// --------------------------------------------------------------------------
623//
624// Retrieves from MExtractedSignalCam:
625// - first used LoGain FADC slice
626//
627// Retrieves from MGeomCam:
628// - number of pixels
629// - number of pixel areas
630// - number of sectors
631//
632// For all TOrdCollection's (including the averaged ones), the following steps are performed:
633//
634// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
635// - MExtractedSignalPix::GetExtractedSignalHiGain();
636// - MExtractedSignalPix::GetExtractedSignalLoGain();
637//
638// 2) Set number of saturated slices (MHCalibrationChargePix::AddSaturated()) with:
639// - MExtractedSignalPix::IsHiGainSaturated();
640// - MExtractedSignalPix::GetNumLoGainSaturated();
641//
642// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
643// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
644// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
645//
646Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
647{
648
649 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
650 if (!signal)
651 {
652 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
653 return kFALSE;
654 }
655
656 const UInt_t npixels = fGeom->GetNumPixels();
657 const UInt_t nareas = fGeom->GetNumAreas();
658 const UInt_t nsectors = fGeom->GetNumSectors();
659
660 const Int_t hifirst = fSignal->GetFirstUsedSliceHiGain();
661 const Int_t hilast = fSignal->GetLastUsedSliceHiGain();
662 const Int_t lofirst = fSignal->GetFirstUsedSliceLoGain();
663 const Int_t lolast = fSignal->GetLastUsedSliceLoGain();
664
665 fSumhiarea .Reset();
666 fSumloarea .Reset();
667 fTimehiarea .Reset();
668 fTimeloarea .Reset();
669 fSumhisector.Reset();
670 fSumlosector.Reset();
671 fTimehisector.Reset();
672 fTimelosector.Reset();
673
674 fSathiarea .Reset();
675 fSatloarea .Reset();
676 fSathisector.Reset();
677 fSatlosector.Reset();
678
679 for (UInt_t i=0; i<npixels; i++)
680 {
681 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
682 if (histhi.IsExcluded())
683 continue;
684
685 const Int_t aidx = (*fGeom)[i].GetAidx();
686 const Int_t sector = (*fGeom)[i].GetSector();
687
688 const MExtractedSignalPix &pix = (*signal)[i];
689
690 const Int_t sathi = pix.IsHiGainSaturated() ? 1 : 0;
691 histhi.AddSaturated(sathi);
692 if (sathi)
693 {
694 fSathiarea[aidx]++;
695 fSathisector[sector]++;
696 }
697
698 if (pix.IsHiGainValid())
699 {
700 const Float_t sumhi = pix.GetExtractedSignalHiGain();
701 if (IsOscillations())
702 histhi.FillHistAndArray(sumhi);
703 else
704 histhi.FillHist(sumhi);
705
706 fSumhiarea[aidx] += sumhi;
707 fSumhisector[sector] += sumhi;
708 }
709
710 if (!IsLoGain())
711 continue;
712
713 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
714
715 const Int_t satlo = pix.IsLoGainSaturated() ? 1 : 0;
716 histlo.AddSaturated(satlo);
717 if (satlo)
718 {
719 fSatloarea[aidx]++;
720 fSatlosector[sector]++;
721 }
722
723 /*
724 // Previously:
725 fSatloarea[aidx] += pix.GetNumHiGainSaturated();
726 fSatlosector[sector] += pix.GetNumHiGainSaturated();
727 */
728
729 if (pix.IsLoGainValid())
730 {
731 const Float_t sumlo = pix.GetExtractedSignalLoGain();
732
733 if (IsOscillations())
734 histlo.FillHistAndArray(sumlo);
735 else
736 histlo.FillHist(sumlo);
737
738 fSumloarea[aidx] += sumlo;
739 fSumlosector[sector] += sumlo;
740 }
741 }
742
743 const Int_t num = fRawEvt->GetNumPixels();
744 const Int_t nhi = fRunHeader->GetNumSamplesHiGain();
745 for (int pixid=0; pixid<num; pixid++)
746 {
747 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
748 if (histhi.IsExcluded())
749 continue;
750
751 const Float_t timehi = fRawEvt->GetMaxPos(pixid, hifirst, hilast)+hifirst;
752 histhi.FillAbsTime(timehi);
753
754 const Int_t aidx = (*fGeom)[pixid].GetAidx();
755 const Int_t sector = (*fGeom)[pixid].GetSector();
756
757 fTimehiarea [aidx] += timehi;
758 fTimehisector[sector] += timehi;
759
760 if (IsLoGain())
761 {
762 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
763
764 const Float_t timelo = fRawEvt->GetMaxPos(pixid, nhi+lofirst, nhi+lolast)+lofirst;
765 histlo.FillAbsTime(timelo);
766
767 fTimeloarea[aidx] += timelo;
768 fTimelosector[sector] += timelo;
769 }
770 }
771
772 for (UInt_t j=0; j<nareas; j++)
773 {
774
775 const Int_t npix = fAverageAreaNum[j];
776
777 if (npix == 0)
778 continue;
779
780 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
781
782 if (IsOscillations())
783 hipix.FillHistAndArray(fSumhiarea [j]/npix);
784 else
785 hipix.FillHist(fSumhiarea[j]/npix);
786
787 hipix.AddSaturated ((Float_t)fSathiarea [j]/npix > 0.5 ? 1 : 0);
788 hipix.FillAbsTime (fTimehiarea[j]/npix);
789
790 if (IsLoGain())
791 {
792 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
793 if (IsOscillations())
794 lopix.FillHistAndArray(fSumloarea [j]/npix);
795 else
796 lopix.FillHist(fSumloarea [j]/npix);
797 lopix.AddSaturated ((Float_t)fSatloarea [j]/npix > 0.5 ? 1 : 0);
798 lopix.FillAbsTime (fTimeloarea[j]/npix);
799 }
800 }
801
802 for (UInt_t j=0; j<nsectors; j++)
803 {
804
805 const Int_t npix = fAverageSectorNum[j];
806
807 if (npix == 0)
808 continue;
809
810 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
811
812 if (IsOscillations())
813 hipix.FillHistAndArray(fSumhisector [j]/npix);
814 else
815 hipix.FillHist(fSumhisector [j]/npix);
816
817 hipix.AddSaturated ((Float_t)fSathisector[j]/npix > 0.5 ? 1 : 0);
818 hipix.FillAbsTime (fTimehisector[j]/npix);
819
820 if (IsLoGain())
821 {
822 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
823
824 if (IsOscillations())
825 lopix.FillHistAndArray(fSumlosector [j]/npix);
826 else
827 lopix.FillHist(fSumlosector [j]/npix);
828
829 lopix.AddSaturated ((Float_t)fSatlosector[j]/npix > 0.5 ? 1 : 0);
830 lopix.FillAbsTime (fTimelosector[j]/npix);
831 }
832 }
833
834 return kTRUE;
835}
836
837// --------------------------------------------------------------------------
838//
839// For all TOrdCollection's (including the averaged ones), the following steps are performed:
840//
841// 1) Returns if the pixel is excluded.
842// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
843// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
844// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
845// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
846// otherwise the Hi-Gain ones.
847// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
848// with the flags:
849// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
850// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
851// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
852// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
853//
854Bool_t MHCalibrationChargeCam::FinalizeHists()
855{
856
857 *fLog << endl;
858
859 TH1F *h = NULL;
860
861 const Int_t hifirst = fSignal->GetFirstUsedSliceHiGain();
862 const Int_t hilast = fSignal->GetLastUsedSliceHiGain();
863 const Int_t lofirst = fSignal->GetFirstUsedSliceLoGain();
864 const Int_t lolast = fSignal->GetLastUsedSliceLoGain();
865
866 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
867 {
868
869 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
870 if (histhi.IsExcluded())
871 continue;
872
873 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
874
875 const Int_t numsat = histhi.GetSaturated();
876
877 pix.SetNumSaturated(numsat);
878
879 if (numsat > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
880 {
881 pix.SetHiGainSaturation();
882 if (IsOscillations())
883 histhi.CreateFourierSpectrum();
884 continue;
885 }
886
887 MBadPixelsPix &bad = (*fBadPixels)[i];
888
889 h = histhi.GetHGausHist();
890
891 Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
892 if (overflow > fOverflowLimit*histhi.GetHGausHist()->GetEntries())
893 {
894 *fLog << warn
895 << "HiGain Hist-overflow " << overflow
896 << " times in: " << histhi.GetName() << " (w/o saturation!) " << endl;
897 bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
898 }
899
900 overflow = h->GetBinContent(0);
901 if (overflow > fOverflowLimit*histhi.GetHGausHist()->GetEntries())
902 {
903 *fLog << warn
904 << "HiGain Hist-underflow " << overflow
905 << " times in pix: " << histhi.GetName() << " (w/o saturation!) " << endl;
906 bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
907 }
908
909 FinalizeAbsTimes(histhi, pix, bad, hifirst, hilast);
910 }
911
912 if (IsLoGain())
913 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
914 {
915
916 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
917
918 if (histlo.IsExcluded())
919 continue;
920
921 MBadPixelsPix &bad = (*fBadPixels)[i];
922 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
923 {
924 *fLog << warn << "Pixel " << setw(4) << i << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
925 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
926 if (IsOscillations())
927 histlo.CreateFourierSpectrum();
928 continue;
929 }
930
931 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i] ;
932 if (!pix.IsHiGainSaturation())
933 continue;
934
935 h = histlo.GetHGausHist();
936
937 Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
938 if (overflow > fOverflowLimit*histlo.GetHGausHist()->GetEntries())
939 {
940 *fLog << warn
941 << "LoGain Hist-overflow " << overflow
942 << " times in: " << histlo.GetName() << " (w/o saturation!) " << endl;
943 bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
944 }
945
946 overflow = h->GetBinContent(0);
947 if (overflow > fOverflowLimit*histlo.GetHGausHist()->GetEntries())
948 {
949 *fLog << warn
950 << "LoGain Hist-underflow " << overflow
951 << " times in: " << histlo.GetName() << " (w/o saturation!) " << endl;
952 bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
953 }
954
955 FinalizeAbsTimes(histlo, pix, bad, lofirst, lolast);
956 }
957
958 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
959 {
960
961 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
962 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
963
964 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
965 {
966 pix.SetHiGainSaturation();
967 if (IsOscillations())
968 histhi.CreateFourierSpectrum();
969 continue;
970 }
971
972 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
973 FinalizeAbsTimes(histhi, pix, bad, hifirst, hilast);
974 }
975
976 if (IsLoGain())
977 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
978 {
979
980 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
981
982 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
983 {
984 *fLog << warn << "Area " << setw(4) << j << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
985 if (IsOscillations())
986 histlo.CreateFourierSpectrum();
987 continue;
988 }
989
990 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageArea(j);
991 if (pix.IsHiGainSaturation())
992 {
993 MBadPixelsPix &bad = fCam->GetAverageBadArea(j);
994 FinalizeAbsTimes(histlo, pix, bad, lofirst, lolast);
995 }
996
997 }
998
999 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
1000 {
1001
1002 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
1003 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
1004
1005 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
1006 {
1007 pix.SetHiGainSaturation();
1008 if (IsOscillations())
1009 histhi.CreateFourierSpectrum();
1010 continue;
1011 }
1012
1013 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
1014 FinalizeAbsTimes(histhi, pix, bad, hifirst, hilast);
1015 }
1016
1017 if (IsLoGain())
1018 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
1019 {
1020
1021 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
1022 MBadPixelsPix &bad = fCam->GetAverageBadSector(j);
1023
1024 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
1025 {
1026 *fLog << warn << "Sector " << setw(4) << j << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
1027 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
1028 if (IsOscillations())
1029 histlo.CreateFourierSpectrum();
1030 continue;
1031 }
1032
1033 MCalibrationChargePix &pix = (MCalibrationChargePix&)fCam->GetAverageSector(j);
1034 if (pix.IsHiGainSaturation())
1035 FinalizeAbsTimes(histlo, pix, bad, lofirst, lolast);
1036 }
1037
1038 //
1039 // Perform the fitting for the High Gain (done in MHCalibrationCam)
1040 //
1041 FitHiGainArrays(*fCam, *fBadPixels,
1042 MBadPixelsPix::kHiGainNotFitted,
1043 MBadPixelsPix::kHiGainOscillating);
1044
1045 //
1046 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
1047 //
1048 if (IsLoGain())
1049 FitLoGainArrays(*fCam, *fBadPixels,
1050 MBadPixelsPix::kLoGainNotFitted,
1051 MBadPixelsPix::kLoGainOscillating);
1052
1053 //
1054 // Check for pixels which have the high-gain saturated, but the low-gain
1055 // switch not applied in sufficient cases. Have to exclude these pixels,
1056 // although they not abnormal. They simply cannot be calibrated...
1057 //
1058 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
1059 {
1060
1061 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
1062 if (histlo.IsExcluded())
1063 continue;
1064
1065 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*fCam)[i];
1066 if (!pix.IsHiGainSaturation())
1067 continue;
1068
1069 //
1070 // Now,treat only low-gain calibrated pixels:
1071 //
1072 const Double_t lim = fNumLoGainBlackoutLimit*histlo.GetHGausHist()->GetEntries();
1073 if (histlo.GetBlackout() <= lim)
1074 continue;
1075
1076 MBadPixelsPix &bad = (*fBadPixels)[i];
1077 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1078 bad.SetUncalibrated(MBadPixelsPix::kLoGainBlackout);
1079 }
1080
1081 return kTRUE;
1082}
1083
1084// --------------------------------------------------------------------------------
1085//
1086// Fill the absolute time results into MCalibrationChargePix
1087//
1088// Check absolute time validity:
1089// - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
1090// - Mean arrival time is at least fUpperLimit slices from the upper edge
1091//
1092void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
1093 Int_t first, Int_t last)
1094{
1095 const Float_t mean = hist.GetAbsTimeMean();
1096 const Float_t rms = hist.GetAbsTimeRms();
1097
1098 pix.SetAbsTimeMean(mean);
1099 pix.SetAbsTimeRms(rms);
1100
1101 const Float_t lowerlimit = (Float_t)first+1;// + fTimeLowerLimit;
1102 const Float_t upperlimit = (Float_t)last -1;// - fTimeUpperLimit;
1103
1104 // FIXME: instead of checking whether the maximum is in the first or
1105 // last extracted slice we should check whether the extractor
1106 // was able to properly extract the signal!!!
1107
1108 if (mean<lowerlimit)
1109 {
1110 *fLog << warn << hist.GetName() << ": Mean Arr.Time: "
1111 << Form("%4.1f < %4.1f (first used fadc sl+1)", mean, lowerlimit) << endl;
1112 bad.SetUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin);
1113 }
1114
1115 if (mean>upperlimit)
1116 {
1117 *fLog << warn << hist.GetName() << ": Mean Arr.Time: "
1118 << Form("%4.1f > %4.1f (last used fadc sl-1)", mean, upperlimit) << endl;
1119 bad.SetUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins);
1120 }
1121}
1122
1123// --------------------------------------------------------------------------
1124//
1125// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
1126// - MBadPixelsPix::kLoGainSaturation
1127//
1128// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
1129// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
1130// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
1131// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
1132// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
1133//
1134void MHCalibrationChargeCam::FinalizeBadPixels()
1135{
1136
1137 for (Int_t i=0; i<fBadPixels->GetSize(); i++)
1138 {
1139
1140 MBadPixelsPix &bad = (*fBadPixels)[i];
1141 MCalibrationPix &pix = (*fCam)[i];
1142
1143 if (bad.IsUncalibrated(MBadPixelsPix::kHiGainNotFitted))
1144 if (!pix.IsHiGainSaturation())
1145 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1146
1147 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainNotFitted))
1148 if (pix.IsHiGainSaturation())
1149 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1150
1151 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainSaturation))
1152 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1153
1154 if (IsOscillations())
1155 {
1156 if (bad.IsUncalibrated(MBadPixelsPix::kHiGainOscillating))
1157 if (!pix.IsHiGainSaturation())
1158 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1159
1160 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainOscillating))
1161 if (pix.IsHiGainSaturation())
1162 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1163 }
1164 }
1165}
1166
1167// --------------------------------------------------------------------------
1168//
1169// Calls MHCalibrationPix::DrawClone() for pixel idx
1170//
1171void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
1172{
1173 (*this)[idx].DrawClone();
1174}
1175
1176
1177// -----------------------------------------------------------------------------
1178//
1179// Default draw:
1180//
1181// Displays the averaged areas, both High Gain and Low Gain
1182//
1183// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1184//
1185void MHCalibrationChargeCam::Draw(const Option_t *opt)
1186{
1187
1188 const Int_t nareas = fAverageHiGainAreas->GetSize();
1189 if (nareas == 0)
1190 return;
1191
1192 TString option(opt);
1193 option.ToLower();
1194
1195 if (!option.Contains("datacheck"))
1196 {
1197 MHCalibrationCam::Draw(opt);
1198 return;
1199 }
1200
1201 //
1202 // From here on , the datacheck - Draw
1203 //
1204 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
1205 pad->SetBorderMode(0);
1206 pad->Divide(1,nareas);
1207
1208 //
1209 // Loop over inner and outer pixels
1210 //
1211 for (Int_t i=0; i<nareas;i++)
1212 {
1213 pad->cd(i+1);
1214
1215 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(i);
1216
1217 //
1218 // Ask for Hi-Gain saturation
1219 //
1220 if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
1221 {
1222 if (i>=fAverageLoGainAreas->GetSize())
1223 continue;
1224
1225 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(i);
1226 DrawDataCheckPixel(lopix,i ? fOuterRefCharge : fInnerRefCharge);
1227 }
1228 else
1229 DrawDataCheckPixel(hipix,i ? fOuterRefCharge : fInnerRefCharge);
1230 }
1231}
1232
1233// -----------------------------------------------------------------------------
1234//
1235// Draw the average pixel for the datacheck:
1236//
1237// Displays the averaged areas, both High Gain and Low Gain
1238//
1239// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1240//
1241void MHCalibrationChargeCam::DrawDataCheckPixel(MHCalibrationChargePix &pix, const Float_t refline)
1242{
1243 if (pix.IsEmpty())
1244 return;
1245
1246 TVirtualPad *pad = gPad;
1247 pad->Divide(1,2, 1e-10, 1e-10);
1248 pad->cd(1);
1249
1250 gPad->SetBorderMode(0);
1251 gPad->SetTicks();
1252 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
1253 gPad->SetLogy();
1254
1255 TH1F *hist = pix.GetHGausHist();
1256 TAxis *xaxe = hist->GetXaxis();
1257 TAxis *yaxe = hist->GetYaxis();
1258 xaxe->CenterTitle();
1259 yaxe->CenterTitle();
1260 xaxe->SetTitleSize(0.07);
1261 yaxe->SetTitleSize(0.07);
1262 xaxe->SetTitleOffset(0.7);
1263 yaxe->SetTitleOffset(0.55);
1264 xaxe->SetLabelSize(0.06);
1265 yaxe->SetLabelSize(0.06);
1266 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
1267 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
1268 xaxe->SetRange(hist->GetMaximumBin()-50, hist->GetMaximumBin()+50);
1269
1270 gStyle->SetOptFit();
1271
1272 hist->Draw();
1273
1274 TF1 *fit = pix.GetFGausFit();
1275 if (fit)
1276 {
1277 switch (fColor)
1278 {
1279 case MCalibrationCam::kGREEN:
1280 fit->SetLineColor(kGreen);
1281 break;
1282 case MCalibrationCam::kBLUE:
1283 fit->SetLineColor(kBlue);
1284 break;
1285 case MCalibrationCam::kUV:
1286 fit->SetLineColor(51);
1287 break;
1288 case MCalibrationCam::kCT1:
1289 fit->SetLineColor(6);
1290 break;
1291 default:
1292 fit->SetLineColor(kRed);
1293 }
1294 fit->Draw("same");
1295 }
1296
1297 DisplayRefLines(hist,refline);
1298
1299 pad->cd(2);
1300 gPad->SetBorderMode(0);
1301 gPad->SetTicks();
1302
1303 pix.CreateGraphEvents();
1304 TGraph *gr = pix.GetGraphEvents();
1305 if (gr)
1306 {
1307 TH1F *null2 = gr->GetHistogram();
1308
1309 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
1310 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
1311 null2->SetStats(kFALSE);
1312
1313 //
1314 // set the labels bigger
1315 //
1316 TAxis *xaxe2 = null2->GetXaxis();
1317 TAxis *yaxe2 = null2->GetYaxis();
1318 xaxe2->CenterTitle();
1319 yaxe2->CenterTitle();
1320 xaxe2->SetTitleSize(0.07);
1321 yaxe2->SetTitleSize(0.07);
1322 xaxe2->SetTitleOffset(0.7);
1323 yaxe2->SetTitleOffset(0.55);
1324 xaxe2->SetLabelSize(0.06);
1325 yaxe2->SetLabelSize(0.06);
1326 xaxe2->SetRangeUser(0, pix.GetEvents()->GetSize()/pix.GetEventFrequency());
1327 }
1328 pix.DrawEvents();
1329
1330 // newpad->cd(3);
1331 // pix.DrawPowerSpectrum(*newpad,4);
1332}
1333
1334void MHCalibrationChargeCam::DisplayRefLines(const TH1F *hist, const Float_t refline) const
1335{
1336 TLine *line = new TLine(refline, 0, refline, hist->GetMaximum());
1337 line->SetLineColor(51);
1338 line->SetLineStyle(2);
1339 line->SetLineWidth(3);
1340 line->SetBit(kCanDelete);
1341 line->Draw();
1342
1343 TLegend *leg = new TLegend(0.8,0.01,0.99,0.15);
1344 leg->AddEntry(line, "10 Leds UV", "l");
1345 leg->SetBit(kCanDelete);
1346 leg->Draw();
1347}
1348
1349Int_t MHCalibrationChargeCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1350{
1351
1352 Int_t rc = MHCalibrationCam::ReadEnv(env,prefix,print);
1353 if (rc==kERROR)
1354 return kERROR;
1355
1356 if (IsEnvDefined(env, prefix, "HiGainNbins", print))
1357 {
1358 fNbins = GetEnvValue(env, prefix, "HiGainNbins", fNbins);
1359 rc = kTRUE;
1360 }
1361
1362 if (IsEnvDefined(env, prefix, "HiGainFirst", print))
1363 {
1364 fFirst = GetEnvValue(env, prefix, "HiGainFirst", fFirst);
1365 rc = kTRUE;
1366 }
1367
1368 if (IsEnvDefined(env, prefix, "HiGainLast", print))
1369 {
1370 fLast = GetEnvValue(env, prefix, "HiGainLast", fLast);
1371 rc = kTRUE;
1372 }
1373
1374 if (IsEnvDefined(env, prefix, "LoGainNbins", print))
1375 {
1376 fLoGainNbins = GetEnvValue(env, prefix, "LoGainNbins", fLoGainNbins);
1377 rc = kTRUE;
1378 }
1379
1380 if (IsEnvDefined(env, prefix, "LoGainFirst", print))
1381 {
1382 fLoGainFirst = GetEnvValue(env, prefix, "LoGainFirst", fLoGainFirst);
1383 rc = kTRUE;
1384 }
1385
1386 if (IsEnvDefined(env, prefix, "LoGainLast", print))
1387 {
1388 fLoGainLast = GetEnvValue(env, prefix, "LoGainLast", fLoGainLast);
1389 rc = kTRUE;
1390 }
1391
1392 if (IsEnvDefined(env, prefix, "TimeLowerLimit", print))
1393 {
1394 SetTimeLowerLimit(GetEnvValue(env, prefix, "TimeLowerLimit", fTimeLowerLimit));
1395 rc = kTRUE;
1396 }
1397
1398 if (IsEnvDefined(env, prefix, "TimeUpperLimit", print))
1399 {
1400 SetTimeUpperLimit(GetEnvValue(env, prefix, "TimeUpperLimit", fTimeUpperLimit));
1401 rc = kTRUE;
1402 }
1403
1404 if (IsEnvDefined(env, prefix, "ReferenceFile", print))
1405 {
1406 SetReferenceFile(GetEnvValue(env,prefix,"ReferenceFile",fReferenceFile.Data()));
1407 rc = kTRUE;
1408 }
1409
1410 if (IsEnvDefined(env, prefix, "NumHiGainSaturationLimit", print))
1411 {
1412 SetNumHiGainSaturationLimit(GetEnvValue(env, prefix, "NumHiGainSaturationLimit", fNumHiGainSaturationLimit));
1413 rc = kTRUE;
1414 }
1415
1416 if (IsEnvDefined(env, prefix, "NumLoGainSaturationLimit", print))
1417 {
1418 SetNumLoGainSaturationLimit(GetEnvValue(env, prefix, "NumLoGainSaturationLimit", fNumLoGainSaturationLimit));
1419 rc = kTRUE;
1420 }
1421
1422 if (IsEnvDefined(env, prefix, "NumLoGainBlackoutLimit", print))
1423 {
1424 SetNumLoGainBlackoutLimit(GetEnvValue(env, prefix, "NumLoGainBlackoutLimit", fNumLoGainBlackoutLimit));
1425 rc = kTRUE;
1426 }
1427
1428 TEnv refenv(fReferenceFile);
1429
1430 fInnerRefCharge = refenv.GetValue("InnerRefCharge",fInnerRefCharge);
1431 fOuterRefCharge = refenv.GetValue("OuterRefCharge",fOuterRefCharge);
1432
1433 return rc;
1434}
Note: See TracBrowser for help on using the repository browser.