source: trunk/MagicSoft/Mars/mcalib/MHCalibrationChargeBlindPix.cc@ 3953

Last change on this file since 3953 was 3953, checked in by gaug, 21 years ago
*** empty log message ***
File size: 30.4 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// FADC counts (Blind Pixel) = 8*10^-3 * GAIN * N_phe
67//
68//////////////////////////////////////////////////////////////////////////////
69#include "MHCalibrationChargeBlindPix.h"
70
71#include <TStyle.h>
72#include <TCanvas.h>
73#include <TPaveText.h>
74
75#include <TVector.h>
76#include <TF1.h>
77#include <TH1.h>
78#include <TRandom.h>
79
80#include "MLog.h"
81#include "MLogManip.h"
82
83#include "MParList.h"
84
85#include "MRawEvtData.h"
86#include "MRawEvtPixelIter.h"
87
88#include "MExtractedSignalBlindPixel.h"
89#include "MCalibrationChargeBlindPix.h"
90
91ClassImp(MHCalibrationChargeBlindPix);
92
93using namespace std;
94
95const Double_t MHCalibrationChargeBlindPix::gkElectronicAmp = 0.008;
96const Double_t MHCalibrationChargeBlindPix::gkElectronicAmpErr = 0.002;
97
98const Int_t MHCalibrationChargeBlindPix::fgChargeNbins = 5300;
99const Axis_t MHCalibrationChargeBlindPix::fgChargeFirst = -100.5;
100const Axis_t MHCalibrationChargeBlindPix::fgChargeLast = 5199.5;
101const Float_t MHCalibrationChargeBlindPix::fgSinglePheCut = 200.;
102const Float_t MHCalibrationChargeBlindPix::fgNumSinglePheLimit = 50.;
103// --------------------------------------------------------------------------
104//
105// Default Constructor.
106//
107// Sets:
108// - the default number for fNbins (fgChargeNbins)
109// - the default number for fFirst (fgChargeFirst)
110// - the default number for fLast (fgChargeLast)
111// - the default number for fSinglePheCut (fgSingePheCut)
112// - the default number for fNumSinglePheLimit (fgNumSinglePheLimit)
113// - the default number of bins after stripping (30)
114//
115// - the default name of the fHGausHist ("HCalibrationChargeBlindPix")
116// - the default title of the fHGausHist ("Distribution of Summed FADC slices Blind Pixel ")
117// - the default x-axis title for fHGausHist ("Sum FADC Slices")
118// - the default y-axis title for fHGausHist ("Nr. of events")
119//
120// Initializes:
121// - all pointers to NULL
122// - fASinglePheFADCSlices(0);
123// - fAPedestalFADCSlices(0);
124//
125// Calls:
126// - Clear()
127//
128MHCalibrationChargeBlindPix::MHCalibrationChargeBlindPix(const char *name, const char *title)
129 : fBlindPix(NULL), fSignal(NULL), fRawEvt(NULL),
130 fSinglePheFit(NULL),
131 fFitLegend(NULL),
132 fHSinglePheFADCSlices(NULL), fHPedestalFADCSlices(NULL)
133{
134
135 fName = name ? name : "MHCalibrationChargeBlindPix";
136 fTitle = title ? title : "Statistics of the FADC sums of Blind Pixel calibration events";
137
138 SetNbins( fgChargeNbins );
139 SetFirst( fgChargeFirst );
140 SetLast ( fgChargeLast );
141
142 SetSinglePheCut();
143 SetNumSinglePheLimit();
144
145 SetBinsAfterStripping(30);
146
147 fHGausHist.SetName("HCalibrationChargeBlindPix");
148 fHGausHist.SetTitle("Distribution of Summed FADC slices Blind Pixel");
149 fHGausHist.SetXTitle("Sum FADC Slices");
150 fHGausHist.SetYTitle("Nr. of events");
151
152 fASinglePheFADCSlices(0);
153 fAPedestalFADCSlices(0);
154
155 Clear();
156}
157
158// --------------------------------------------------------------------------
159//
160// Default Destructor.
161//
162// Deletes (if Pointer is not NULL):
163//
164// - fSinglePheFit
165// - fFitLegend
166// - fHSinglePheFADCSlices
167// - fHPedestalFADCSlices
168//
169MHCalibrationChargeBlindPix::~MHCalibrationChargeBlindPix()
170{
171
172 if (fSinglePheFit)
173 delete fSinglePheFit;
174
175 if (fFitLegend)
176 delete fFitLegend;
177
178 if (fHSinglePheFADCSlices)
179 delete fHSinglePheFADCSlices;
180
181 if (fHPedestalFADCSlices)
182 delete fHPedestalFADCSlices;
183
184}
185
186// --------------------------------------------------------------------------
187//
188// Sets:
189// - all variables to 0., except the fit result variables to -999.
190// - all flags to kFALSE
191// - all pointers to NULL
192// - the default fit function (kEPoisson5)
193//
194// Deletes:
195// - all pointers unequal NULL
196//
197// Calls:
198// - MHCalibrationChargePix::Clear()
199//
200void MHCalibrationChargeBlindPix::Clear(Option_t *o)
201{
202
203 fLambda = -999.;
204 fMu0 = -999.;
205 fMu1 = -999.;
206 fSigma0 = -999.;
207 fSigma1 = -999.;
208
209 fLambdaErr = -999.;
210 fMu0Err = -999.;
211 fMu1Err = -999.;
212 fSigma0Err = -999.;
213 fSigma1Err = -999.;
214
215 fLambdaCheck = -999.;
216 fLambdaCheckErr = -999.;
217
218 fFitFunc = kEPoisson5;
219
220 fNumSinglePhes = 0;
221 fNumPedestals = 0;
222
223 fChisquare = 0.;
224 fNDF = 0 ;
225 fProb = 0.;
226
227 SetSinglePheFitOK ( kFALSE );
228 SetPedestalFitOK ( kFALSE );
229
230 if (fFitLegend)
231 {
232 delete fFitLegend;
233 fFitLegend = NULL;
234 }
235
236 if (fSinglePheFit)
237 {
238 delete fSinglePheFit;
239 fSinglePheFit = NULL;
240 }
241
242 if (fHSinglePheFADCSlices)
243 {
244 delete fHSinglePheFADCSlices;
245 fHSinglePheFADCSlices = NULL;
246 }
247
248 if (fHPedestalFADCSlices)
249 {
250 delete fHPedestalFADCSlices;
251 fHPedestalFADCSlices = NULL;
252 }
253
254
255 MHCalibrationChargePix::Clear();
256 return;
257}
258
259// --------------------------------------------------------------------------
260//
261// Set bit kSinglePheFitOK from outside
262//
263void MHCalibrationChargeBlindPix::SetSinglePheFitOK (const Bool_t b)
264{
265 b ? SETBIT(fFlags,kSinglePheFitOK) : CLRBIT(fFlags,kSinglePheFitOK);
266}
267
268// --------------------------------------------------------------------------
269//
270// Set bit kPedestalFitOK from outside
271//
272void MHCalibrationChargeBlindPix::SetPedestalFitOK(const Bool_t b)
273{
274 b ? SETBIT(fFlags,kPedestalFitOK) : CLRBIT(fFlags,kPedestalFitOK);
275}
276
277// --------------------------------------------------------------------------
278//
279// Ask for status of bit kSinglePheFitOK
280//
281const Bool_t MHCalibrationChargeBlindPix::IsSinglePheFitOK() const
282{
283 return TESTBIT(fFlags,kSinglePheFitOK);
284}
285
286// --------------------------------------------------------------------------
287//
288// Ask for status of bit kPedestalFitOK
289//
290const Bool_t MHCalibrationChargeBlindPix::IsPedestalFitOK() const
291{
292 return TESTBIT(fFlags,kPedestalFitOK);
293}
294
295// --------------------------------------------------------------------------
296//
297// Gets the pointers to:
298// - MRawEvtData
299// - MExtractedSignalBlindPixel
300//
301// Calls:
302// - MHGausHist::InitBins()
303//
304Bool_t MHCalibrationChargeBlindPix::SetupFill(const MParList *pList)
305{
306
307 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
308 if (!fRawEvt)
309 {
310 *fLog << err << "MRawEvtData not found... aborting." << endl;
311 return kFALSE;
312 }
313
314 fSignal = (MExtractedSignalBlindPixel*)pList->FindObject("MExtractedSignalBlindPixel");
315 if (!fSignal)
316 {
317 *fLog << err << "MExtractedSignalBlindPixel not found... aborting " << endl;
318 return kFALSE;
319 }
320
321 InitBins();
322
323 return kTRUE;
324}
325
326// --------------------------------------------------------------------------
327//
328// Gets or creates the pointers to:
329// - MCalibrationChargeBlindPix
330//
331Bool_t MHCalibrationChargeBlindPix::ReInit(MParList *pList)
332{
333
334 fBlindPix = (MCalibrationChargeBlindPix*)pList->FindCreateObj("MCalibrationChargeBlindPix");
335 if (!fBlindPix)
336 return kFALSE;
337
338 return kTRUE;
339}
340
341// --------------------------------------------------------------------------
342//
343// Retrieves from MExtractedSignalBlindPixel:
344// - number of FADC samples
345// - extracted signal
346// - blind Pixel ID
347//
348// Resizes (if necessary):
349// - fASinglePheFADCSlices to sum of HiGain and LoGain samples
350// - fAPedestalFADCSlices to sum of HiGain and LoGain samples
351//
352// Fills the following histograms:
353// - MHGausEvents::FillHistAndArray(signal)
354//
355// Creates MRawEvtPixelIter, jumps to blind pixel ID,
356// fills the vectors fASinglePheFADCSlices and fAPedestalFADCSlices
357// with the full FADC slices, depending on the size of the signal w.r.t. fSinglePheCut
358//
359Bool_t MHCalibrationChargeBlindPix::Fill(const MParContainer *par, const Stat_t w)
360{
361
362 const Int_t samples = (Int_t)fRawEvt->GetNumHiGainSamples()+(Int_t)fRawEvt->GetNumLoGainSamples();
363
364 if (fASinglePheFADCSlices.GetNrows() != samples)
365 {
366 fASinglePheFADCSlices.ResizeTo(samples);
367 fAPedestalFADCSlices.ResizeTo(samples);
368 }
369
370 Float_t slices = (Float_t)fSignal->GetNumFADCSamples();
371
372 if (slices == 0.)
373 {
374 *fLog << err << "Number of used signal slices in MExtractedSignalBlindPix is zero ... abort."
375 << endl;
376 return kFALSE;
377 }
378
379 //
380 // Signal extraction and histogram filling
381 //
382 const Float_t signal = (Float_t)fSignal->GetExtractedSignal();
383 FillHistAndArray(signal);
384
385 //
386 // IN order to study the single-phe posistion, we extract the slices
387 //
388 const Int_t blindpixIdx = fSignal->GetBlindPixelIdx();
389
390 MRawEvtPixelIter pixel(fRawEvt);
391 pixel.Jump(blindpixIdx);
392
393 if (signal > fSinglePheCut)
394 FillSinglePheFADCSlices(pixel);
395 else
396 FillPedestalFADCSlices(pixel);
397
398 return kTRUE;
399}
400
401// --------------------------------------------------------------------------
402//
403// Returns kFALSE, if empty
404//
405// - Creates the fourier spectrum and sets bit MHGausEvents::IsFourierSpectrumOK()
406// - Retrieves the pedestals from MExtractedSignalBlindPixel
407// - Normalizes fASinglePheFADCSlices and fAPedestalFADCSlices
408// - Executes FitPedestal()
409// - Executes FitSinglePhe()
410// - Retrieves fit results and stores them in MCalibrationChargeBlindPix
411//
412Bool_t MHCalibrationChargeBlindPix::Finalize()
413{
414
415 if (IsEmpty())
416 {
417 *fLog << err << GetDescriptor() << ": My histogram has not been filled !! " << endl;
418 return kFALSE;
419 }
420
421 CreateFourierSpectrum();
422 fBlindPix->SetOscillating ( !IsFourierSpectrumOK() );
423
424 fBlindPix->SetValid(kTRUE);
425
426 fMeanPedestal = fSignal->GetPed();
427 fMeanPedestalErr = fSignal->GetPedErr();
428 fSigmaPedestal = fSignal->GetPedRms();
429 fSigmaPedestalErr = fSignal->GetPedRmsErr();
430
431 if (fNumSinglePhes > 1)
432 for (Int_t i=0;i<fASinglePheFADCSlices.GetNrows();i++)
433 fASinglePheFADCSlices[i] = fASinglePheFADCSlices[i]/fNumSinglePhes;
434 if (fNumPedestals > 1)
435 for (Int_t i=0;i<fAPedestalFADCSlices.GetNrows();i++)
436 fAPedestalFADCSlices[i] = fAPedestalFADCSlices[i]/fNumPedestals;
437
438 FitPedestal();
439
440 if (FitSinglePhe())
441 fBlindPix->SetSinglePheFitOK();
442 else
443 fBlindPix->SetValid(kFALSE);
444
445 fBlindPix->SetLambda ( fLambda );
446 fBlindPix->SetLambdaVar ( fLambdaErr*fLambdaErr );
447 fBlindPix->SetMu0 ( fMu0 );
448 fBlindPix->SetMu0Err ( fMu0Err );
449 fBlindPix->SetMu1 ( fMu1 );
450 fBlindPix->SetMu1Err ( fMu1Err );
451 fBlindPix->SetSigma0 ( fSigma0 );
452 fBlindPix->SetSigma0Err ( fSigma0Err );
453 fBlindPix->SetSigma1 ( fSigma1 );
454 fBlindPix->SetSigma1Err ( fSigma1Err );
455 fBlindPix->SetProb ( fProb );
456
457 fBlindPix->SetLambdaCheck ( fLambdaCheck );
458 fBlindPix->SetLambdaCheckErr ( fLambdaCheckErr );
459
460 return kTRUE;
461}
462
463
464// --------------------------------------------------------------------------
465//
466// Checks again for the size and fills fASinglePheFADCSlices with the FADC slice entries
467//
468void MHCalibrationChargeBlindPix::FillSinglePheFADCSlices(const MRawEvtPixelIter &iter)
469{
470
471 const Int_t n = iter.GetNumHiGainSamples() + iter.GetNumLoGainSamples();
472
473 if (fASinglePheFADCSlices.GetNrows() < n)
474 fASinglePheFADCSlices.ResizeTo(n);
475
476 Int_t i=0;
477
478 Byte_t *start = iter.GetHiGainSamples();
479 Byte_t *end = start + iter.GetNumHiGainSamples();
480
481 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
482 fASinglePheFADCSlices(i) = fASinglePheFADCSlices(i) + (Float_t)*ptr;
483
484 start = iter.GetLoGainSamples();
485 end = start + iter.GetNumLoGainSamples();
486
487 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
488 fASinglePheFADCSlices(i) = fASinglePheFADCSlices(i) + (Float_t)*ptr;
489
490 fNumSinglePhes++;
491}
492
493// --------------------------------------------------------------------------
494//
495// Checks again for the size and fills fAPedestalFADCSlices with the FADC slice entries
496//
497void MHCalibrationChargeBlindPix::FillPedestalFADCSlices(const MRawEvtPixelIter &iter)
498{
499
500 const Int_t n = iter.GetNumHiGainSamples() + iter.GetNumLoGainSamples();
501
502 if (fAPedestalFADCSlices.GetNrows() < n)
503 fAPedestalFADCSlices.ResizeTo(n);
504
505 Int_t i = 0;
506 Byte_t *start = iter.GetHiGainSamples();
507 Byte_t *end = start + iter.GetNumHiGainSamples();
508
509 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
510 fAPedestalFADCSlices(i) = fAPedestalFADCSlices(i)+ (Float_t)*ptr;
511
512 start = iter.GetLoGainSamples();
513 end = start + iter.GetNumLoGainSamples();
514
515 for (Byte_t *ptr = start; ptr < end; ptr++, i++)
516 fAPedestalFADCSlices(i) = fAPedestalFADCSlices(i)+ (Float_t)*ptr;
517
518 fNumPedestals++;
519}
520
521
522// --------------------------------------------------------------------------
523//
524// Task to simulate single phe spectrum with the given parameters
525//
526Bool_t MHCalibrationChargeBlindPix::SimulateSinglePhe(Double_t lambda, Double_t mu0, Double_t mu1, Double_t sigma0, Double_t sigma1)
527{
528
529 gRandom->SetSeed();
530
531 if (fHGausHist.GetIntegral() != 0)
532 {
533 *fLog << err << "Histogram " << fHGausHist.GetTitle() << " is already filled. " << endl;
534 *fLog << err << "Create new class MHCalibrationBlindPixel for simulation! " << endl;
535 return kFALSE;
536 }
537
538 if (!InitFit())
539 return kFALSE;
540
541 for (Int_t i=0;i<10000; i++)
542 fHGausHist.Fill(fSinglePheFit->GetRandom());
543
544 return kTRUE;
545}
546
547// --------------------------------------------------------------------------
548//
549// - Get the ranges from the stripped histogram
550// - choose reasonable start values for the fit
551// - initialize the fit function depending on fFitFunc
552// - initialize parameter names and limits depending on fFitFunc
553//
554Bool_t MHCalibrationChargeBlindPix::InitFit()
555{
556
557 //
558 // Get the fitting ranges
559 //
560 Axis_t rmin = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetFirst());
561 Axis_t rmax = fHGausHist.GetBinCenter(fHGausHist.GetXaxis()->GetLast());
562
563 if (rmin < 0.)
564 rmin = 0.;
565
566 //
567 // First guesses for the fit (should be as close to reality as possible,
568 // otherwise the fit goes gaga because of high number of dimensions ...
569 //
570 const Stat_t entries = fHGausHist.Integral("width");
571 const Double_t lambda_guess = 0.1;
572 const Double_t maximum_bin = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
573 const Double_t norm = entries/TMath::Sqrt(TMath::TwoPi());
574
575 //
576 // Initialize the fit function
577 //
578 switch (fFitFunc)
579 {
580 case kEPoisson4:
581 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto4,rmin,rmax,6);
582 break;
583 case kEPoisson5:
584 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto5,rmin,rmax,6);
585 break;
586 case kEPoisson6:
587 fSinglePheFit = new TF1("SinglePheFit",&fPoissonKto6,rmin,rmax,6);
588 break;
589 case kEPolya:
590 fSinglePheFit = new TF1("SinglePheFit",&fPolya,rmin,rmax,8);
591 break;
592 case kEMichele:
593 break;
594
595 default:
596 *fLog << warn << "WARNING: Could not find Fit Function for Blind Pixel " << endl;
597 return kFALSE;
598 break;
599 }
600
601 if (!fSinglePheFit)
602 {
603 *fLog << warn << dbginf << "WARNING: Could not create fit function for Single Phe fit" << endl;
604 return kFALSE;
605 }
606
607 const Double_t mu_0_guess = maximum_bin;
608 const Double_t si_0_guess = 40.;
609 const Double_t mu_1_guess = mu_0_guess + 100.;
610 const Double_t si_1_guess = si_0_guess + si_0_guess;
611 // Michele
612// const Double_t lambda_1cat_guess = 0.5;
613 // const Double_t lambda_1dyn_guess = 0.5;
614 // const Double_t mu_1cat_guess = mu_0_guess + 50.;
615 // const Double_t mu_1dyn_guess = mu_0_guess + 20.;
616 // const Double_t si_1cat_guess = si_0_guess + si_0_guess;
617 // const Double_t si_1dyn_guess = si_0_guess;
618 // Polya
619 const Double_t excessPoisson_guess = 0.5;
620 const Double_t delta1_guess = 8.;
621 const Double_t delta2_guess = 5.;
622 const Double_t electronicAmp_guess = gkElectronicAmp;
623 const Double_t electronicAmp_limit = gkElectronicAmpErr;
624
625 //
626 // Initialize boundaries and start parameters
627 //
628 switch (fFitFunc)
629 {
630
631 case kEPoisson4:
632 fSinglePheFit->SetParNames( "#lambda", "#mu_{0}", "#mu_{1}", "#sigma_{0}", "#sigma_{1}","Area");
633 fSinglePheFit->SetParameters(lambda_guess,fMeanPedestal,mu_1_guess,fSigmaPedestal,si_1_guess,norm);
634
635 fSinglePheFit->SetParLimits(0,0.,0.5);
636 fSinglePheFit->SetParLimits(1,
637 fMeanPedestal-5.*fMeanPedestalErr,
638 fMeanPedestal+5.*fMeanPedestalErr);
639 fSinglePheFit->SetParLimits(2,rmin,rmax);
640 fSinglePheFit->SetParLimits(3,
641 fSigmaPedestal-5.*fSigmaPedestalErr,
642 fSigmaPedestal+5.*fSigmaPedestalErr);
643 fSinglePheFit->SetParLimits(4,0.,(rmax-rmin));
644 fSinglePheFit->SetParLimits(5,norm-(0.5*norm),norm+(0.5*norm));
645 break;
646 case kEPoisson5:
647 case kEPoisson6:
648 fSinglePheFit->SetParameters(lambda_guess,mu_0_guess,mu_1_guess,si_0_guess,si_1_guess,norm);
649 fSinglePheFit->SetParNames("#lambda","#mu_{0}","#mu_{1}","#sigma_{0}","#sigma_{1}","Area");
650 fSinglePheFit->SetParLimits(0,0.,1.);
651 fSinglePheFit->SetParLimits(1,rmin,(rmax-rmin)/1.5);
652 fSinglePheFit->SetParLimits(2,(rmax-rmin)/2.,(rmax-0.05*(rmax-rmin)));
653 fSinglePheFit->SetParLimits(3,1.0,(rmax-rmin)/2.0);
654 fSinglePheFit->SetParLimits(4,1.0,(rmax-rmin)/2.5);
655 fSinglePheFit->SetParLimits(5,norm-0.1,norm+0.1);
656 break;
657
658 case kEPolya:
659 fSinglePheFit->SetParameters(lambda_guess, excessPoisson_guess,
660 delta1_guess,delta2_guess,
661 electronicAmp_guess,
662 fSigmaPedestal,
663 norm,
664 fMeanPedestal);
665 fSinglePheFit->SetParNames("#lambda","b_{tot}",
666 "#delta_{1}","#delta_{2}",
667 "amp_{e}","#sigma_{0}",
668 "Area", "#mu_{0}");
669 fSinglePheFit->SetParLimits(0,0.,1.);
670 fSinglePheFit->SetParLimits(1,0.,1.);
671 fSinglePheFit->SetParLimits(2,6.,12.);
672 fSinglePheFit->SetParLimits(3,3.,8.);
673 fSinglePheFit->SetParLimits(4,electronicAmp_guess-electronicAmp_limit,
674 electronicAmp_guess+electronicAmp_limit);
675 fSinglePheFit->SetParLimits(5,
676 fSigmaPedestal-3.*fSigmaPedestalErr,
677 fSigmaPedestal+3.*fSigmaPedestalErr);
678 fSinglePheFit->SetParLimits(6,norm-0.1,norm+0.1);
679 fSinglePheFit->SetParLimits(7,
680 fMeanPedestal-3.*fMeanPedestalErr,
681 fMeanPedestal+3.*fMeanPedestalErr);
682 break;
683 case kEMichele:
684 break;
685
686 default:
687 *fLog << warn << "WARNING: Could not find Fit Function for Blind Pixel " << endl;
688 return kFALSE;
689 break;
690 }
691
692 fSinglePheFit->SetRange(rmin,rmax);
693
694 return kTRUE;
695}
696
697// --------------------------------------------------------------------------
698//
699// - Retrieve the parameters depending on fFitFunc
700// - Retrieve probability, Chisquare and NDF
701//
702void MHCalibrationChargeBlindPix::ExitFit()
703{
704
705
706 //
707 // Finalize
708 //
709 switch (fFitFunc)
710 {
711
712 case kEPoisson4:
713 case kEPoisson5:
714 case kEPoisson6:
715 case kEPoisson7:
716 fLambda = fSinglePheFit->GetParameter(0);
717 fMu0 = fSinglePheFit->GetParameter(1);
718 fMu1 = fSinglePheFit->GetParameter(2);
719 fSigma0 = fSinglePheFit->GetParameter(3);
720 fSigma1 = fSinglePheFit->GetParameter(4);
721
722 fLambdaErr = fSinglePheFit->GetParError(0);
723 fMu0Err = fSinglePheFit->GetParError(1);
724 fMu1Err = fSinglePheFit->GetParError(2);
725 fSigma0Err = fSinglePheFit->GetParError(3);
726 fSigma1Err = fSinglePheFit->GetParError(4);
727 break;
728 case kEPolya:
729 fLambda = fSinglePheFit->GetParameter(0);
730 fMu0 = fSinglePheFit->GetParameter(7);
731 fMu1 = 0.;
732 fSigma0 = fSinglePheFit->GetParameter(5);
733 fSigma1 = 0.;
734
735 fLambdaErr = fSinglePheFit->GetParError(0);
736 fMu0Err = fSinglePheFit->GetParError(7);
737 fMu1Err = 0.;
738 fSigma0Err = fSinglePheFit->GetParError(5);
739 fSigma1Err = 0.;
740 default:
741 break;
742 }
743
744 fProb = fSinglePheFit->GetProb();
745 fChisquare = fSinglePheFit->GetChisquare();
746 fNDF = fSinglePheFit->GetNDF();
747
748 *fLog << all << "Results of the Blind Pixel Fit: " << endl;
749 *fLog << all << "Chisquare: " << fChisquare << endl;
750 *fLog << all << "DoF: " << fNDF << endl;
751 *fLog << all << "Probability: " << fProb << endl;
752
753}
754
755// --------------------------------------------------------------------------
756//
757// - Executes InitFit()
758// - Fits the fHGausHist with fSinglePheFit
759// - Executes ExitFit()
760//
761// The fit result is accepted under condition:
762// 1) The results are not nan's
763// 2) The NDF is not smaller than fNDFLimit (5)
764// 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
765// 4) at least fNumSinglePheLimit events are in the single Photo-electron peak
766//
767Bool_t MHCalibrationChargeBlindPix::FitSinglePhe(Option_t *opt)
768{
769
770 if (!InitFit())
771 return kFALSE;
772
773 fHGausHist.Fit(fSinglePheFit,opt);
774
775 ExitFit();
776
777 //
778 // The fit result is accepted under condition:
779 // 1) The results are not nan's
780 // 2) The NDF is not smaller than fNDFLimit (5)
781 // 3) The Probability is greater than fProbLimit (default 0.001 == 99.9%)
782 // 4) at least fNumSinglePheLimit events are in the single Photo-electron peak
783 //
784 if ( TMath::IsNaN(fLambda)
785 || TMath::IsNaN(fLambdaErr)
786 || TMath::IsNaN(fProb)
787 || TMath::IsNaN(fMu0)
788 || TMath::IsNaN(fMu0Err)
789 || TMath::IsNaN(fMu1)
790 || TMath::IsNaN(fMu1Err)
791 || TMath::IsNaN(fSigma0)
792 || TMath::IsNaN(fSigma0Err)
793 || TMath::IsNaN(fSigma1)
794 || TMath::IsNaN(fSigma1Err)
795 || fNDF < fNDFLimit
796 || fProb < fProbLimit )
797 return kFALSE;
798
799 const Stat_t entries = fHGausHist.Integral("width");
800 const Float_t numSinglePhe = TMath::Exp(-1.0*fLambda)*fLambda*entries;
801
802 if (numSinglePhe < fNumSinglePheLimit)
803 {
804 *fLog << warn << "WARNING - Statistics is too low: Only " << numSinglePhe
805 << " in the Single Photo-Electron peak " << endl;
806 return kFALSE;
807 }
808 else
809 *fLog << all << numSinglePhe << " in Single Photo-Electron peak " << endl;
810
811 SetSinglePheFitOK();
812 return kTRUE;
813}
814
815// --------------------------------------------------------------------------
816//
817// - Retrieves limits for the fit
818// - Fits the fHGausHist with Gauss
819// - Retrieves the results to fLambdaCheck and fLambdaCheckErr
820// - Sets a flag IsPedestalFitOK()
821//
822void MHCalibrationChargeBlindPix::FitPedestal (Option_t *opt)
823{
824
825 // Perform the cross-check fitting only the pedestal:
826 const Axis_t rmin = 0.;
827 const Axis_t rmax = fHGausHist.GetBinCenter(fHGausHist.GetMaximumBin());
828
829 FitGaus(opt, rmin, rmax);
830
831 const Stat_t entries = fHGausHist.Integral("width");
832 const Double_t fitarea = fFGausFit->GetParameter(0);
833 const Double_t pedarea = fitarea * TMath::Sqrt(TMath::TwoPi()) * fFGausFit->GetParameter(2);
834
835 fLambdaCheck = TMath::Log(entries/pedarea);
836 fLambdaCheckErr = fFGausFit->GetParError(0)/fFGausFit->GetParameter(0)
837 + fFGausFit->GetParError(2)/fFGausFit->GetParameter(2);
838
839
840 SetPedestalFitOK();
841 return;
842}
843
844
845// -------------------------------------------------------------------------
846//
847// Draw a legend with the fit results
848//
849void MHCalibrationChargeBlindPix::DrawLegend()
850{
851
852 if (!fFitLegend)
853 {
854 fFitLegend = new TPaveText(0.05,0.05,0.95,0.95);
855 fFitLegend->SetLabel(Form("%s%s", "Results of the single PhE Fit (",
856 (fFitFunc = kEPoisson4) ? "Poisson(k=4))" :
857 (fFitFunc = kEPoisson5) ? "Poisson(k=5))" :
858 (fFitFunc = kEPoisson6) ? "Poisson(k=4))" :
859 (fFitFunc = kEPolya ) ? "Polya(k=4))" :
860 (fFitFunc = kEMichele ) ? "Michele)" : " none )" ));
861 fFitLegend->SetTextSize(0.05);
862 }
863 else
864 fFitLegend->Clear();
865
866 const TString line1 =
867 Form("Mean: #lambda = %2.2f #pm %2.2f",fLambda,fLambdaErr);
868 TText *t1 = fFitLegend->AddText(line1.Data());
869 t1->SetBit(kCanDelete);
870
871 const TString line6 =
872 Form("Mean #lambda (check) = %2.2f #pm %2.2f",fLambdaCheck,fLambdaCheckErr);
873 TText *t2 = fFitLegend->AddText(line6.Data());
874 t2->SetBit(kCanDelete);
875
876 const TString line2 =
877 Form("Pedestal: #mu_{0} = %2.2f #pm %2.2f",fMu0,fMu0Err);
878 TText *t3 = fFitLegend->AddText(line2.Data());
879 t3->SetBit(kCanDelete);
880
881 const TString line3 =
882 Form("Width Pedestal: #sigma_{0} = %2.2f #pm %2.2f",fSigma0,fSigma0Err);
883 TText *t4 = fFitLegend->AddText(line3.Data());
884 t4->SetBit(kCanDelete);
885
886 const TString line4 =
887 Form("1^{st} Phe-peak: #mu_{1} = %2.2f #pm %2.2f",fMu1,fMu1Err);
888 TText *t5 = fFitLegend->AddText(line4.Data());
889 t5->SetBit(kCanDelete);
890
891 const TString line5 =
892 Form("Width 1^{st} Phe-peak: #sigma_{1} = %2.2f #pm %2.2f",fSigma1,fSigma1Err);
893 TText *t6 = fFitLegend->AddText(line5.Data());
894 t6->SetBit(kCanDelete);
895
896 const TString line7 =
897 Form("#chi^{2} / N_{dof}: %4.2f / %3i",fChisquare,fNDF);
898 TText *t7 = fFitLegend->AddText(line7.Data());
899 t7->SetBit(kCanDelete);
900
901 const TString line8 =
902 Form("Probability: %4.2f ",fProb);
903 TText *t8 = fFitLegend->AddText(line8.Data());
904 t8->SetBit(kCanDelete);
905
906 if (IsSinglePheFitOK())
907 {
908 TText *t = fFitLegend->AddText(0.,0.,"Result of the Fit: OK");
909 t->SetBit(kCanDelete);
910 }
911 else
912 {
913 TText *t = fFitLegend->AddText("Result of the Fit: NOT OK");
914 t->SetBit(kCanDelete);
915 }
916
917 fFitLegend->SetFillColor(IsSinglePheFitOK() ? 80 : 2);
918 fFitLegend->Draw();
919
920 return;
921}
922
923
924// -------------------------------------------------------------------------
925//
926// Draw the histogram
927//
928// The following options can be chosen:
929//
930// "": displays the fHGausHist, the legend and fASinglePheFADCSlices and fAPedestalFADCSlices
931// "all": executes additionally MHGausEvents::Draw(), with option "fourierevents"
932//
933void MHCalibrationChargeBlindPix::Draw(Option_t *opt)
934{
935
936 TString option(opt);
937 option.ToLower();
938
939 Int_t win = 1;
940
941 TVirtualPad *oldpad = gPad ? gPad : MH::MakeDefCanvas(this,900, 600);
942 TVirtualPad *pad = NULL;
943
944 oldpad->SetBorderMode(0);
945
946 if (option.Contains("all"))
947 {
948 option.ReplaceAll("all","");
949 oldpad->Divide(2,1);
950 win = 2;
951 oldpad->cd(1);
952 TVirtualPad *newpad = gPad;
953 pad = newpad;
954 pad->Divide(2,2);
955 pad->cd(1);
956 }
957 else
958 {
959 pad = oldpad;
960 pad->Divide(2,2);
961 pad->cd(1);
962 }
963
964 if (!IsEmpty())
965 gPad->SetLogy();
966
967 gPad->SetTicks();
968
969 fHGausHist.Draw(opt);
970 if (fFGausFit)
971 {
972 fFGausFit->SetLineColor(kBlue);
973 fFGausFit->Draw("same");
974 }
975 if (fSinglePheFit)
976 {
977 fSinglePheFit->SetLineColor(IsSinglePheFitOK() ? kGreen : kRed);
978 fSinglePheFit->Draw("same");
979 }
980
981 pad->cd(2);
982 DrawLegend();
983
984 pad->cd(3);
985
986 if (fASinglePheFADCSlices.GetNrows() != 0)
987 {
988 if (fHSinglePheFADCSlices)
989 delete fHSinglePheFADCSlices;
990 fHSinglePheFADCSlices = new TH1F(fASinglePheFADCSlices);
991 fHSinglePheFADCSlices->SetName("SinglePheFADCSlices");
992 fHSinglePheFADCSlices->SetTitle(Form("%s%f","Assumed Single Phe FADC Slices, Sum > ",fSinglePheCut));
993 fHSinglePheFADCSlices->SetXTitle("FADC slice number");
994 fHSinglePheFADCSlices->SetYTitle("FADC counts");
995 fHSinglePheFADCSlices->Draw(opt);
996 }
997
998 pad->cd(4);
999 if (fAPedestalFADCSlices.GetNrows() != 0)
1000 {
1001
1002 if (fHPedestalFADCSlices)
1003 delete fHPedestalFADCSlices;
1004
1005 fHPedestalFADCSlices = new TH1F(fAPedestalFADCSlices);
1006 fHPedestalFADCSlices->SetName("PedestalFADCSlices");
1007 fHPedestalFADCSlices->SetTitle(Form("%s%f","Pedestal FADC Slices, Sum < ",fSinglePheCut));
1008 fHPedestalFADCSlices->SetXTitle("FADC slice number");
1009 fHPedestalFADCSlices->SetYTitle("FADC counts");
1010 fHPedestalFADCSlices->Draw(opt);
1011 }
1012
1013 if (win < 2)
1014 return;
1015
1016 oldpad->cd(2);
1017 MHGausEvents::Draw("fourierevents");
1018}
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
Note: See TracBrowser for help on using the repository browser.