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

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