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

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