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

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