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

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