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

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