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

Last change on this file since 8286 was 8279, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 49.3 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MHCalibrationChargeCam.cc,v 1.49 2007-01-30 14:16:45 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-2006
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/////////////////////////////////////////////////////////////////////////////
117#include "MHCalibrationChargeCam.h"
118#include "MHCalibrationCam.h"
119
120#include "MLog.h"
121#include "MLogManip.h"
122
123#include "MParList.h"
124
125#include "MHCalibrationChargePix.h"
126#include "MHCalibrationPix.h"
127
128#include "MCalibrationIntensityCam.h"
129#include "MCalibrationChargeCam.h"
130#include "MCalibrationChargePix.h"
131
132#include "MGeomCam.h"
133#include "MGeomPix.h"
134
135#include "MBadPixelsIntensityCam.h"
136#include "MBadPixelsCam.h"
137#include "MBadPixelsPix.h"
138
139#include "MRawEvtData.h"
140#include "MRawRunHeader.h"
141#include "MRawEvtPixelIter.h"
142
143#include "MExtractedSignalCam.h"
144#include "MExtractedSignalPix.h"
145
146#include "MArrayI.h"
147#include "MArrayD.h"
148
149#include <TOrdCollection.h>
150#include <TPad.h>
151#include <TVirtualPad.h>
152#include <TCanvas.h>
153#include <TStyle.h>
154#include <TF1.h>
155#include <TLatex.h>
156#include <TLegend.h>
157#include <TGraph.h>
158#include <TEnv.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)
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 SetNbins(fgChargeHiGainNbins);
236 SetFirst(fgChargeHiGainFirst);
237 SetLast (fgChargeHiGainLast );
238
239 SetProbLimit(fgProbLimit);
240
241 SetLoGainNbins(fgChargeLoGainNbins);
242 SetLoGainFirst(fgChargeLoGainFirst);
243 SetLoGainLast (fgChargeLoGainLast );
244
245 SetHistName (gsHistName .Data());
246 SetHistTitle (gsHistTitle .Data());
247 SetHistXTitle(gsHistXTitle.Data());
248 SetHistYTitle(gsHistYTitle.Data());
249
250 SetAbsHistName (gsAbsHistName .Data());
251 SetAbsHistTitle (gsAbsHistTitle .Data());
252 SetAbsHistXTitle(gsAbsHistXTitle.Data());
253 SetAbsHistYTitle(gsAbsHistYTitle.Data());
254
255 SetReferenceFile();
256
257 fInnerRefCharge = 278.;
258 fOuterRefCharge = 282.;
259}
260
261// --------------------------------------------------------------------------
262//
263// Creates new MHCalibrationChargeCam only with the averaged areas:
264// the rest has to be retrieved directly, e.g. via:
265// MHCalibrationChargeCam *cam = MParList::FindObject("MHCalibrationChargeCam");
266// - cam->GetAverageSector(5).DrawClone();
267// - (*cam)[100].DrawClone()
268//
269TObject *MHCalibrationChargeCam::Clone(const char *) const
270{
271
272 MHCalibrationChargeCam *cam = new MHCalibrationChargeCam();
273
274 //
275 // Copy the data members
276 //
277 cam->fColor = fColor;
278 cam->fRunNumbers = fRunNumbers;
279 cam->fPulserFrequency = fPulserFrequency;
280 cam->fFlags = fFlags;
281 cam->fNbins = fNbins;
282 cam->fFirst = fFirst;
283 cam->fLast = fLast;
284 cam->fLoGainNbins = fLoGainNbins;
285 cam->fLoGainFirst = fLoGainFirst;
286 cam->fLoGainLast = fLoGainLast;
287 cam->fReferenceFile = fReferenceFile;
288
289 //
290 // Copy the MArrays
291 //
292 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
293 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
294 cam->fAverageAreaSat = fAverageAreaSat;
295 cam->fAverageAreaSigma = fAverageAreaSigma;
296 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
297 cam->fAverageAreaNum = fAverageAreaNum;
298 cam->fAverageSectorNum = fAverageSectorNum;
299
300 if (!IsAverageing())
301 return cam;
302
303 const Int_t navhi = fAverageHiGainAreas->GetSize();
304
305 for (int i=0; i<navhi; i++)
306 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
307
308 if (IsLoGain())
309 {
310
311 const Int_t navlo = fAverageLoGainAreas->GetSize();
312 for (int i=0; i<navlo; i++)
313 cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
314
315 }
316
317 return cam;
318}
319
320// --------------------------------------------------------------------------
321//
322// Gets the pointers to:
323// - MRawEvtData
324//
325Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
326{
327
328 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
329 if (!fRawEvt)
330 {
331 *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
332 return kFALSE;
333 }
334
335 return kTRUE;
336}
337
338// --------------------------------------------------------------------------
339//
340// Gets or creates the pointers to:
341// - MExtractedSignalCam
342// - MCalibrationChargeCam or MCalibrationIntensityChargeCam
343// - MBadPixelsCam
344//
345// Initializes the number of used FADC slices from MExtractedSignalCam
346// into MCalibrationChargeCam and test for changes in that variable
347//
348// Calls:
349// - InitHiGainArrays()
350// - InitLoGainArrays()
351//
352// Sets:
353// - fSumhiarea to nareas
354// - fSumloarea to nareas
355// - fTimehiarea to nareas
356// - fTimeloarea to nareas
357// - fSumhisector to nsectors
358// - fSumlosector to nsectors
359// - fTimehisector to nsectors
360// - fTimelosector to nsectors
361// - fSathiarea to nareas
362// - fSatloarea to nareas
363// - fSathisector to nsectors
364// - fSatlosector to nsectors
365//
366Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
367{
368
369 MExtractedSignalCam *signal =
370 (MExtractedSignalCam*)pList->FindObject(AddSerialNumber("MExtractedSignalCam"));
371 if (!signal)
372 {
373 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
374 return kFALSE;
375 }
376
377 if (!InitCams(pList,"Charge"))
378 return kFALSE;
379
380 fFirstHiGain = signal->GetFirstUsedSliceHiGain();
381 fLastHiGain = signal->GetLastUsedSliceHiGain();
382 fFirstLoGain = signal->GetFirstUsedSliceLoGain();
383 fLastLoGain = signal->GetLastUsedSliceLoGain();
384
385 const Int_t npixels = fGeom->GetNumPixels();
386 const Int_t nsectors = fGeom->GetNumSectors();
387 const Int_t nareas = fGeom->GetNumAreas();
388
389 //
390 // In case of the intense blue, double the range
391 //
392 if (fGeom->InheritsFrom("MGeomCamMagic"))
393 if ( fColor == MCalibrationCam::kBLUE)
394 SetLoGainLast(2.*fLoGainLast - fLoGainFirst);
395
396 InitHiGainArrays(npixels,nareas,nsectors);
397 InitLoGainArrays(npixels,nareas,nsectors);
398
399 fSumhiarea .Set(nareas);
400 fSumloarea .Set(nareas);
401 fTimehiarea .Set(nareas);
402 fTimeloarea .Set(nareas);
403 fSumhisector.Set(nsectors);
404 fSumlosector.Set(nsectors);
405 fTimehisector.Set(nsectors);
406 fTimelosector.Set(nsectors);
407
408 fSathiarea .Set(nareas);
409 fSatloarea .Set(nareas);
410 fSathisector.Set(nsectors);
411 fSatlosector.Set(nsectors);
412
413 return kTRUE;
414}
415
416// --------------------------------------------------------------------------
417//
418// Retrieve:
419// - fRunHeader->GetNumSamplesHiGain();
420//
421// Initializes the High Gain Arrays:
422//
423// - For every entry in the expanded arrays:
424// * Initialize an MHCalibrationChargePix
425// * Set Binning from fNbins, fFirst and fLast
426// * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
427// * Set Histgram names and titles from fHistName and fHistTitle
428// * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
429// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
430// * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
431// * Call InitHists
432//
433//
434void MHCalibrationChargeCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
435{
436
437 MBadPixelsCam *badcam = fIntensBad
438 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
439
440 TH1F *h;
441
442 const Int_t higainsamples = fRunHeader->GetNumSamplesHiGain();
443
444 if (fHiGainArray->GetSize()==0)
445 {
446 for (Int_t i=0; i<npixels; i++)
447 {
448 fHiGainArray->AddAt(new MHCalibrationChargePix(Form("%sHiGainPix%04d",fHistName.Data(),i),
449 Form("%s High Gain Pixel%04d",fHistTitle.Data(),i)),i);
450
451 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)[i];
452
453 pix.SetNbins(fNbins);
454 pix.SetFirst(fFirst);
455 pix.SetLast (fLast);
456
457 pix.SetProbLimit(fProbLimit);
458
459 pix.SetAbsTimeNbins(higainsamples);
460 pix.SetAbsTimeFirst(-0.5);
461 pix.SetAbsTimeLast(higainsamples-0.5);
462
463 InitHists(pix,(*badcam)[i],i);
464
465 h = pix.GetHAbsTime();
466
467 h->SetName (Form("H%sHiGainPix%04d",fAbsHistName.Data(),i));
468 h->SetTitle(Form("%s High Gain Pixel %04d",fAbsHistTitle.Data(),i));
469 h->SetXTitle(fAbsHistXTitle.Data());
470 h->SetYTitle(fAbsHistYTitle.Data());
471 }
472 }
473
474
475 if (fAverageHiGainAreas->GetSize()==0)
476 {
477 for (Int_t j=0; j<nareas; j++)
478 {
479 fAverageHiGainAreas->AddAt(new MHCalibrationChargePix(Form("%sHiGainArea%d",fHistName.Data(),j),
480 Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
481
482 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
483
484 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
485 pix.SetFirst(fFirst);
486 pix.SetLast (fLast);
487
488 pix.SetAbsTimeNbins(higainsamples);
489 pix.SetAbsTimeFirst(-0.5);
490 pix.SetAbsTimeLast(higainsamples-0.5);
491
492 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
493
494 h = pix.GetHAbsTime();
495
496 h->SetName (Form("H%sHiGainArea%d",fAbsHistName.Data(),j));
497 h->SetTitle(Form("%s%s%d",fAbsHistTitle.Data(),
498 " averaged on event-by-event basis High Gain Area Idx ",j));
499 h->SetXTitle(fAbsHistXTitle.Data());
500 h->SetYTitle(fAbsHistYTitle.Data());
501 }
502 }
503
504 if (fAverageHiGainSectors->GetSize()==0)
505 {
506 for (Int_t j=0; j<nsectors; j++)
507 {
508 fAverageHiGainSectors->AddAt(new MHCalibrationChargePix(Form("%sHiGainSector%02d",fHistName.Data(),j),
509 Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
510
511 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
512
513 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
514 pix.SetFirst(fFirst);
515 pix.SetLast (fLast);
516
517 pix.SetAbsTimeNbins(higainsamples);
518 pix.SetAbsTimeFirst(-0.5);
519 pix.SetAbsTimeLast(higainsamples-0.5);
520
521 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
522
523 h = pix.GetHAbsTime();
524
525 h->SetName (Form("H%sHiGainSector%02d",fAbsHistName.Data(),j));
526 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
527 " averaged on event-by-event basis High Gain Area Sector ",j));
528 h->SetXTitle(fAbsHistXTitle.Data());
529 h->SetYTitle(fAbsHistYTitle.Data());
530 }
531 }
532}
533
534//--------------------------------------------------------------------------------------
535//
536// Return, if IsLoGain() is kFALSE
537//
538// Retrieve:
539// - fRunHeader->GetNumSamplesHiGain();
540//
541// Initializes the Low Gain Arrays:
542//
543// - For every entry in the expanded arrays:
544// * Initialize an MHCalibrationChargePix
545// * Set Binning from fNbins, fFirst and fLast
546// * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
547// * Set Histgram names and titles from fHistName and fHistTitle
548// * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
549// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
550// * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
551// * Call InitHists
552//
553void MHCalibrationChargeCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
554{
555
556 if (!IsLoGain())
557 return;
558
559
560 MBadPixelsCam *badcam = fIntensBad
561 ? (MBadPixelsCam*) fIntensBad->GetCam() : fBadPixels;
562
563 const Int_t logainsamples = fRunHeader->GetNumSamplesLoGain();
564
565 TH1F *h;
566
567 if (fLoGainArray->GetSize()==0 )
568 {
569 for (Int_t i=0; i<npixels; i++)
570 {
571 fLoGainArray->AddAt(new MHCalibrationChargePix(Form("%sLoGainPix%04d",fHistName.Data(),i),
572 Form("%s Low Gain Pixel %04d",fHistTitle.Data(),i)),i);
573
574 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)(i);
575
576 pix.SetNbins(fLoGainNbins);
577 pix.SetFirst(fLoGainFirst);
578 pix.SetLast (fLoGainLast);
579
580 pix.SetProbLimit(fProbLimit);
581
582 pix.SetAbsTimeNbins(logainsamples);
583 pix.SetAbsTimeFirst(-0.5);
584 pix.SetAbsTimeLast(logainsamples-0.5);
585 pix.SetPickupLimit(fgLoGainPickupLimit);
586 pix.SetBlackoutLimit(fgLoGainBlackoutLimit);
587
588 InitHists(pix,(*badcam)[i],i);
589
590 h = pix.GetHAbsTime();
591
592 h->SetName (Form("H%sLoGainPix%04d",fAbsHistName.Data(),i));
593 h->SetTitle(Form("%s Low Gain Pixel %04d",fAbsHistTitle.Data(),i));
594 h->SetXTitle(fAbsHistXTitle.Data());
595 h->SetYTitle(fAbsHistYTitle.Data());
596 }
597 }
598
599 if (fAverageLoGainAreas->GetSize()==0)
600 {
601 for (Int_t j=0; j<nareas; j++)
602 {
603 fAverageLoGainAreas->AddAt(new MHCalibrationChargePix(Form("%sLoGainArea%d",fHistName.Data(),j),
604 Form("%s Low Gain Area Idx %d",fHistTitle.Data(),j)),j);
605
606 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
607
608 pix.SetNbins(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
609 pix.SetFirst(fLoGainFirst);
610 pix.SetLast (fLoGainLast);
611
612 pix.SetAbsTimeNbins(logainsamples);
613 pix.SetAbsTimeFirst(-0.5);
614 pix.SetAbsTimeLast(logainsamples-0.5);
615
616 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
617
618 h = pix.GetHAbsTime();
619
620 h->SetName (Form("H%sLoGainArea%02d",fAbsHistName.Data(),j));
621 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
622 " averaged on event-by-event basis Low Gain Area Idx ",j));
623 h->SetXTitle(fAbsHistXTitle.Data());
624 h->SetYTitle(fAbsHistYTitle.Data());
625 }
626 }
627
628
629 if (fAverageLoGainSectors->GetSize()==0 && IsLoGain())
630 {
631 for (Int_t j=0; j<nsectors; j++)
632 {
633 fAverageLoGainSectors->AddAt(new MHCalibrationChargePix(Form("%sLoGainSector%02d",fHistName.Data(),j),
634 Form("%s Low Gain Sector %02d",fHistTitle.Data(),j)),j);
635
636 MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
637
638 pix.SetNbins(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
639 pix.SetFirst(fLoGainFirst);
640 pix.SetLast (fLoGainLast);
641
642 pix.SetAbsTimeNbins(logainsamples);
643 pix.SetAbsTimeFirst(-0.5);
644 pix.SetAbsTimeLast(logainsamples-0.5);
645
646 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
647
648 h = pix.GetHAbsTime();
649
650 h->SetName (Form("H%sLoGainSector%02d",fAbsHistName.Data(),j));
651 h->SetTitle(Form("%s%s%02d",fAbsHistTitle.Data(),
652 " averaged on event-by-event basis Low Gain Area Sector ",j));
653 h->SetXTitle(fAbsHistXTitle.Data());
654 h->SetYTitle(fAbsHistYTitle.Data());
655 }
656 }
657}
658
659
660// --------------------------------------------------------------------------
661//
662// Retrieves from MExtractedSignalCam:
663// - first used LoGain FADC slice
664//
665// Retrieves from MGeomCam:
666// - number of pixels
667// - number of pixel areas
668// - number of sectors
669//
670// For all TOrdCollection's (including the averaged ones), the following steps are performed:
671//
672// 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
673// - MExtractedSignalPix::GetExtractedSignalHiGain();
674// - MExtractedSignalPix::GetExtractedSignalLoGain();
675//
676// 2) Set number of saturated slices (MHCalibrationChargePix::AddSaturated()) with:
677// - MExtractedSignalPix::IsHiGainSaturated();
678// - MExtractedSignalPix::GetNumLoGainSaturated();
679//
680// 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
681// - MRawEvtPixelIter::GetIdxMaxHiGainSample();
682// - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
683//
684Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
685{
686
687 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
688 if (!signal)
689 {
690 *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
691 return kFALSE;
692 }
693
694 const UInt_t npixels = fGeom->GetNumPixels();
695 const UInt_t nareas = fGeom->GetNumAreas();
696 const UInt_t nsectors = fGeom->GetNumSectors();
697 const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
698
699 fSumhiarea .Reset();
700 fSumloarea .Reset();
701 fTimehiarea .Reset();
702 fTimeloarea .Reset();
703 fSumhisector.Reset();
704 fSumlosector.Reset();
705 fTimehisector.Reset();
706 fTimelosector.Reset();
707
708 fSathiarea .Reset();
709 fSatloarea .Reset();
710 fSathisector.Reset();
711 fSatlosector.Reset();
712
713 for (UInt_t i=0; i<npixels; i++)
714 {
715 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
716 if (histhi.IsExcluded())
717 continue;
718
719 const Int_t aidx = (*fGeom)[i].GetAidx();
720 const Int_t sector = (*fGeom)[i].GetSector();
721
722 const MExtractedSignalPix &pix = (*signal)[i];
723
724 const Int_t sathi = pix.IsHiGainSaturated() ? 1 : 0;
725 histhi.AddSaturated(sathi);
726 if (sathi)
727 {
728 fSathiarea[aidx]++;
729 fSathisector[sector]++;
730 }
731
732 if (pix.IsHiGainValid())
733 {
734 const Float_t sumhi = pix.GetExtractedSignalHiGain();
735 if (IsOscillations())
736 histhi.FillHistAndArray(sumhi);
737 else
738 histhi.FillHist(sumhi);
739
740 fSumhiarea[aidx] += sumhi;
741 fSumhisector[sector] += sumhi;
742 }
743
744 if (!IsLoGain())
745 continue;
746
747 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
748
749 const Int_t satlo = pix.IsLoGainSaturated() ? 1 : 0;
750 histlo.AddSaturated(satlo);
751 if (satlo)
752 {
753 fSatloarea[aidx]++;
754 fSatlosector[sector]++;
755 }
756
757 /*
758 // Previously:
759 fSatloarea[aidx] += pix.GetNumHiGainSaturated();
760 fSatlosector[sector] += pix.GetNumHiGainSaturated();
761 */
762
763 if (pix.IsLoGainValid())
764 {
765 const Float_t sumlo = pix.GetExtractedSignalLoGain();
766
767 if (IsOscillations())
768 histlo.FillHistAndArray(sumlo);
769 else
770 histlo.FillHist(sumlo);
771
772 fSumloarea[aidx] += sumlo;
773 fSumlosector[sector] += sumlo;
774 }
775 }
776
777 MRawEvtPixelIter pixel(fRawEvt);
778 while (pixel.Next())
779 {
780 const UInt_t pixid = pixel.GetPixelId();
781
782 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
783
784 if (histhi.IsExcluded())
785 continue;
786
787 const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
788
789 histhi.FillAbsTime(timehi);
790
791 const Int_t aidx = (*fGeom)[pixid].GetAidx();
792 const Int_t sector = (*fGeom)[pixid].GetSector();
793
794 fTimehiarea [aidx] += timehi;
795 fTimehisector[sector] += timehi;
796
797 if (IsLoGain())
798 {
799 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
800
801 const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
802 histlo.FillAbsTime(timelo);
803
804 fTimeloarea[aidx] += timelo;
805 fTimelosector[sector] += timelo;
806 }
807 }
808
809 for (UInt_t j=0; j<nareas; j++)
810 {
811
812 const Int_t npix = fAverageAreaNum[j];
813
814 if (npix == 0)
815 continue;
816
817 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
818
819
820 if (IsOscillations())
821 hipix.FillHistAndArray(fSumhiarea [j]/npix);
822 else
823 hipix.FillHist(fSumhiarea[j]/npix);
824
825 hipix.AddSaturated ((Float_t)fSathiarea [j]/npix > 0.5 ? 1 : 0);
826 hipix.FillAbsTime (fTimehiarea[j]/npix);
827
828 if (IsLoGain())
829 {
830 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
831 if (IsOscillations())
832 lopix.FillHistAndArray(fSumloarea [j]/npix);
833 else
834 lopix.FillHist(fSumloarea [j]/npix);
835 lopix.AddSaturated ((Float_t)fSatloarea [j]/npix > 0.5 ? 1 : 0);
836 lopix.FillAbsTime (fTimeloarea[j]/npix);
837 }
838 }
839
840 for (UInt_t j=0; j<nsectors; j++)
841 {
842
843 const Int_t npix = fAverageSectorNum[j];
844
845 if (npix == 0)
846 continue;
847
848 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
849
850 if (IsOscillations())
851 hipix.FillHistAndArray(fSumhisector [j]/npix);
852 else
853 hipix.FillHist(fSumhisector [j]/npix);
854
855 hipix.AddSaturated ((Float_t)fSathisector[j]/npix > 0.5 ? 1 : 0);
856 hipix.FillAbsTime (fTimehisector[j]/npix);
857
858 if (IsLoGain())
859 {
860 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
861
862 if (IsOscillations())
863 lopix.FillHistAndArray(fSumlosector [j]/npix);
864 else
865 lopix.FillHist(fSumlosector [j]/npix);
866
867 lopix.AddSaturated ((Float_t)fSatlosector[j]/npix > 0.5 ? 1 : 0);
868 lopix.FillAbsTime (fTimelosector[j]/npix);
869 }
870 }
871
872 return kTRUE;
873}
874
875// --------------------------------------------------------------------------
876//
877// For all TOrdCollection's (including the averaged ones), the following steps are performed:
878//
879// 1) Returns if the pixel is excluded.
880// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
881// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
882// 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
883// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
884// otherwise the Hi-Gain ones.
885// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
886// with the flags:
887// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
888// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
889// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
890// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
891//
892Bool_t MHCalibrationChargeCam::FinalizeHists()
893{
894
895 *fLog << endl;
896
897 TH1F *h = NULL;
898
899 MCalibrationCam *chargecam = fIntensCam ? fIntensCam->GetCam() : fCam;
900 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
901
902
903 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
904 {
905
906 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
907
908 if (histhi.IsExcluded())
909 continue;
910
911 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i] ;
912
913 const Int_t numsat = histhi.GetSaturated();
914
915 pix.SetNumSaturated(numsat);
916
917 if (numsat > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
918 {
919 pix.SetHiGainSaturation();
920 if (IsOscillations())
921 histhi.CreateFourierSpectrum();
922 continue;
923 }
924
925 MBadPixelsPix &bad = (*badcam)[i];
926
927 h = histhi.GetHGausHist();
928
929 Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
930 if (overflow > fOverflowLimit*histhi.GetHGausHist()->GetEntries())
931 {
932 *fLog << warn
933 << "HiGain Hist-overflow " << overflow
934 << " times in: " << histhi.GetName() << " (w/o saturation!) " << endl;
935 bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
936 }
937
938 overflow = h->GetBinContent(0);
939 if (overflow > fOverflowLimit*histhi.GetHGausHist()->GetEntries())
940 {
941 *fLog << warn
942 << "HiGain Hist-underflow " << overflow
943 << " times in pix: " << histhi.GetName() << " (w/o saturation!) " << endl;
944 bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
945 }
946
947 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
948 }
949
950 if (IsLoGain())
951 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
952 {
953
954 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
955 MBadPixelsPix &bad = (*badcam)[i];
956
957 if (histlo.IsExcluded())
958 continue;
959
960 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
961 {
962 *fLog << warn << "Pixel " << setw(4) << i << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
963 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
964 if (IsOscillations())
965 histlo.CreateFourierSpectrum();
966 continue;
967 }
968
969 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i] ;
970 if (!pix.IsHiGainSaturation())
971 continue;
972
973 h = histlo.GetHGausHist();
974
975 Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
976 if (overflow > fOverflowLimit*histlo.GetHGausHist()->GetEntries())
977 {
978 *fLog << warn
979 << "LoGain Hist-overflow " << overflow
980 << " times in: " << histlo.GetName() << " (w/o saturation!) " << endl;
981 bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
982 }
983
984 overflow = h->GetBinContent(0);
985 if (overflow > fOverflowLimit*histlo.GetHGausHist()->GetEntries())
986 {
987 *fLog << warn
988 << "LoGain Hist-underflow " << overflow
989 << " times in: " << histlo.GetName() << " (w/o saturation!) " << endl;
990 bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
991 }
992
993 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
994 }
995
996 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
997 {
998
999 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
1000 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(j);
1001
1002 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
1003 {
1004 pix.SetHiGainSaturation();
1005 if (IsOscillations())
1006 histhi.CreateFourierSpectrum();
1007 continue;
1008 }
1009
1010 MBadPixelsPix &bad = chargecam->GetAverageBadArea(j);
1011 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
1012 }
1013
1014 if (IsLoGain())
1015 for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
1016 {
1017
1018 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
1019
1020 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
1021 {
1022 *fLog << warn << "Area " << setw(4) << j << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
1023 if (IsOscillations())
1024 histlo.CreateFourierSpectrum();
1025 continue;
1026 }
1027
1028 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(j);
1029
1030 if (pix.IsHiGainSaturation())
1031 {
1032 MBadPixelsPix &bad = chargecam->GetAverageBadArea(j);
1033 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
1034 }
1035
1036 }
1037
1038 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
1039 {
1040
1041 MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
1042 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(j);
1043
1044 if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
1045 {
1046 pix.SetHiGainSaturation();
1047 if (IsOscillations())
1048 histhi.CreateFourierSpectrum();
1049 continue;
1050 }
1051
1052 MBadPixelsPix &bad = chargecam->GetAverageBadSector(j);
1053 FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
1054 }
1055
1056 if (IsLoGain())
1057 for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
1058 {
1059
1060 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
1061 MBadPixelsPix &bad = chargecam->GetAverageBadSector(j);
1062
1063 if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
1064 {
1065 *fLog << warn << "Sector " << setw(4) << j << ": More than " << Form("%.1f%%", fNumLoGainSaturationLimit*100) << " lo-gains saturated." << endl;
1066 bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
1067 if (IsOscillations())
1068 histlo.CreateFourierSpectrum();
1069 continue;
1070 }
1071
1072 MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(j);
1073
1074 if (pix.IsHiGainSaturation())
1075 FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
1076 }
1077
1078 //
1079 // Perform the fitting for the High Gain (done in MHCalibrationCam)
1080 //
1081 FitHiGainArrays(*chargecam, *badcam,
1082 MBadPixelsPix::kHiGainNotFitted,
1083 MBadPixelsPix::kHiGainOscillating);
1084
1085 //
1086 // Perform the fitting for the Low Gain (done in MHCalibrationCam)
1087 //
1088 if (IsLoGain())
1089 FitLoGainArrays(*chargecam, *badcam,
1090 MBadPixelsPix::kLoGainNotFitted,
1091 MBadPixelsPix::kLoGainOscillating);
1092
1093 //
1094 // Check for pixels which have the high-gain saturated, but the low-gain
1095 // switch not applied in sufficient cases. Have to exclude these pixels,
1096 // although they not abnormal. They simply cannot be calibrated...
1097 //
1098 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
1099 {
1100
1101 MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
1102 if (histlo.IsExcluded())
1103 continue;
1104
1105 MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i] ;
1106 if (!pix.IsHiGainSaturation())
1107 continue;
1108
1109 //
1110 // Now,treat only low-gain calibrated pixels:
1111 //
1112 const Double_t lim = fNumLoGainBlackoutLimit*histlo.GetHGausHist()->GetEntries();
1113 if (histlo.GetBlackout() <= lim)
1114 continue;
1115
1116 MBadPixelsPix &bad = (*badcam)[i];
1117 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1118 bad.SetUncalibrated(MBadPixelsPix::kLoGainBlackout);
1119 }
1120
1121 return kTRUE;
1122}
1123
1124// --------------------------------------------------------------------------------
1125//
1126// Fill the absolute time results into MCalibrationChargePix
1127//
1128// Check absolute time validity:
1129// - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
1130// - Mean arrival time is at least fUpperLimit slices from the upper edge
1131//
1132void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
1133 Byte_t first, Byte_t last)
1134{
1135 const Float_t mean = hist.GetAbsTimeMean();
1136 const Float_t rms = hist.GetAbsTimeRms();
1137
1138 pix.SetAbsTimeMean(mean);
1139 pix.SetAbsTimeRms(rms);
1140
1141 const Float_t lowerlimit = (Float_t)first;// + fTimeLowerLimit;
1142 const Float_t upperlimit = (Float_t)last ;// - fTimeUpperLimit;
1143
1144 // FIXME: instead of checking whether the maximum is in the first or
1145 // last extracted slice we should check whether the extractor
1146 // was able to properly extract the signal!!!
1147
1148 if (mean<lowerlimit)
1149 {
1150 *fLog << warn << hist.GetName() << ": Mean Arr.Time: "
1151 << Form("%4.1f < %4.1f, (%3.1f + %3.1)", mean, lowerlimit, fTimeLowerLimit, first) << endl;
1152 bad.SetUncalibrated(MBadPixelsPix::kMeanTimeInFirstBin);
1153 }
1154
1155 if (mean>upperlimit)
1156 {
1157 *fLog << warn << hist.GetName() << ": Mean Arr.Time: "
1158 << Form("%4.1f > %4.1f, (%3.1f - %3.1)", mean, upperlimit, fTimeUpperLimit, last) << endl;
1159 bad.SetUncalibrated(MBadPixelsPix::kMeanTimeInLast2Bins);
1160 }
1161}
1162
1163// --------------------------------------------------------------------------
1164//
1165// Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
1166// - MBadPixelsPix::kLoGainSaturation
1167//
1168// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
1169// - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
1170// - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
1171// - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
1172// - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
1173//
1174void MHCalibrationChargeCam::FinalizeBadPixels()
1175{
1176
1177 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
1178 MCalibrationCam *chargecam = fIntensCam ? fIntensCam->GetCam() : fCam;
1179
1180 for (Int_t i=0; i<badcam->GetSize(); i++)
1181 {
1182
1183 MBadPixelsPix &bad = (*badcam)[i];
1184 MCalibrationPix &pix = (*chargecam)[i];
1185
1186 if (bad.IsUncalibrated(MBadPixelsPix::kHiGainNotFitted))
1187 if (!pix.IsHiGainSaturation())
1188 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1189
1190 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainNotFitted))
1191 if (pix.IsHiGainSaturation())
1192 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1193
1194 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainSaturation))
1195 bad.SetUnsuitable(MBadPixelsPix::kUnsuitableRun);
1196
1197 if (IsOscillations())
1198 {
1199 if (bad.IsUncalibrated(MBadPixelsPix::kHiGainOscillating))
1200 if (!pix.IsHiGainSaturation())
1201 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1202
1203 if (bad.IsUncalibrated(MBadPixelsPix::kLoGainOscillating))
1204 if (pix.IsHiGainSaturation())
1205 bad.SetUnsuitable(MBadPixelsPix::kUnreliableRun);
1206 }
1207 }
1208}
1209
1210// --------------------------------------------------------------------------
1211//
1212// Calls MHCalibrationPix::DrawClone() for pixel idx
1213//
1214void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
1215{
1216 (*this)[idx].DrawClone();
1217}
1218
1219
1220// -----------------------------------------------------------------------------
1221//
1222// Default draw:
1223//
1224// Displays the averaged areas, both High Gain and Low Gain
1225//
1226// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1227//
1228void MHCalibrationChargeCam::Draw(const Option_t *opt)
1229{
1230
1231 const Int_t nareas = fAverageHiGainAreas->GetSize();
1232 if (nareas == 0)
1233 return;
1234
1235 TString option(opt);
1236 option.ToLower();
1237
1238 if (!option.Contains("datacheck"))
1239 {
1240 MHCalibrationCam::Draw(opt);
1241 return;
1242 }
1243
1244 //
1245 // From here on , the datacheck - Draw
1246 //
1247 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
1248 pad->SetBorderMode(0);
1249 pad->Divide(1,nareas);
1250
1251 //
1252 // Loop over inner and outer pixels
1253 //
1254 for (Int_t i=0; i<nareas;i++)
1255 {
1256
1257 pad->cd(i+1);
1258
1259 MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(i);
1260 //
1261 // Ask for Hi-Gain saturation
1262 //
1263 if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
1264 {
1265 MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(i);
1266 DrawDataCheckPixel(lopix,i ? fOuterRefCharge : fInnerRefCharge);
1267 }
1268 else
1269 DrawDataCheckPixel(hipix,i ? fOuterRefCharge : fInnerRefCharge);
1270
1271 }
1272}
1273
1274// -----------------------------------------------------------------------------
1275//
1276// Draw the average pixel for the datacheck:
1277//
1278// Displays the averaged areas, both High Gain and Low Gain
1279//
1280// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
1281//
1282void MHCalibrationChargeCam::DrawDataCheckPixel(MHCalibrationChargePix &pix, const Float_t refline)
1283{
1284 if (pix.IsEmpty())
1285 return;
1286
1287 TVirtualPad *newpad = gPad;
1288 newpad->Divide(1,2);
1289 newpad->cd(1);
1290
1291 gPad->SetTicks();
1292 if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
1293 gPad->SetLogy();
1294
1295 TH1F *hist = pix.GetHGausHist();
1296
1297 TH1F *null = new TH1F("Null",hist->GetTitle(),100,
1298 pix.GetFirst() > 0. ? pix.GetFirst() : 0.,
1299 pix.GetLast() > pix.GetFirst()
1300 ? ( pix.GetLast() > 450.
1301 ? ( fColor == MCalibrationCam::kBLUE ? 800. : 450. )
1302 : pix.GetLast() )
1303 : pix.GetFirst()*2.);
1304
1305 null->SetMaximum(1.1*hist->GetMaximum());
1306 null->SetDirectory(NULL);
1307 null->SetBit(kCanDelete);
1308 null->SetStats(kFALSE);
1309 //
1310 // set the labels bigger
1311 //
1312 TAxis *xaxe = null->GetXaxis();
1313 TAxis *yaxe = null->GetYaxis();
1314 xaxe->CenterTitle();
1315 yaxe->CenterTitle();
1316 xaxe->SetTitleSize(0.07);
1317 yaxe->SetTitleSize(0.07);
1318 xaxe->SetTitleOffset(0.7);
1319 yaxe->SetTitleOffset(0.55);
1320 xaxe->SetLabelSize(0.06);
1321 yaxe->SetLabelSize(0.06);
1322 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
1323 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
1324
1325 null->Draw();
1326 hist->Draw("same");
1327
1328 gStyle->SetOptFit();
1329
1330 TF1 *fit = pix.GetFGausFit();
1331
1332 if (fit)
1333 {
1334 switch ( fColor )
1335 {
1336 case MCalibrationCam::kGREEN:
1337 fit->SetLineColor(kGreen);
1338 break;
1339 case MCalibrationCam::kBLUE:
1340 fit->SetLineColor(kBlue);
1341 break;
1342 case MCalibrationCam::kUV:
1343 fit->SetLineColor(106);
1344 break;
1345 case MCalibrationCam::kCT1:
1346 fit->SetLineColor(006);
1347 break;
1348 default:
1349 fit->SetLineColor(kRed);
1350 }
1351 fit->Draw("same");
1352 }
1353
1354 DisplayRefLines(null,refline);
1355
1356 newpad->cd(2);
1357 gPad->SetTicks();
1358
1359 TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
1360
1361 null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
1362 null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
1363 null2->SetDirectory(NULL);
1364 null2->SetBit(kCanDelete);
1365 null2->SetStats(kFALSE);
1366 //
1367 // set the labels bigger
1368 //
1369 TAxis *xaxe2 = null2->GetXaxis();
1370 TAxis *yaxe2 = null2->GetYaxis();
1371 xaxe2->CenterTitle();
1372 yaxe2->CenterTitle();
1373 xaxe2->SetTitleSize(0.07);
1374 yaxe2->SetTitleSize(0.07);
1375 xaxe2->SetTitleOffset(0.7);
1376 yaxe2->SetTitleOffset(0.55);
1377 xaxe2->SetLabelSize(0.06);
1378 yaxe2->SetLabelSize(0.06);
1379
1380 pix.CreateGraphEvents();
1381 TGraph *gr = pix.GetGraphEvents();
1382 if (gr)
1383 {
1384 xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
1385 yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
1386 }
1387
1388 null2->Draw();
1389
1390 pix.DrawEvents("same");
1391
1392 // newpad->cd(3);
1393 // pix.DrawPowerSpectrum(*newpad,4);
1394
1395 return;
1396
1397}
1398
1399
1400void MHCalibrationChargeCam::DisplayRefLines(const TH1F *hist, const Float_t refline) const
1401{
1402
1403 TGraph *uv10 = new TGraph(2);
1404 uv10->SetPoint(0,refline,0.);
1405 uv10->SetPoint(1,refline,hist->GetMaximum());
1406 uv10->SetBit(kCanDelete);
1407 uv10->SetLineColor(106);
1408 uv10->SetLineStyle(2);
1409 uv10->SetLineWidth(3);
1410 uv10->Draw("L");
1411
1412 TLegend *leg = new TLegend(0.8,0.55,0.99,0.99);
1413 leg->SetBit(kCanDelete);
1414 leg->AddEntry(uv10,"10 Leds UV","l");
1415
1416 leg->Draw();
1417}
1418
1419Int_t MHCalibrationChargeCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1420{
1421
1422 Bool_t rc = kFALSE;
1423
1424 if (MHCalibrationCam::ReadEnv(env,prefix,print))
1425 rc = kTRUE;
1426
1427 if (IsEnvDefined(env, prefix, "HiGainNbins", print))
1428 {
1429 SetNbins(GetEnvValue(env, prefix, "HiGainNbins", fNbins));
1430 rc = kTRUE;
1431 }
1432
1433 if (IsEnvDefined(env, prefix, "HiGainFirst", print))
1434 {
1435 SetFirst(GetEnvValue(env, prefix, "HiGainFirst", fFirst));
1436 rc = kTRUE;
1437 }
1438
1439 if (IsEnvDefined(env, prefix, "HiGainLast", print))
1440 {
1441 SetLast(GetEnvValue(env, prefix, "HiGainLast", fLast));
1442 rc = kTRUE;
1443 }
1444
1445 if (IsEnvDefined(env, prefix, "LoGainNbins", print))
1446 {
1447 SetLoGainNbins(GetEnvValue(env, prefix, "LoGainNbins", fLoGainNbins));
1448 rc = kTRUE;
1449 }
1450
1451 if (IsEnvDefined(env, prefix, "LoGainFirst", print))
1452 {
1453 SetLoGainFirst(GetEnvValue(env, prefix, "LoGainFirst", fLoGainFirst));
1454 rc = kTRUE;
1455 }
1456
1457 if (IsEnvDefined(env, prefix, "LoGainLast", print))
1458 {
1459 SetLoGainLast(GetEnvValue(env, prefix, "LoGainLast", fLoGainLast));
1460 rc = kTRUE;
1461 }
1462
1463 if (IsEnvDefined(env, prefix, "TimeLowerLimit", print))
1464 {
1465 SetTimeLowerLimit(GetEnvValue(env, prefix, "TimeLowerLimit", fTimeLowerLimit));
1466 rc = kTRUE;
1467 }
1468
1469 if (IsEnvDefined(env, prefix, "TimeUpperLimit", print))
1470 {
1471 SetTimeUpperLimit(GetEnvValue(env, prefix, "TimeUpperLimit", fTimeUpperLimit));
1472 rc = kTRUE;
1473 }
1474
1475 if (IsEnvDefined(env, prefix, "ReferenceFile", print))
1476 {
1477 SetReferenceFile(GetEnvValue(env,prefix,"ReferenceFile",fReferenceFile.Data()));
1478 rc = kTRUE;
1479 }
1480
1481 if (IsEnvDefined(env, prefix, "NumHiGainSaturationLimit", print))
1482 {
1483 SetNumHiGainSaturationLimit(GetEnvValue(env, prefix, "NumHiGainSaturationLimit", fNumHiGainSaturationLimit));
1484 rc = kTRUE;
1485 }
1486
1487 if (IsEnvDefined(env, prefix, "NumLoGainSaturationLimit", print))
1488 {
1489 SetNumLoGainSaturationLimit(GetEnvValue(env, prefix, "NumLoGainSaturationLimit", fNumLoGainSaturationLimit));
1490 rc = kTRUE;
1491 }
1492
1493 if (IsEnvDefined(env, prefix, "NumLoGainBlackoutLimit", print))
1494 {
1495 SetNumLoGainBlackoutLimit(GetEnvValue(env, prefix, "NumLoGainBlackoutLimit", fNumLoGainBlackoutLimit));
1496 rc = kTRUE;
1497 }
1498
1499
1500 TEnv refenv(fReferenceFile);
1501
1502 fInnerRefCharge = refenv.GetValue("InnerRefCharge",fInnerRefCharge);
1503 fOuterRefCharge = refenv.GetValue("OuterRefCharge",fOuterRefCharge);
1504
1505 return rc;
1506}
Note: See TracBrowser for help on using the repository browser.