source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationChargeBlindPix.cc@ 4964

Last change on this file since 4964 was 4963, checked in by gaug, 20 years ago
*** empty log message ***
File size: 38.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//
27// MHCalibrationChargeBlindPix
28//
29// Histogram class for the charge calibration of the Blind Pixel.
30// Stores and fits the charges and stores the averaged assumed pedestal and
31// single-phe FADC slice entries. Charges are taken from MExtractedSignalBlindPix.
32// Performs the Single Photo-electron fit to extract the Poisson mean and its errors.
33//
34// Different fits can be chosen with the function ChangeFitFunc().
35//
36// The fit result is accepted under the condition that:
37// 1) the Probability is greater than fProbLimit (default 0.001 == 99.7%)
38// 2) at least fNumSinglePheLimit events are found in the single Photo-electron peak
39//
40// The single FADC slice entries are averaged and stored in fASinglePheFADCSlices, if
41// their sum exceeds fSinglePheCut, otherwise in fAPedestalFADCSlices.
42//
43// Used numbers are the following:
44//
45// Electronic conversion factor:
46// Assume, we have N_e electrons at the anode,
47// thus a charge of N_e*e (e = electron charge) Coulomb.
48//
49// This charge is AC coupled and runs into a R_pre = 50 Ohm resistency.
50// The corresponding current is amplified by a gain factor G_pre = 400
51// (the precision of this value still has to be checked !!!) and again AC coupled to
52// the output.
53// The corresponding signal goes through the whole transmission and
54// amplification chain and is digitized in the FADCs.
55// The conversion Signal Area to FADC counts (Conv_trans) has been measured
56// by David and Oscar to be approx. 3.9 pVs^-1
57//
58// Thus: Conversion FADC counts to Number of Electrons at Anode:
59// FADC counts = (1/Conv_tran) * G_pre * R_pre * e * N_e = 8 * 10^-4 N_e.
60//
61// Also: FADC counts = 8*10^-4 * GAIN * N_phe
62//
63// In the blind pixel, there is an additional pre-amplifier with an amplification of
64// about 10. Therefore, we have for the blind pixel:
65//
66//
67// FADC counts (Blind Pixel) = 8*10^-3 * GAIN * N_phe
68//
69//////////////////////////////////////////////////////////////////////////////
70#include "MHCalibrationChargeBlindPix.h"
71#include "MExtractBlindPixel.h"
72
73#include <TStyle.h>
74#include <TCanvas.h>
75#include <TPaveText.h>
76#include <TPaveStats.h>
77#include <TGraph.h>
78
79#include <TVector.h>
80#include <TF1.h>
81#include <TH1.h>
82#include <TH2D.h>
83#include <TRandom.h>
84
85#include "MLog.h"
86#include "MLogManip.h"
87
88#include "MParList.h"
89
90#include "MRawEvtData.h"
91#include "MRawEvtPixelIter.h"
92
93#include "MExtractedSignalBlindPixel.h"
94#include "MCalibrationChargeBlindPix.h"
95
96ClassImp(MHCalibrationChargeBlindPix);
97
98using namespace std;
99
100const Double_t MHCalibrationChargeBlindPix::gkElectronicAmp = 0.008;
101const Double_t MHCalibrationChargeBlindPix::gkElectronicAmpErr = 0.002;
102const Float_t MHCalibrationChargeBlindPix::gkSignalInitializer = -9999.;
103
104const Int_t MHCalibrationChargeBlindPix::fgChargeNbins = 128;
105const Axis_t MHCalibrationChargeBlindPix::fgChargeFirst = -0.5;
106const Axis_t MHCalibrationChargeBlindPix::fgChargeLast = 511.5;
107const Float_t MHCalibrationChargeBlindPix::fgSinglePheCut = 20.;
108const Float_t MHCalibrationChargeBlindPix::fgNumSinglePheLimit = 50.;
109// --------------------------------------------------------------------------
110//
111// Default Constructor.
112//
113// Sets:
114// - the default number for fNbins (fgChargeNbins)
115// - the default number for fFirst (fgChargeFirst)
116// - the default number for fLast (fgChargeLast)
117// - the default number for fSinglePheCut (fgSingePheCut)
118// - the default number for fNumSinglePheLimit (fgNumSinglePheLimit)
119// - the default number of bins after stripping (30)
120//
121// - the default name of the fHGausHist ("HCalibrationChargeBlindPix")
122// - the default title of the fHGausHist ("Distribution of Summed FADC slices Blind Pixel ")
123// - the default x-axis title for fHGausHist ("Sum FADC Slices")
124// - the default y-axis title for fHGausHist ("Nr. of events")
125//
126// Initializes:
127// - all pointers to NULL
128// - fASinglePheFADCSlices(0);
129// - fAPedestalFADCSlices(0);
130// - fPixId to 0
131//
132// Calls:
133// - Clear()
134//
135MHCalibrationChargeBlindPix::MHCalibrationChargeBlindPix(const char *name, const char *title)
136 : fBlindPix(NULL), fSignal(NULL), fRawEvt(NULL),
137 fSinglePheFit(NULL),
138 fFitLegend(NULL),
139 fHSinglePheFADCSlices(NULL), fHPedestalFADCSlices(NULL)
140{
141
142 fName = name ? name : "MHCalibrationChargeBlindPix";
143 fTitle = title ? title : "Statistics of the FADC sums of Blind Pixel calibration events";
144
145 SetNbins( fgChargeNbins );
146 SetFirst( fgChargeFirst );
147 SetLast ( fgChargeLast );
148
149 fASinglePheFADCSlices.ResizeTo(1);
150 fAPedestalFADCSlices.ResizeTo(1);
151
152 SetSinglePheCut();
153 SetNumSinglePheLimit();
154 SetProbLimit(0.001);
155 SetBinsAfterStripping(0);
156
157 fHGausHist.SetName("HCalibrationChargeBlindPix");
158 fHGausHist.SetTitle("Distribution of Summed FADC slices Blind Pixel");
159 fHGausHist.SetXTitle("Signal Amplitude");
160 fHGausHist.SetYTitle("Nr. of events");
161
162 fPixId = 0;
163
164 Clear();
165}
166
167// --------------------------------------------------------------------------
168//
169// Default Destructor.
170//
171// Deletes (if Pointer is not NULL):
172//
173// - fSinglePheFit
174// - fFitLegend
175// - fHSinglePheFADCSlices
176// - fHPedestalFADCSlices
177//
178MHCalibrationChargeBlindPix::~MHCalibrationChargeBlindPix()
179{
180
181 if (fSinglePheFit)
182 delete fSinglePheFit;
183
184 if (fFitLegend)
185 delete fFitLegend;
186
187 if (fHSinglePheFADCSlices)
188 delete fHSinglePheFADCSlices;
189
190 if (fHPedestalFADCSlices)
191 delete fHPedestalFADCSlices;
192
193}
194
195// --------------------------------------------------------------------------
196//
197// Sets:
198// - all variables to 0., except the fit result variables to gkSignalInitializer
199// - all flags to kFALSE
200// - all pointers to NULL
201// - the default fit function (kEPoisson5)
202//
203// Deletes:
204// - all pointers unequal NULL
205//
206// Calls:
207// - MHCalibrationChargePix::Clear()
208//
209void MHCalibrationChargeBlindPix::Clear(Option_t *o)
210{
211
212 fLambda = gkSignalInitializer;
213 fMu0 = gkSignalInitializer;
214 fMu1 = gkSignalInitializer;
215 fSigma0 = gkSignalInitializer;
216 fSigma1 = gkSignalInitializer;
217 fLambdaErr = gkSignalInitializer;
218 fMu0Err = gkSignalInitializer;
219 fMu1Err = gkSignalInitializer;
220 fSigma0Err = gkSignalInitializer;
221 fSigma1Err = gkSignalInitializer;
222
223 fLambdaCheck = gkSignalInitializer;
224 fLambdaCheckErr = gkSignalInitializer;
225
226 // fFitFunc = kEMichele;
227 fFitFunc = kEPoisson5;
228
229 fNumSinglePhes = 0;
230 fNumPedestals = 0;
231
232 fChisquare = 0.;
233 fNDF = 0 ;
234 fProb = 0.;
235
236 SetSinglePheFitOK ( kFALSE );
237 SetPedestalFitOK ( kFALSE );
238
239 if (fFitLegend)
240 {
241 delete fFitLegend;
242 fFitLegend = NULL;
243 }
244
245 if (fSinglePheFit)
246 {
247 delete fSinglePheFit;
248 fSinglePheFit = NULL;
249 }
250
251 if (fHSinglePheFADCSlices)
252 {
253 delete fHSinglePheFADCSlices;
254 fHSinglePheFADCSlices = NULL;
255 }
256
257 if (fHPedestalFADCSlices)
258 {
259 delete fHPedestalFADCSlices;
260 fHPedestalFADCSlices = NULL;
261 }
262
263 MHCalibrationPix::Clear();
264 return;
265}
266
267#if 0
268// --------------------------------------------------------------------------
269//
270// ATTENTION: This nasty Clone function is necessary since the ROOT cloning
271// lead to crashes on SOME machines (unfortunately not mine...).
272// This function is a workaround in order to achieve the correct
273// DrawClone() behaviour.
274//
275TObject *MHCalibrationChargeBlindPix::Clone(const char *name) const
276{
277
278 MHCalibrationChargeBlindPix &pix =
279 *new MHCalibrationChargeBlindPix(name ? name : fName.Data(),fTitle.Data());
280
281 //
282 // Copy MHGausEvents data members
283 //
284 pix.fBinsAfterStripping = fBinsAfterStripping;
285 pix.fCurrentSize = fCurrentSize;
286 pix.fFlags = fFlags;
287 pix.fPowerProbabilityBins = fPowerProbabilityBins;
288
289 if (fHPowerProbability)
290 pix.fHPowerProbability=(TH1I*)fHPowerProbability->Clone();
291
292 if (fPowerSpectrum)
293 pix.fPowerSpectrum = new TArrayF(*fPowerSpectrum);
294
295 pix.fEvents = fEvents;
296
297 if (fFGausFit)
298 pix.fFGausFit=(TF1*)fFGausFit->Clone();
299 if (fFExpFit)
300 pix.fFExpFit=(TF1*)fFExpFit->Clone();
301
302 pix.fFirst = fFirst;
303
304 if (fGraphEvents)
305 pix.fGraphEvents=(TGraph*)fGraphEvents->Clone();
306 if (fGraphPowerSpectrum)
307 pix.fGraphPowerSpectrum=(TGraph*)fGraphPowerSpectrum->Clone();
308
309 pix.fHGausHist = fHGausHist;
310
311 pix.fLast = fLast;
312 pix.fMean = fMean;
313 pix.fMeanErr = fMeanErr;
314 pix.fNbins = fNbins;
315 pix.fNDFLimit = fNDFLimit;
316 pix.fSigma = fSigma;
317 pix.fSigmaErr = fSigmaErr;
318 pix.fProb = fProb;
319 pix.fProbLimit = fProbLimit;
320
321 //
322 // Copy MHCalibrationPix data members
323 //
324 pix.fEventFrequency = fEventFrequency;
325 pix.fBlackoutLimit = fBlackoutLimit;
326 pix.fSaturated = fSaturated;
327 pix.fPickupLimit = fPickupLimit;
328 pix.fPixId = fPixId;
329
330 //
331 // Copy MHCalibrationChargeBlindPix data members
332 //
333 pix.fSinglePheCut = fSinglePheCut;
334 pix.fNumSinglePheLimit = fNumSinglePheLimit;
335
336 fASinglePheFADCSlices.Copy(pix.fASinglePheFADCSlices);
337 fAPedestalFADCSlices.Copy( pix.fAPedestalFADCSlices);
338
339 pix.fSinglePheFit = (TF1*)fSinglePheFit->Clone();
340
341 pix.fNumSinglePhes = fNumSinglePhes;
342 pix.fNumPedestals = fNumPedestals;
343 pix.fLambda = fLambda;
344 pix.fLambdaCheck = fLambdaCheck;
345 pix.fMu0 = fMu0;
346 pix.fMu1 = fMu1;
347 pix.fSigma0 = fSigma0;
348 pix.fSigma1 = fSigma1;
349 pix.fLambdaErr = fLambdaErr;
350 pix.fLambdaCheckErr = fLambdaCheckErr;
351 pix.fMu0Err = fMu0Err;
352 pix.fMu1Err = fMu1Err;
353 pix.fSigma0Err = fSigma0Err;
354 pix.fSigma1Err = fSigma1Err;
355 pix.fChisquare = fChisquare;
356 pix.fNDF = fNDF;
357 pix.fProb = fProb;
358 pix.fMeanPedestal = fMeanPedestal;
359 pix.fSigmaPedestal = fSigmaPedestal;
360 pix.fMeanPedestalErr = fMeanPedestalErr;
361 pix.fSigmaPedestalErr = fSigmaPedestalErr;
362
363 return &pix;
364}
365#endif
366
367/*
368// --------------------------------------------------------------------------
369//
370// Our own clone function is necessary since root 3.01/06 or Mars 0.4
371// I don't know the reason.
372//
373// Creates new MHCalibrationCam
374//
375TObject *MHCalibrationChargeBlindPix::Clone(const char *) const
376{
377
378 MHCalibrationChargeBlindPix *pix = new MHCalibrationChargeBlindPix();
379 this->Copy(*pix);
380
381 this->fHGausHist.Copy(pix->fHGausHist);
382 this->fSinglePheFit->Copy(*(pix->fSinglePheFit));
383 this->fHSinglePheFADCSlices->Copy(*(pix->fHSinglePheFADCSlices));
384 this->fHPedestalFADCSlices->Copy(*(pix->fHPedestalFADCSlices));
385
386
387 return pix;
388}
389*/
390
391// --------------------------------------------------------------------------
392//
393// Set bit kSinglePheFitOK from outside
394//
395void MHCalibrationChargeBlindPix::SetSinglePheFitOK (const Bool_t b )
396{
397 b ? SETBIT(fFlags,kSinglePheFitOK) : CLRBIT(fFlags,kSinglePheFitOK);
398}
399
400// --------------------------------------------------------------------------
401//
402// Set bit kPedestalFitOK from outside
403//
404void MHCalibrationChargeBlindPix::SetPedestalFitOK(const Bool_t b)
405{
406 b ? SETBIT(fFlags,kPedestalFitOK) : CLRBIT(fFlags,kPedestalFitOK);
407}
408
409// --------------------------------------------------------------------------
410//
411// Ask for status of bit kSinglePheFitOK
412//
413const Bool_t MHCalibrationChargeBlindPix::IsSinglePheFitOK() const
414{
415 return TESTBIT(fFlags,kSinglePheFitOK);
416}
417
418// --------------------------------------------------------------------------
419//
420// Ask for status of bit kPedestalFitOK
421//
422const Bool_t MHCalibrationChargeBlindPix::IsPedestalFitOK() const
423{
424 return TESTBIT(fFlags,kPedestalFitOK);
425}
426
427// --------------------------------------------------------------------------
428//
429// Gets the pointers to:
430// - MRawEvtData
431// - MExtractedSignalBlindPixel
432//
433Bool_t MHCalibrationChargeBlindPix::SetupFill(const MParList *pList)
434{
435
436 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
437 if (!fRawEvt)
438 {
439 *fLog << err << "MRawEvtData not found... aborting." << endl;
440 return kFALSE;
441 }
442
443 fSignal = (MExtractedSignalBlindPixel*)pList->FindObject("MExtractedSignalBlindPixel");
444 if (!fSignal)
445 {
446 *fLog << err << "MExtractedSignalBlindPixel not found... aborting " << endl;
447 return kFALSE;
448 }
449
450 return kTRUE;
451}
452
453// --------------------------------------------------------------------------
454//
455// Gets or creates the pointers to:
456// - MCalibrationChargeBlindPix
457//
458// Calls:
459// - MHGausHist::InitBins()
460//
461Bool_t MHCalibrationChargeBlindPix::ReInit(MParList *pList)
462{
463
464 fBlindPix = (MCalibrationChargeBlindPix*)pList->FindCreateObj("MCalibrationChargeBlindPix");
465 if (!fBlindPix)
466 return kFALSE;
467
468 const Int_t samples = fSignal->GetNumFADCSamples();
469 const Int_t integ = fSignal->IsExtractionType( MExtractBlindPixel::kIntegral );
470
471 //
472 // Modify the histogram size in case, integrals have been used
473 //
474 if ( fLast < samples*integ*fgChargeLast )
475 {
476 SetLast ( samples * (fgChargeLast+0.5) - 0.5 );
477 SetSinglePheCut( samples * fgSinglePheCut );
478 }
479
480 MHCalibrationPix::InitBins();
481
482 return kTRUE;
483}
484
485// --------------------------------------------------------------------------
486//
487// Retrieves from MExtractedSignalBlindPixel:
488// - number of FADC samples
489// - extracted signal
490// - blind Pixel ID
491//
492// Resizes (if necessary):
493// - fASinglePheFADCSlices to sum of HiGain and LoGain samples
494// - fAPedestalFADCSlices to sum of HiGain and LoGain samples
495//
496// Fills the following histograms:
497// - MHGausEvents::FillHistAndArray(signal)
498//
499// Creates MRawEvtPixelIter, jumps to blind pixel ID,
500// fills the vectors fASinglePheFADCSlices and fAPedestalFADCSlices
501// with the full FADC slices, depending on the size of the signal w.r.t. fSinglePheCut
502//
503Bool_t MHCalibrationChargeBlindPix::Fill(const MParContainer *par, const Stat_t w)
504{
505
506 const Int_t samples = (Int_t)fRawEvt->GetNumHiGainSamples()
507 +(Int_t)fRawEvt->GetNumLoGainSamples();
508
509 if (!fASinglePheFADCSlices.IsValid())
510 {
511 fASinglePheFADCSlices.ResizeTo(samples);
512 fAPedestalFADCSlices.ResizeTo(samples);
513 }
514
515 if (fASinglePheFADCSlices.GetNrows() != samples)
516 {
517 fASinglePheFADCSlices.ResizeTo(samples);
518 fAPedestalFADCSlices.ResizeTo(samples);
519 }
520
521 Float_t slices = (Float_t)fSignal->GetNumFADCSamples();
522
523 if (slices == 0.)
524 {
525 *fLog << err
526 << "Number of used signal slices in MExtractedSignalBlindPix "
527 << "is zero ... abort."
528 << endl;
529 return kFALSE;
530 }
531
532 //
533 // Signal extraction and histogram filling
534 //
535 const Float_t signal = fSignal->GetExtractedSignal(fPixId);
536
537 if (signal > -0.5)
538 FillHist(signal);
539 else
540 return kTRUE;
541
542 //
543 // In order to study the single-phe posistion, we extract the slices
544 //
545 const Int_t blindpixIdx = fSignal->GetBlindPixelIdx(fPixId);
546
547 MRawEvtPixelIter pixel(fRawEvt);
548 pixel.Jump(blindpixIdx);
549
550 if (signal > fSinglePheCut)
551 FillSinglePheFADCSlices(pixel);
552 else
553 FillPedestalFADCSlices(pixel);
554
555 return kTRUE;
556}
557
558// --------------------------------------------------------------------------
559//
560// Returns kFALSE, if empty
561//
562// - Creates the fourier spectrum and sets bit MHGausEvents::IsFourierSpectrumOK()
563// - Retrieves the pedestals from MExtractedSignalBlindPixel
564// - Normalizes fASinglePheFADCSlices and fAPedestalFADCSlices
565// - Executes FitPedestal()
566// - Executes FitSinglePhe()
567// - Retrieves fit results and stores them in MCalibrationChargeBlindPix
568//
569Bool_t MHCalibrationChargeBlindPix::Finalize()
570{
571
572 if (IsEmpty())
573 {
574 *fLog << err << GetDescriptor() << " ID: " << fPixId
575 << " My histogram has not been filled !! " << endl;
576 return kFALSE;
577 }
578
579 fMeanPedestal = fSignal->GetPed();
580 fMeanPedestalErr = fSignal->GetPedErr();
581 fSigmaPedestal = fSignal->GetPedRms();
582 fSigmaPedestalErr = fSignal->GetPedRmsErr();
583
584 if (fNumSinglePhes > 1)
585 for (Int_t i=0;i<fASinglePheFADCSlices.GetNrows();i++)
586 fASinglePheFADCSlices[i] = fASinglePheFADCSlices[i]/fNumSinglePhes;
587 if (fNumPedestals > 1)
588 for (Int_t i=0;i<fAPedestalFADCSlices.GetNrows();i++)
589 fAPedestalFADCSlices[i] = fAPedestalFADCSlices[i]/fNumPedestals;
590
591 FitPedestal();
592
593 fBlindPix->SetValid(kTRUE);
594
595 if (FitSinglePhe())
596 fBlindPix->SetSinglePheFitOK();
597 else
598 fBlindPix->SetValid(IsPedestalFitOK());
599
600 fBlindPix->SetLambda ( fLambda );
601 fBlindPix->SetLambdaVar ( fLambdaErr*fLambdaErr );
602 fBlindPix->SetMu0 ( fMu0 );
603 fBlindPix->SetMu0Err ( fMu0Err );
604 fBlindPix->SetMu1 ( fMu1 );
605 fBlindPix->SetMu1Err ( fMu1Err );
606 fBlindPix->SetSigma0 ( fSigma0 );
607 fBlindPix->SetSigma0Err ( fSigma0Err );
608 fBlindPix->SetSigma1 ( fSigma1 );
609 fBlindPix->SetSigma1Err ( fSigma1Err );
610 fBlindPix->SetProb ( fProb );
611
612 fBlindPix->SetLambdaCheck ( fLambdaCheck );
613 fBlindPix->SetLambdaCheckErr ( fLambdaCheckErr );
614
615 return kTRUE;
616}
617
618
619// --------------------------------------------------------------------------
620//
621// Checks again for the size and fills fASinglePheFADCSlices with the FADC slice entries
622//
623void MHCalibrationChargeBlindPix::FillSinglePheFADCSlices(const MRawEvtPixelIter &iter)
624{
625
626 const Int_t n = iter.GetNumHiGainSamples() + iter.GetNumLoGainSamples();
627
628 if (fASinglePheFADCSlices.GetNrows() < n)
629 fASinglePheFADCSlices.ResizeTo(n);
630
631 Int_t i=0;
632
633 Byte_t *start = iter.GetHiGainSamples();
634 Byte_t *end = start + iter.GetNumHiGainSamples();
635
636 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
637 fASinglePheFADCSlices(i) = fASinglePheFADCSlices(i) + (Float_t)*ptr;
638
639 start = iter.GetLoGainSamples();
640 end = start + iter.GetNumLoGainSamples();
641
642 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
643 fASinglePheFADCSlices(i) = fASinglePheFADCSlices(i) + (Float_t)*ptr;
644
645 fNumSinglePhes++;
646}
647
648// --------------------------------------------------------------------------
649//
650// Checks again for the size and fills fAPedestalFADCSlices with the FADC slice entries
651//
652void MHCalibrationChargeBlindPix::FillPedestalFADCSlices(const MRawEvtPixelIter &iter)
653{
654
655 const Int_t n = iter.GetNumHiGainSamples() + iter.GetNumLoGainSamples();
656
657 if (fAPedestalFADCSlices.GetNrows() < n)
658 fAPedestalFADCSlices.ResizeTo(n);
659
660 Int_t i = 0;
661 Byte_t *start = iter.GetHiGainSamples();
662 Byte_t *end = start + iter.GetNumHiGainSamples();
663
664 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
665 fAPedestalFADCSlices(i) = fAPedestalFADCSlices(i)+ (Float_t)*ptr;
666
667 start = iter.GetLoGainSamples();
668 end = start + iter.GetNumLoGainSamples();
669
670 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
671 fAPedestalFADCSlices(i) = fAPedestalFADCSlices(i)+ (Float_t)*ptr;
672
673 fNumPedestals++;
674}
675
676
677// --------------------------------------------------------------------------
678//
679// Task to simulate single phe spectrum with the given parameters
680//
681Bool_t MHCalibrationChargeBlindPix::SimulateSinglePhe(Double_t lambda, Double_t mu0, Double_t mu1, Double_t sigma0, Double_t sigma1)
682{
683
684 gRandom->SetSeed();
685
686 if (fHGausHist.GetIntegral() != 0)
687 {
688 *fLog << err << "Histogram " << fHGausHist.GetTitle() << " is already filled. " << endl;
689 *fLog << err << "Create new class MHCalibrationBlindPixel for simulation! " << endl;
690 return kFALSE;
691 }
692
693 if (!InitFit())
694 return kFALSE;
695
696 for (Int_t i=0;i<10000; i++)
697 fHGausHist.Fill(fSinglePheFit->GetRandom());
698
699 return kTRUE;
700}
701
702// --------------------------------------------------------------------------
703//
704// - Get the ranges from the stripped histogram
705// - choose reasonable start values for the fit
706// - initialize the fit function depending on fFitFunc
707// - initialize parameter names and limits depending on fFitFunc
708//
709Bool_t MHCalibrationChargeBlindPix::InitFit()
710{
711
712 //
713 // Get the fitting ranges
714 //
715 Axis_t rmin = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst());
716 Axis_t rmax = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast());
717
718 if (rmin < 0.)
719 rmin = 0.;
720
721 //
722 // First guesses for the fit (should be as close to reality as possible,
723 // otherwise the fit goes gaga because of high number of dimensions ...
724 //
725 const Stat_t entries = fHGausHist.Integral("width");
726 const Double_t lambda_guess = 0.5;
727 //const Double_t maximum_bin = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
728 const Double_t norm = entries/TMath::Sqrt(TMath::TwoPi());
729
730 //
731 // Initialize the fit function
732 //
733 switch (fFitFunc)
734 {
735 case kEPoisson4:
736 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto4,rmin,rmax,6);
737 rmin += 6.5;
738 break;
739 case kEPoisson5:
740 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto5,rmin,rmax,6);
741 rmin = 0.;
742 break;
743 case kEPoisson6:
744 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto6,rmin,rmax,6);
745 break;
746 case kEPolya:
747 fSinglePheFit = new TF1("SinglePheFit",&fPolya,rmin,rmax,8);
748 break;
749 case kEMichele:
750 fSinglePheFit = new TF1("SinglePheFit",&fFitFuncMichele,rmin,rmax,9);
751 break;
752 default:
753 *fLog << warn << "WARNING: Could not find Fit Function for Blind Pixel " << endl;
754 return kFALSE;
755 break;
756 }
757
758 if (!fSinglePheFit)
759 {
760 *fLog << warn << dbginf << "WARNING: Could not create fit function for Single Phe fit" << endl;
761 return kFALSE;
762 }
763
764 const Double_t mu_0_guess = 13.5;
765 const Double_t si_0_guess = 2.5;
766 const Double_t mu_1_guess = 30.;
767 const Double_t si_1_guess = si_0_guess + si_0_guess;
768 // Michele
769 const Double_t lambda_1cat_guess = 1.00;
770 const Double_t lambda_1dyn_guess = lambda_1cat_guess/10.;
771 const Double_t mu_1cat_guess = 50.;
772 const Double_t mu_1dyn_guess = 17.;
773 const Double_t si_1cat_guess = si_0_guess + si_0_guess;
774 const Double_t si_1dyn_guess = si_0_guess + si_0_guess/2.;
775 // Polya
776 const Double_t excessPoisson_guess = 0.5;
777 const Double_t delta1_guess = 8.;
778 const Double_t delta2_guess = 5.;
779 const Double_t electronicAmp_guess = gkElectronicAmp;
780 const Double_t electronicAmp_limit = gkElectronicAmpErr;
781
782 //
783 // Initialize boundaries and start parameters
784 //
785 switch (fFitFunc)
786 {
787
788 case kEPoisson4:
789 fSinglePheFit->SetParNames( "#lambda", "#mu_{0}", "#mu_{1}", "#sigma_{0}", "#sigma_{1}","Area");
790 fSinglePheFit->SetParameters(lambda_guess,mu_0_guess,mu_1_guess,si_0_guess,si_1_guess,norm);
791 fSinglePheFit->SetParLimits(0,0.,2.);
792 fSinglePheFit->SetParLimits(1,10.,17.);
793 fSinglePheFit->SetParLimits(2,17.,50.);
794 fSinglePheFit->SetParLimits(3,1.,5.);
795 fSinglePheFit->SetParLimits(4,5.,30.);
796 fSinglePheFit->SetParLimits(5,norm-(0.5*norm),norm+(0.7*norm));
797 break;
798 case kEPoisson5:
799 case kEPoisson6:
800 fSinglePheFit->SetParNames("#lambda","#mu_{0}","#mu_{1}","#sigma_{0}","#sigma_{1}","Area");
801 fSinglePheFit->SetParameters(lambda_guess,mu_0_guess,800.,si_0_guess,500.,norm);
802 fSinglePheFit->SetParLimits(0,0.,2.);
803 fSinglePheFit->SetParLimits(1,0.,100.);
804 fSinglePheFit->SetParLimits(2,300.,1500.);
805 fSinglePheFit->SetParLimits(3,30.,250.);
806 fSinglePheFit->SetParLimits(4,100.,1000.);
807 fSinglePheFit->SetParLimits(5,norm/1.5,norm*1.5);
808 break;
809
810 case kEPolya:
811 fSinglePheFit->SetParameters(lambda_guess, excessPoisson_guess,
812 delta1_guess,delta2_guess,
813 electronicAmp_guess,
814 fSigmaPedestal,
815 norm,
816 fMeanPedestal);
817 fSinglePheFit->SetParNames("#lambda","b_{tot}",
818 "#delta_{1}","#delta_{2}",
819 "amp_{e}","#sigma_{0}",
820 "Area", "#mu_{0}");
821 fSinglePheFit->SetParLimits(0,0.,1.);
822 fSinglePheFit->SetParLimits(1,0.,1.);
823 fSinglePheFit->SetParLimits(2,6.,12.);
824 fSinglePheFit->SetParLimits(3,3.,8.);
825 fSinglePheFit->SetParLimits(4,electronicAmp_guess-electronicAmp_limit,
826 electronicAmp_guess+electronicAmp_limit);
827 fSinglePheFit->SetParLimits(5,
828 fSigmaPedestal-3.*fSigmaPedestalErr,
829 fSigmaPedestal+3.*fSigmaPedestalErr);
830 fSinglePheFit->SetParLimits(6,norm-0.1,norm+0.1);
831 fSinglePheFit->SetParLimits(7,
832 fMeanPedestal-3.*fMeanPedestalErr,
833 fMeanPedestal+3.*fMeanPedestalErr);
834 break;
835 case kEMichele:
836 fSinglePheFit->SetParNames("#lambda_{cat}","#lambda_{dyn}",
837 "#mu_{0}","#mu_{1cat}","#mu_{1dyn}",
838 "#sigma_{0}","#sigma_{1cat}","#sigma_{1dyn}",
839 "Area");
840 fSinglePheFit->SetParameters(lambda_1cat_guess, lambda_1dyn_guess,
841 mu_0_guess, mu_1cat_guess,mu_1dyn_guess,
842 si_0_guess, si_1cat_guess,si_1dyn_guess,
843 norm);
844 fSinglePheFit->SetParLimits(0,0.01,2.0);
845 fSinglePheFit->SetParLimits(1,0.,0.5);
846 fSinglePheFit->SetParLimits(2,10.,16.);
847 fSinglePheFit->SetParLimits(3,25.,50.);
848 fSinglePheFit->SetParLimits(4,16.,18.5);
849 fSinglePheFit->SetParLimits(5,1.,5.);
850 fSinglePheFit->SetParLimits(6,10.,50.);
851 fSinglePheFit->SetParLimits(7,5.,10.);
852 fSinglePheFit->SetParLimits(8,norm/2.,norm*2.5);
853 break;
854
855 default:
856 *fLog << warn << "WARNING: Could not find Fit Function for Blind Pixel " << endl;
857 return kFALSE;
858 break;
859 }
860
861 fSinglePheFit->SetRange(rmin,rmax);
862
863 return kTRUE;
864}
865
866// --------------------------------------------------------------------------
867//
868// - Retrieve the parameters depending on fFitFunc
869// - Retrieve probability, Chisquare and NDF
870//
871void MHCalibrationChargeBlindPix::ExitFit()
872{
873
874
875 //
876 // Finalize
877 //
878 switch (fFitFunc)
879 {
880
881 case kEPoisson4:
882 case kEPoisson5:
883 case kEPoisson6:
884 case kEPoisson7:
885 fLambda = fSinglePheFit->GetParameter(0);
886 fMu0 = fSinglePheFit->GetParameter(1);
887 fMu1 = fSinglePheFit->GetParameter(2);
888 fSigma0 = fSinglePheFit->GetParameter(3);
889 fSigma1 = fSinglePheFit->GetParameter(4);
890
891 fLambdaErr = fSinglePheFit->GetParError(0);
892 fMu0Err = fSinglePheFit->GetParError(1);
893 fMu1Err = fSinglePheFit->GetParError(2);
894 fSigma0Err = fSinglePheFit->GetParError(3);
895 fSigma1Err = fSinglePheFit->GetParError(4);
896 break;
897 case kEPolya:
898 fLambda = fSinglePheFit->GetParameter(0);
899 fMu0 = fSinglePheFit->GetParameter(7);
900 fMu1 = 0.;
901 fSigma0 = fSinglePheFit->GetParameter(5);
902 fSigma1 = 0.;
903
904 fLambdaErr = fSinglePheFit->GetParError(0);
905 fMu0Err = fSinglePheFit->GetParError(7);
906 fMu1Err = 0.;
907 fSigma0Err = fSinglePheFit->GetParError(5);
908 fSigma1Err = 0.;
909 case kEMichele:
910 fLambda = fSinglePheFit->GetParameter(0);
911 fMu0 = fSinglePheFit->GetParameter(2);
912 fMu1 = fSinglePheFit->GetParameter(3);
913 fSigma0 = fSinglePheFit->GetParameter(5);
914 fSigma1 = fSinglePheFit->GetParameter(6);
915
916 fLambdaErr = fSinglePheFit->GetParError(0);
917 fMu0Err = fSinglePheFit->GetParError(2);
918 fMu1Err = fSinglePheFit->GetParError(3);
919 fSigma0Err = fSinglePheFit->GetParError(5);
920 fSigma1Err = fSinglePheFit->GetParError(6);
921 break;
922 default:
923 break;
924 }
925
926 fProb = fSinglePheFit->GetProb();
927 fChisquare = fSinglePheFit->GetChisquare();
928 fNDF = fSinglePheFit->GetNDF();
929
930 *fLog << all << "Results of the Blind Pixel Fit: " << endl;
931 *fLog << all << "Chisquare: " << fChisquare << endl;
932 *fLog << all << "DoF: " << fNDF << endl;
933 *fLog << all << "Probability: " << fProb << endl;
934
935}
936
937// --------------------------------------------------------------------------
938//
939// - Executes InitFit()
940// - Fits the fHGausHist with fSinglePheFit
941// - Executes ExitFit()
942//
943// The fit result is accepted under condition:
944// 1) The results are not nan's
945// 2) The NDF is not smaller than fNDFLimit (5)
946// 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
947// 4) at least fNumSinglePheLimit events are in the single Photo-electron peak
948//
949Bool_t MHCalibrationChargeBlindPix::FitSinglePhe(Option_t *opt)
950{
951
952 if (!InitFit())
953 return kFALSE;
954
955 fHGausHist.Fit(fSinglePheFit,opt);
956
957 ExitFit();
958
959 //
960 // The fit result is accepted under condition:
961 // 1) The results are not nan's
962 // 2) The NDF is not smaller than fNDFLimit (5)
963 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
964 // 4) at least fNumSinglePheLimit events are in the single Photo-electron peak
965 //
966 if ( TMath::IsNaN(fLambda)
967 || TMath::IsNaN(fLambdaErr)
968 || TMath::IsNaN(fProb)
969 || TMath::IsNaN(fMu0)
970 || TMath::IsNaN(fMu0Err)
971 || TMath::IsNaN(fMu1)
972 || TMath::IsNaN(fMu1Err)
973 || TMath::IsNaN(fSigma0)
974 || TMath::IsNaN(fSigma0Err)
975 || TMath::IsNaN(fSigma1)
976 || TMath::IsNaN(fSigma1Err)
977 || fNDF < fNDFLimit
978 || fProb < fProbLimit )
979 return kFALSE;
980
981 const Stat_t entries = fHGausHist.Integral("width");
982 const Float_t numSinglePhe = TMath::Exp(-1.0*fLambda)*fLambda*entries;
983
984 if (numSinglePhe < fNumSinglePheLimit)
985 {
986 *fLog << warn << "WARNING - Statistics is too low: Only " << numSinglePhe
987 << " in the Single Photo-Electron peak " << endl;
988 return kFALSE;
989 }
990 else
991 *fLog << all << numSinglePhe << " in Single Photo-Electron peak " << endl;
992
993 SetSinglePheFitOK();
994 return kTRUE;
995}
996
997// --------------------------------------------------------------------------
998//
999// - Retrieves limits for the fit
1000// - Fits the fHGausHist with Gauss
1001// - Retrieves the results to fLambdaCheck and fLambdaCheckErr
1002// - Sets a flag IsPedestalFitOK()
1003//
1004void MHCalibrationChargeBlindPix::FitPedestal (Option_t *opt)
1005{
1006
1007 // Perform the cross-check fitting only the pedestal:
1008 const Axis_t rmin = 0.;
1009 // const Axis_t rmax = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
1010 const Axis_t rmax = fSinglePheCut;
1011
1012 FitGaus(opt, rmin, rmax);
1013
1014 const Stat_t entries = fHGausHist.Integral("width");
1015 const Double_t pedarea = fFGausFit->Integral(0.,fSinglePheCut);
1016
1017 fLambdaCheck = TMath::Log(entries/pedarea);
1018 // estimate the error by the error of the obtained area from the Gauss-function:
1019 fLambdaCheckErr = fFGausFit->GetParError(0)/fFGausFit->GetParameter(0);
1020
1021 SetPedestalFitOK(IsGausFitOK());
1022 return;
1023}
1024
1025
1026// -------------------------------------------------------------------------
1027//
1028// Draw a legend with the fit results
1029//
1030void MHCalibrationChargeBlindPix::DrawLegend(Option_t *opt)
1031{
1032
1033 TString option(opt);
1034
1035 if (!fFitLegend)
1036 {
1037 fFitLegend = new TPaveText(0.05,0.05,0.95,0.95);
1038 fFitLegend->SetLabel(Form("%s%s", "Results of the single PhE Fit (",
1039 (fFitFunc == kEPoisson4) ? "Poisson(k=4))" :
1040 (fFitFunc == kEPoisson5) ? "Poisson(k=5))" :
1041 (fFitFunc == kEPoisson6) ? "Poisson(k=6))" :
1042 (fFitFunc == kEPolya ) ? "Polya(k=4))" :
1043 (fFitFunc == kEMichele ) ? "Michele)" : " none )" ));
1044 fFitLegend->SetTextSize(0.05);
1045 }
1046 else
1047 fFitLegend->Clear();
1048
1049 const TString line1 =
1050 Form("Mean: #lambda = %2.2f #pm %2.2f",fLambda,fLambdaErr);
1051 TText *t1 = fFitLegend->AddText(line1.Data());
1052 t1->SetBit(kCanDelete);
1053
1054 const TString line6 =
1055 Form("Mean #lambda_{check} = %2.2f #pm %2.2f",fLambdaCheck,fLambdaCheckErr);
1056 TText *t2 = fFitLegend->AddText(line6.Data());
1057 t2->SetBit(kCanDelete);
1058
1059 if (option.Contains("datacheck"))
1060 {
1061 if (fLambda + 3.*fLambdaErr < fLambdaCheck - 3.*fLambdaCheckErr
1062 ||
1063 fLambda - 3.*fLambdaErr > fLambdaCheck + 3.*fLambdaCheckErr )
1064 {
1065 TText *t = fFitLegend->AddText("#lambda and #lambda_{check} more than 3#sigma apart!");
1066 t->SetBit(kCanDelete);
1067 }
1068 }
1069 else
1070 {
1071
1072 const TString line2 =
1073 Form("Pedestal: #mu_{0} = %2.2f #pm %2.2f",fMu0,fMu0Err);
1074 TText *t3 = fFitLegend->AddText(line2.Data());
1075 t3->SetBit(kCanDelete);
1076
1077 const TString line3 =
1078 Form("Width Pedestal: #sigma_{0} = %2.2f #pm %2.2f",fSigma0,fSigma0Err);
1079 TText *t4 = fFitLegend->AddText(line3.Data());
1080 t4->SetBit(kCanDelete);
1081
1082 const TString line4 =
1083 Form("1^{st} Phe-peak: #mu_{1} = %2.2f #pm %2.2f",fMu1,fMu1Err);
1084 TText *t5 = fFitLegend->AddText(line4.Data());
1085 t5->SetBit(kCanDelete);
1086
1087 const TString line5 =
1088 Form("Width 1^{st} Phe-peak: #sigma_{1} = %2.2f #pm %2.2f",fSigma1,fSigma1Err);
1089 TText *t6 = fFitLegend->AddText(line5.Data());
1090 t6->SetBit(kCanDelete);
1091 }
1092
1093 const TString line7 =
1094 Form("#chi^{2} / N_{dof}: %4.2f / %3i",fChisquare,fNDF);
1095 TText *t7 = fFitLegend->AddText(line7.Data());
1096 t7->SetBit(kCanDelete);
1097
1098 const TString line8 =
1099 Form("Probability: %6.4f ",fProb);
1100 TText *t8 = fFitLegend->AddText(line8.Data());
1101 t8->SetBit(kCanDelete);
1102
1103 if (IsSinglePheFitOK())
1104 {
1105 TText *t = fFitLegend->AddText("Result of the Fit: OK");
1106 t->SetBit(kCanDelete);
1107 }
1108 else
1109 {
1110 TText *t = fFitLegend->AddText("Result of the Fit: NOT OK");
1111 t->SetBit(kCanDelete);
1112 }
1113
1114 fFitLegend->SetFillColor(IsSinglePheFitOK() ? 80 : 2);
1115 fFitLegend->Draw();
1116
1117 return;
1118}
1119
1120
1121// -------------------------------------------------------------------------
1122//
1123// Draw the histogram
1124//
1125// The following options can be chosen:
1126//
1127// "": displays the fHGausHist, the fits, the legend and fASinglePheFADCSlices and fAPedestalFADCSlices
1128// "all": executes additionally MHGausEvents::Draw(), with option "fourierevents"
1129// "datacheck" display the fHGausHist, the fits and the legend
1130//
1131void MHCalibrationChargeBlindPix::Draw(Option_t *opt)
1132{
1133
1134 TString option(opt);
1135 option.ToLower();
1136
1137 Int_t win = 1;
1138
1139 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,900, 600);
1140 TVirtualPad *pad = NULL;
1141
1142 if (option.Contains("all"))
1143 {
1144 option.ReplaceAll("all","");
1145 oldpad->Divide(2,1);
1146 win = 2;
1147 oldpad->cd(1);
1148 TVirtualPad *newpad = gPad;
1149 pad = newpad;
1150 pad->Divide(2,2);
1151 pad->cd(1);
1152 }
1153 else if (option.Contains("datacheck"))
1154 {
1155 pad = oldpad;
1156 pad->Divide(1,2);
1157 pad->cd(1);
1158 fHGausHist.SetStats(0);
1159 }
1160 else
1161 {
1162 pad = oldpad;
1163 pad->Divide(2,2);
1164 pad->cd(1);
1165 }
1166
1167 if (!IsEmpty() && !IsOnlyOverflow() && !IsOnlyUnderflow())
1168 gPad->SetLogy();
1169
1170 gPad->SetTicks();
1171
1172 fHGausHist.Draw();
1173 if (fFGausFit )
1174 {
1175 fFGausFit->SetLineColor(kBlue);
1176 fFGausFit->Draw("same");
1177 if (!option.Contains("datacheck"))
1178 {
1179 TLine *line = new TLine(fSinglePheCut, 0., fSinglePheCut, 10.);
1180 line->SetBit(kCanDelete);
1181 line->SetLineColor(kBlue);
1182 line->SetLineWidth(3);
1183 line->DrawLine(fSinglePheCut, 0., fSinglePheCut, 2.);
1184 }
1185 }
1186
1187 if (fSinglePheFit)
1188 {
1189 fSinglePheFit->SetFillStyle(0);
1190 fSinglePheFit->SetLineWidth(3);
1191 fSinglePheFit->SetLineColor(IsSinglePheFitOK() ? kGreen : kRed);
1192 fSinglePheFit->Draw("same");
1193 }
1194
1195 pad->cd(2);
1196 DrawLegend(option.Data());
1197
1198 if (option.Contains("datacheck"))
1199 return;
1200
1201 pad->cd(3);
1202
1203 if (fASinglePheFADCSlices.GetNrows()!=1)
1204 {
1205 if (fHSinglePheFADCSlices)
1206 delete fHSinglePheFADCSlices;
1207
1208 fHSinglePheFADCSlices = new TH1F(fASinglePheFADCSlices);
1209 fHSinglePheFADCSlices->SetName("SinglePheFADCSlices");
1210 fHSinglePheFADCSlices->SetTitle(Form("%s%4.1f","Assumed Single Phe FADC Slices, Sum > ",fSinglePheCut));
1211 fHSinglePheFADCSlices->SetXTitle("FADC slice number");
1212 fHSinglePheFADCSlices->SetYTitle("FADC counts");
1213 const Int_t nbins = fHSinglePheFADCSlices->GetNbinsX();
1214 TH2D *nulls = new TH2D("Nulls",fHSinglePheFADCSlices->GetTitle(),nbins,0.,
1215 fHSinglePheFADCSlices->GetXaxis()->GetBinCenter(nbins),
1216 100,0.,50.);
1217 nulls->SetDirectory(NULL);
1218 nulls->SetBit(kCanDelete);
1219 nulls->GetXaxis()->SetTitle(fHSinglePheFADCSlices->GetXaxis()->GetTitle());
1220 nulls->GetYaxis()->SetTitle(fHSinglePheFADCSlices->GetYaxis()->GetTitle());
1221 nulls->GetXaxis()->CenterTitle();
1222 nulls->GetYaxis()->CenterTitle();
1223 nulls->SetStats(0);
1224 nulls->Draw();
1225 fHSinglePheFADCSlices->Draw("same");
1226 }
1227
1228 pad->cd(4);
1229 if (fAPedestalFADCSlices.GetNrows()!=1)
1230 {
1231
1232 if (fHPedestalFADCSlices)
1233 delete fHPedestalFADCSlices;
1234
1235 fHPedestalFADCSlices = new TH1F(fAPedestalFADCSlices);
1236 fHPedestalFADCSlices->SetName("PedestalFADCSlices");
1237 fHPedestalFADCSlices->SetTitle(Form("%s%4.1f","Pedestal FADC Slices, Sum < ",fSinglePheCut));
1238 fHPedestalFADCSlices->SetXTitle("FADC slice number");
1239 fHPedestalFADCSlices->SetYTitle("FADC counts");
1240 const Int_t nbins = fHPedestalFADCSlices->GetNbinsX();
1241 TH2D *nullp = new TH2D("Nullp",fHPedestalFADCSlices->GetTitle(),nbins,0.,
1242 fHPedestalFADCSlices->GetXaxis()->GetBinCenter(nbins),
1243 100,0.,50.);
1244 nullp->SetDirectory(NULL);
1245 nullp->SetBit(kCanDelete);
1246 nullp->GetXaxis()->SetTitle(fHPedestalFADCSlices->GetXaxis()->GetTitle());
1247 nullp->GetYaxis()->SetTitle(fHPedestalFADCSlices->GetYaxis()->GetTitle());
1248 nullp->GetXaxis()->CenterTitle();
1249 nullp->GetYaxis()->CenterTitle();
1250 nullp->SetStats(0);
1251 nullp->Draw();
1252 fHPedestalFADCSlices->Draw("same");
1253 }
1254
1255 if (win < 2)
1256 return;
1257
1258 oldpad->cd(2);
1259 MHCalibrationPix::Draw("fourierevents");
1260}
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
Note: See TracBrowser for help on using the repository browser.