source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationChargeCam.cc@ 9463

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