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

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