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

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