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

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