]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/macros/FitCDFLocal.C
including switch to set on/off iso-track core removal, cleaning and bug fix
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / macros / FitCDFLocal.C
1 TH1F *psproperMCsecJPSI();
2 void LoadLib();
3 char * inputDistr = "data/result_data_70any.root";
4
5 Double_t *resParamFF = new Double_t[9];
6 Double_t *resParamFS = new Double_t[9];
7 Double_t *resParamSS = new Double_t[9];
8 Double_t *invMassParam = new Double_t[9];
9 Double_t *bkgParam = new Double_t[10];
10 Double_t Fb = 0.1;
11 Double_t Fsig = 0.291701; // 2.4 - 4
12
13 Bool_t range = 0; // 1 signal (2.92 - 3.16) 
14                   // 0 all (2.4 - 4) 
15
16 //select types of candidates used for likelihood fit
17 TString resType = "FF;FS;SS";
18
19 Double_t weightType[] = {0.,0.,0.};
20
21 void FitCDFLocal();
22 AliDielectronBtoJPSItoEleCDFfitFCN *likely_obj = 0x0;
23
24 void  FitCDFLocal(){
25
26   ///////////////////////////////////////////////////////////////////
27   //
28   // Example macro to read local N-Tuples of JPSI 
29   // and bkg candidates and perform log-likelihood 
30   // minimization using the minimization handler class
31   //
32   // Origin: C. Di Giglio
33   //
34   ///////////////////////////////////////////////////////////////////
35
36   TH1F* hCsiMCPithya = new TH1F();
37   TH1F* hCsiMCevtgen = new TH1F();
38   // background parameters
39   bkgParam[0] = 3.07667e-04;
40   bkgParam[1] = 3.44581e-04;
41   bkgParam[2] = 9.31570e-03;
42   bkgParam[3] = 1.76852e+03;
43   bkgParam[4] = 6.57758e+03;
44   bkgParam[5] = 6.10483e+03;
45   bkgParam[6] = 1.63737e+04;
46   bkgParam[7] = 1.33739e-01;
47   bkgParam[8] = 1.02785e+00;
48   bkgParam[9] = 1.07043e+00;
49
50   // resolution parameters FF
51   resParamFF[0] = 1.30781e+03;
52   resParamFF[1] = 0.;
53   resParamFF[2] = 2.65464e+01;
54   resParamFF[3] = 2.09853e+03;
55   resParamFF[4] = 0.;
56   resParamFF[5] = 5.75582e+01;
57   resParamFF[6] = 1.77041e+02;
58   resParamFF[7] = 3.49328e+00;
59   resParamFF[8] = 3.74067e+02;
60
61   // resolution parameters FS
62   resParamFS[0] = 2.88183e+03;
63   resParamFS[1] = 1.65460e+00;
64   resParamFS[2] = 1.11174e+02;
65   resParamFS[3] = 2.10971e+03;
66   resParamFS[4] = 1.11906e+00;
67   resParamFS[5] = 4.87024e+01;
68   resParamFS[6] = 3.56243e+02;
69   resParamFS[7] = 3.22939e+00;
70   resParamFS[8] = 5.94102e+02;
71  
72
73   // resolution parameters SS
74   resParamSS[0] = 5.77864e+03;
75   resParamSS[1] = 1.38674e+00;
76   resParamSS[2] = 2.80921e+02;
77   resParamSS[3] = 1.17984e+04;
78   resParamSS[4] = -5.67447e-01;
79   resParamSS[5] = 1.04591e+02;
80   resParamSS[6] = 1.14000e+03;
81   resParamSS[7] = 1.25609e+00;
82   resParamSS[8] = 1.43915e+03;
83
84   // invariant mass parameters
85   invMassParam[0] = 3.08547e+00;
86   invMassParam[1] = 3.21099e-02;
87   invMassParam[2] = 5.17599e-01;
88   invMassParam[3] = 9.62091e+01;
89   invMassParam[4] = 1.36800e+00;
90   invMassParam[5] = 1.24918e+02;
91   invMassParam[6] = 1.12386e+00;
92   invMassParam[7] = 1.63711e+00;
93   invMassParam[8] = -1.20962e+01; 
94
95   hCsiMCPithya = psproperMCsecJPSI();  
96   Double_t integral = 0;
97   for(int i=1;i<hCsiMCPithya->GetNbinsX()+1; i++) integral += (hCsiMCPithya->GetBinContent(i)*hCsiMCPithya->GetBinWidth(i));
98   hCsiMCPithya->Scale(1./integral);
99
100   Double_t* x=0x0; Double_t* m=0x0; Int_t* type=0; Int_t n=0; 
101   AliDielectronBtoJPSItoEle* aBtoJPSItoEle =new AliDielectronBtoJPSItoEle();
102   Double_t paramInputValues[45] =  { bkgParam[3], // fWeightRes [0]
103                                      bkgParam[4], // Fpos [1]
104                                      bkgParam[5], // FNeg [2]
105                                      bkgParam[6], // FSym [3]
106                                      bkgParam[0], // LamdaPos [4]
107                                      bkgParam[1], // LambdaNeg [5]
108                                      bkgParam[2], // LamdaSym [6]
109                                      Fb, // fB [7]
110                                      Fsig, // FSig [8]
111                                      invMassParam[0], // FMean [9]
112                                      invMassParam[4], // fNexp [10]
113                                      invMassParam[1], // fNsigma [11]
114                                      invMassParam[2], // fAlpha [12]
115                                      invMassParam[3],//fNorm [13]
116                                      invMassParam[5], // fBkgNorm [14]
117                                      invMassParam[6], // fBkgMean [15]
118                                      invMassParam[7], // fBkgSlope [16]
119                                      invMassParam[8], // fBkgConst [17]
120                                      resParamFF[0], // Gaus1Norm [18]
121                                      resParamFF[3], // Gaus2Norm [19]
122                                      resParamFF[1], // Mean1Res[20]
123                                      resParamFF[2], // Sigma1Res[21]
124                                      resParamFF[4], // Mean2Res[22]
125                                      resParamFF[5],  // Sigma2Res[23]
126                                      resParamFF[6], // alfa res [24]
127                                      resParamFF[7], // lambda res [25]
128                                      resParamFF[8], // norm res [26]
129                                      resParamFS[0], // Gaus1Norm [27]
130                                      resParamFS[3], // Gaus2Norm [28]
131                                      resParamFS[1], // Mean1Res[29]
132                                      resParamFS[2], // Sigma1Res[30]
133                                      resParamFS[4], // Mean2Res[31]
134                                      resParamFS[5],  // Sigma2Res[32]
135                                      resParamFS[6], // alfa res [33]
136                                      resParamFS[7], // lambda res [34]
137                                      resParamFS[8], // norm res [35]
138                                      resParamSS[0], // Gaus1Norm [36]
139                                      resParamSS[3], // Gaus2Norm [37]
140                                      resParamSS[1], // Mean1Res[38]
141                                      resParamSS[2], // Sigma1Res[39]
142                                      resParamSS[4], // Mean2Res[40]
143                                      resParamSS[5], // Sigma2Res[41]
144                                      resParamSS[6], // alfa res [42]
145                                      resParamSS[7], // lambda res [43]
146                                      resParamSS[8]  // norm res [44]
147                                      }; // Starting values for parameters 
148  
149   TFile f(inputDistr);
150   TList * list = (TList*)f.Get("resultAOD");
151   TNtuple *nt;
152   if(range == 0) nt=(TNtuple*)list->FindObject("fNtupleJPSItype");
153   if(range == 1) nt=(TNtuple*)list->FindObject("fNtupleJPSItype_signal");
154   nt->ls();
155
156   aBtoJPSItoEle->SetResTypeAnalysis(resType);
157   aBtoJPSItoEle->ReadCandidates(nt,x,m,type,n); // read N-Tuples
158   printf("+++\n+++ Number of total candidates (prim J/psi + secondary J/psi + bkg) ---> %d candidates \n+++\n",n);
159
160   aBtoJPSItoEle->SetFitHandler(x,m,type,n); // Set the fit handler with given values of x, m, # of candidates 
161
162   aBtoJPSItoEle->CloneMCtemplate(hCsiMCPithya);    // clone MC template and copy internally
163                                                    // in this way any model can be setted from outside
164   //aBtoJPSItoEle->CloneMCtemplate(hCsiMCEvtGen);  // clone MC template and copy internally
165                                                    // in this way any model can be setted from outside
166
167   aBtoJPSItoEle->SetCsiMC(); // Pass the MC template to the CDF fit function
168
169
170   AliDielectronBtoJPSItoEleCDFfitHandler* fithandler = aBtoJPSItoEle->GetCDFFitHandler(); // Get the fit handler
171
172   //
173   // Set some fit options through the handler class
174   //
175   fithandler->SetErrorDef(0.5); // tells Minuit that the error interval is the one in which
176                                 // the function differs from the minimum for less than setted value
177
178   fithandler->SetPrintStatus(kTRUE);
179
180   fithandler->SetParamStartValues(paramInputValues);
181   fithandler->SetCrystalBallFunction(kTRUE);
182   
183   fithandler->FixParam(0,kTRUE); // bkg weights
184   fithandler->FixParam(1,kTRUE); 
185   fithandler->FixParam(2,kTRUE);
186   fithandler->FixParam(3,kTRUE);
187   fithandler->FixParam(4,kTRUE); // lPos
188   fithandler->FixParam(5,kTRUE); // lNeg  Fix bkg exponential   
189   fithandler->FixParam(6,kTRUE); // lSym
190   //fithandler->FixParam(8,kTRUE); // Fsig fix
191   
192   fithandler->FixParam(9,kTRUE);
193   fithandler->FixParam(10,kTRUE); // Fix CristalBall param
194   fithandler->FixParam(11,kTRUE);
195   fithandler->FixParam(12,kTRUE);  
196   fithandler->FixParam(13,kTRUE); 
197   fithandler->FixParam(14,kTRUE);
198   fithandler->FixParam(15,kTRUE); // Fix Bkg invMass param
199   fithandler->FixParam(16,kTRUE);
200   fithandler->FixParam(17,kTRUE);
201   fithandler->FixParam(18,kTRUE);
202   fithandler->FixParam(19,kTRUE); // norm2
203   fithandler->FixParam(20,kTRUE);
204   fithandler->FixParam(21,kTRUE);
205   fithandler->FixParam(22,kTRUE); // mean2
206   fithandler->FixParam(23,kTRUE); // sigma2
207   fithandler->FixParam(24,kTRUE); // alfaRes
208   fithandler->FixParam(25,kTRUE); // lambdaRes
209   fithandler->FixParam(26,kTRUE); // ConstRes
210   // resolution func for FS and SS types
211   fithandler->FixParam(27,kTRUE);
212   fithandler->FixParam(28,kTRUE); // norm2
213   fithandler->FixParam(29,kTRUE);
214   fithandler->FixParam(30,kTRUE);
215   fithandler->FixParam(31,kTRUE); // mean2
216   fithandler->FixParam(32,kTRUE); // sigma2
217   fithandler->FixParam(33,kTRUE); // alfaRes
218   fithandler->FixParam(34,kTRUE); // lambdaRes
219   fithandler->FixParam(35,kTRUE); // ConstRes
220   fithandler->FixParam(36,kTRUE);
221   fithandler->FixParam(37,kTRUE); // norm2
222   fithandler->FixParam(38,kTRUE);
223   fithandler->FixParam(39,kTRUE);
224   fithandler->FixParam(40,kTRUE); // mean2
225   fithandler->FixParam(41,kTRUE); // sigma2
226   fithandler->FixParam(42,kTRUE); // alfaRes
227   fithandler->FixParam(43,kTRUE); // lambdaRes
228   fithandler->FixParam(44,kTRUE); // ConstRes
229
230   //invMass funcions normalized between bandLow - bandUp -->
231   Double_t bandLow;
232   Double_t bandUp;
233
234   if(range==1){
235    bandLow = 2.92;
236    bandUp = 3.16;
237   }
238   
239   if(range==0){
240   bandLow = 2.4;
241   bandUp = 4.;
242   }
243   
244   Double_t massLow = (TDatabasePDG::Instance()->GetParticle(443)->Mass()) - bandLow;
245   Double_t massHigh = bandUp - (TDatabasePDG::Instance()->GetParticle(443)->Mass());
246   fithandler->SetMassWndLow(massLow);
247   fithandler->SetMassWndHigh(massHigh);
248
249   likely_obj = fithandler->LikelihoodPointer();
250   likely_obj->SetAllParameters(paramInputValues); 
251   likely_obj->ComputeMassIntegral();
252   likely_obj->PrintStatus(); 
253
254   aBtoJPSItoEle->DoMinimization();
255  
256   Double_t Fsig_fromFit = fithandler->GetParameter(8);
257   Double_t Fb_fromFit = fithandler->GetParameter(7);
258
259   Double_t Fsig_err = fithandler->GetParameterError(8);
260   Double_t Fb_err = fithandler->GetParameterError(7);  
261
262   //fill histos from Ntupla
263   TH1F *fInvMass = new TH1F("fInvMass","Invariant Mass; InvMass[GeV]; Entries/40MeV",40,2.4,2.4+40*.04); // step 40MeV
264   TH1F *fpsproperSignal = new TH1F("psproper_decay_length",Form("psproper_decay_length_distrib(%f < M < %f);X [#mu m];Entries/40#mu m",bandLow,bandUp),150,-3000.,3000.); 
265   
266   Float_t mass , psproper, typeCand;
267   TString arrType[]={"SS","FS","FF"}; 
268   nt->SetBranchAddress("Xdecaytime",&psproper);
269   nt->SetBranchAddress("Mass",&mass);
270   nt->SetBranchAddress("Type",&typeCand);
271   Int_t fNcurrent=0; Double_t nCandSel = 0;
272   Int_t nb = (Int_t)nt->GetEvent(fNcurrent);
273
274   for (Int_t iev=0; iev<(nt->GetEntries()); iev++){
275    if(resType.Contains(arrType[(Int_t)typeCand])){
276    nCandSel += 1;
277    weightType[(Int_t)typeCand] += 1.;
278    fInvMass->Fill(mass);
279    fpsproperSignal->Fill(psproper);}
280    fNcurrent++;
281    nb = (Int_t) nt->GetEvent(fNcurrent);
282    }
283    likely_obj->SetWeightType(weightType[2]/nCandSel,weightType[1]/nCandSel,weightType[0]/nCandSel);
284   
285   // draw psproper total
286   TCanvas *d6 = new TCanvas();
287   d6->SetLogy();
288   Double_t norm1 = ((Double_t)fpsproperSignal->GetEntries())*fpsproperSignal->GetBinWidth(1);
289   TF1 *psproperTot = likely_obj->GetEvaluateCDFDecayTimeTotalDistrAllTypes(-1.e+04, 1.e+04,norm1);
290   TFitResultPtr rPsproper = fpsproperSignal->Fit(psproperTot->GetName(),"S");
291   TLatex *lat=new TLatex;
292   lat->SetNDC(kTRUE);
293   lat->SetTextColor(1);lat->SetTextFont(42);lat->SetTextSize(.035);
294   fpsproperSignal->DrawCopy("E");
295   lat->DrawLatex(0.53, 0.82, Form("#chi^{2}/dof = %4.3f ",(rPsproper->Chi2()/(Double_t)rPsproper->Ndf())));
296   if(range == 0) lat->DrawLatex(0.53, 0.72, Form("F_{Sig}[%1.2f - %1.2f] = %4.3f #pm %4.3f", bandLow, bandUp, Fsig_fromFit, Fsig_err));
297   else lat->DrawLatex(0.53, 0.72, Form("F_{Sig}[%1.2f - %1.2f](fixed from LS) = %4.3f", bandLow, bandUp, Fsig_fromFit));
298   lat->DrawLatex(0.53, 0.62, Form("F_{B} = %4.3f #pm %4.3f", Fb_fromFit, Fb_err)); 
299
300   //prompt jpsi
301   Double_t normPrompt = (psproperTot->GetParameter(0))*Fsig_fromFit*(1 - Fb_fromFit);
302   TF1 *prompt = likely_obj->GetResolutionFuncAllTypes(-1.e+04,1.e+04,normPrompt);
303   prompt->SetLineColor(2);
304   prompt->Draw("same");
305  
306   //secondary Jpsi 
307   Double_t normSec = (psproperTot->GetParameter(0))*Fsig_fromFit*Fb_fromFit;
308   TF1 *templateMC = likely_obj->GetFunBAllTypes(-1.e+04,1.e+04,normSec);
309   templateMC->SetLineColor(6);
310   templateMC->SetFillColor(6);
311   templateMC->SetFillStyle(3005);
312   templateMC->Draw("same");
313  
314   //bkg
315   Double_t normBkg = (psproperTot->GetParameter(0))*(1 - Fsig_fromFit);
316   TF1 *psProperBack = likely_obj->GetEvaluateCDFDecayTimeBkgDistrAllTypes(-1.e+04,1.e+04,normBkg);
317   psProperBack->SetLineColor(3);
318   psProperBack->Draw("same");
319
320   fpsproperSignal->Sumw2();
321   fpsproperSignal->DrawCopy("same"); 
322   
323   //legend
324   TLegend *leg=new TLegend(0.17,0.72,0.42,0.88);
325   leg->SetBorderSize(0); leg->SetFillColor(0); leg->SetTextFont(42);
326   leg->SetFillStyle(0); leg->SetMargin(0.25); //separation symbol-text
327   leg->SetEntrySeparation(0.15);
328   leg->AddEntry(psproperTot, "all","l");
329   leg->AddEntry(prompt, "prompt J/#psi","l");
330   leg->AddEntry(templateMC, "secondary J/#psi","l");
331   leg->AddEntry(psProperBack, "bkg","l");
332   leg->Draw("same");
333
334   // draw invariant mass 
335   Double_t norm =((Double_t)fInvMass->GetEntries())*fInvMass->GetBinWidth(1);
336   TF1 *invMassFunc = likely_obj->GetEvaluateCDFInvMassTotalDistr(bandLow,bandUp,norm);
337   TCanvas *d5 = new TCanvas();
338   Double_t intTot = invMassFunc->Integral(bandLow,bandUp); 
339   printf("intTot (%f-%f)= %f \n",bandLow,bandUp,intTot);
340   TFitResultPtr rMass = fInvMass->Fit(invMassFunc->GetName(),"S");
341   fInvMass->DrawCopy("E"); 
342   lat->DrawLatex(0.53, 0.82, Form("#chi^{2}/dof = %4.2f ",(rMass->Chi2()/(Double_t)rMass->Ndf())));
343   if(range == 0) lat->DrawLatex(0.53, 0.72, Form("F_{Sig}[%1.2f - %1.2f] = %4.3f #pm %4.3f", bandLow, bandUp, Fsig_fromFit, Fsig_err));
344   else lat->DrawLatex(0.53, 0.72, Form("F_{Sig}[%1.2f - %1.2f](fixed from LS) = %4.3f", bandLow, bandUp, Fsig_fromFit));
345
346   TF1 *invMassSig = likely_obj->GetEvaluateCDFInvMassSigDistr(2.4,4,norm*Fsig_fromFit);
347   TF1 *invMassBkg = likely_obj->GetEvaluateCDFInvMassBkgDistr(2.4,4,norm*(1.-Fsig_fromFit));
348
349   invMassSig->SetLineColor(2);
350   invMassBkg->SetLineColor(3);
351
352   invMassSig->Draw("same");
353   invMassBkg->Draw("same");
354   
355   if(range == 0){ 
356   Double_t integSig = invMassSig->Integral(2.92,3.16);
357   Double_t integBkg = invMassBkg->Integral(2.92,3.16);
358   Double_t Fsig_new = integSig/(integSig+integBkg); 
359   printf(" Fsig(rescaled) = %f \n",Fsig_new);
360   }
361   
362   f.Close();
363   return;
364 }
365
366 void LoadLib(){
367   gSystem->Load("libTree.so");
368   gSystem->Load("libGeom.so");
369   gSystem->Load("libPhysics.so");
370   gSystem->Load("libVMC.so");
371   gSystem->Load("libMinuit.so");
372   gSystem->Load("libSTEERBase.so");
373   gSystem->Load("libESD.so");
374   gSystem->Load("libAOD.so");
375   gSystem->Load("libANALYSIS.so");
376   gSystem->Load("libANALYSISalice.so");
377   gSystem->Load("libCORRFW.so");
378   gSystem->Load("libPWG3dielectron.so");
379  }
380
381
382 TH1F *psproperMCsecJPSI(){
383    TH1F *hDecayTimeMCjpsifromB = new TH1F("hDecayTimeMCjpsifromB","B --> J/#Psi MC pseudo proper decay length",150,-3000,3000);
384    hDecayTimeMCjpsifromB->SetBinContent(67,1);
385    hDecayTimeMCjpsifromB->SetBinContent(68,2);
386    hDecayTimeMCjpsifromB->SetBinContent(69,1);
387    hDecayTimeMCjpsifromB->SetBinContent(70,1);
388    hDecayTimeMCjpsifromB->SetBinContent(71,1);
389    hDecayTimeMCjpsifromB->SetBinContent(72,9);
390    hDecayTimeMCjpsifromB->SetBinContent(73,16);
391    hDecayTimeMCjpsifromB->SetBinContent(74,51);
392    hDecayTimeMCjpsifromB->SetBinContent(75,204);
393    hDecayTimeMCjpsifromB->SetBinContent(76,6494);
394    hDecayTimeMCjpsifromB->SetBinContent(77,5249);
395    hDecayTimeMCjpsifromB->SetBinContent(78,4530);
396    hDecayTimeMCjpsifromB->SetBinContent(79,3786);
397    hDecayTimeMCjpsifromB->SetBinContent(80,3376);
398    hDecayTimeMCjpsifromB->SetBinContent(81,2892);
399    hDecayTimeMCjpsifromB->SetBinContent(82,2585);
400    hDecayTimeMCjpsifromB->SetBinContent(83,2261);
401    hDecayTimeMCjpsifromB->SetBinContent(84,1997);
402    hDecayTimeMCjpsifromB->SetBinContent(85,1712);
403    hDecayTimeMCjpsifromB->SetBinContent(86,1521);
404    hDecayTimeMCjpsifromB->SetBinContent(87,1318);
405    hDecayTimeMCjpsifromB->SetBinContent(88,1293);
406    hDecayTimeMCjpsifromB->SetBinContent(89,1085);
407    hDecayTimeMCjpsifromB->SetBinContent(90,936);
408    hDecayTimeMCjpsifromB->SetBinContent(91,872);
409    hDecayTimeMCjpsifromB->SetBinContent(92,768);
410    hDecayTimeMCjpsifromB->SetBinContent(93,686);
411    hDecayTimeMCjpsifromB->SetBinContent(94,594);
412    hDecayTimeMCjpsifromB->SetBinContent(95,545);
413    hDecayTimeMCjpsifromB->SetBinContent(96,545);
414    hDecayTimeMCjpsifromB->SetBinContent(97,462);
415    hDecayTimeMCjpsifromB->SetBinContent(98,405);
416    hDecayTimeMCjpsifromB->SetBinContent(99,376);
417    hDecayTimeMCjpsifromB->SetBinContent(100,367);
418    hDecayTimeMCjpsifromB->SetBinContent(101,355);
419    hDecayTimeMCjpsifromB->SetBinContent(102,252);
420    hDecayTimeMCjpsifromB->SetBinContent(103,250);
421    hDecayTimeMCjpsifromB->SetBinContent(104,207);
422    hDecayTimeMCjpsifromB->SetBinContent(105,204);
423    hDecayTimeMCjpsifromB->SetBinContent(106,182);
424    hDecayTimeMCjpsifromB->SetBinContent(107,168);
425    hDecayTimeMCjpsifromB->SetBinContent(108,125);
426    hDecayTimeMCjpsifromB->SetBinContent(109,142);
427    hDecayTimeMCjpsifromB->SetBinContent(110,116);
428    hDecayTimeMCjpsifromB->SetBinContent(111,132);
429    hDecayTimeMCjpsifromB->SetBinContent(112,100);
430    hDecayTimeMCjpsifromB->SetBinContent(113,115);
431    hDecayTimeMCjpsifromB->SetBinContent(114,93);
432    hDecayTimeMCjpsifromB->SetBinContent(115,85);
433    hDecayTimeMCjpsifromB->SetBinContent(116,96);
434    hDecayTimeMCjpsifromB->SetBinContent(117,73);
435    hDecayTimeMCjpsifromB->SetBinContent(118,76);
436    hDecayTimeMCjpsifromB->SetBinContent(119,56);
437    hDecayTimeMCjpsifromB->SetBinContent(120,53);
438    hDecayTimeMCjpsifromB->SetBinContent(121,46);
439    hDecayTimeMCjpsifromB->SetBinContent(122,50);
440    hDecayTimeMCjpsifromB->SetBinContent(123,42);
441    hDecayTimeMCjpsifromB->SetBinContent(124,40);
442    hDecayTimeMCjpsifromB->SetBinContent(125,36);
443    hDecayTimeMCjpsifromB->SetBinContent(126,43);
444    hDecayTimeMCjpsifromB->SetBinContent(127,31);
445    hDecayTimeMCjpsifromB->SetBinContent(128,20);
446    hDecayTimeMCjpsifromB->SetBinContent(129,18);
447    hDecayTimeMCjpsifromB->SetBinContent(130,23);
448    hDecayTimeMCjpsifromB->SetBinContent(131,19);
449    hDecayTimeMCjpsifromB->SetBinContent(132,18);
450    hDecayTimeMCjpsifromB->SetBinContent(133,29);
451    hDecayTimeMCjpsifromB->SetBinContent(134,20);
452    hDecayTimeMCjpsifromB->SetBinContent(135,26);
453    hDecayTimeMCjpsifromB->SetBinContent(136,13);
454    hDecayTimeMCjpsifromB->SetBinContent(137,12);
455    hDecayTimeMCjpsifromB->SetBinContent(138,12);
456    hDecayTimeMCjpsifromB->SetBinContent(139,14);
457    hDecayTimeMCjpsifromB->SetBinContent(140,8);
458    hDecayTimeMCjpsifromB->SetBinContent(141,14);
459    hDecayTimeMCjpsifromB->SetBinContent(142,9);
460    hDecayTimeMCjpsifromB->SetBinContent(143,14);
461    hDecayTimeMCjpsifromB->SetBinContent(144,13);
462    hDecayTimeMCjpsifromB->SetBinContent(145,13);
463    hDecayTimeMCjpsifromB->SetBinContent(146,11);
464    hDecayTimeMCjpsifromB->SetBinContent(147,5);
465    hDecayTimeMCjpsifromB->SetBinContent(148,11);
466    hDecayTimeMCjpsifromB->SetBinContent(149,7);
467    hDecayTimeMCjpsifromB->SetBinContent(150,6);
468    hDecayTimeMCjpsifromB->SetBinContent(151,136);
469    hDecayTimeMCjpsifromB->SetBinError(67,1);
470    hDecayTimeMCjpsifromB->SetBinError(68,1.414214);
471    hDecayTimeMCjpsifromB->SetBinError(69,1);
472    hDecayTimeMCjpsifromB->SetBinError(70,1);
473    hDecayTimeMCjpsifromB->SetBinError(71,1);
474    hDecayTimeMCjpsifromB->SetBinError(72,3);
475    hDecayTimeMCjpsifromB->SetBinError(73,4);
476    hDecayTimeMCjpsifromB->SetBinError(74,7.141428);
477    hDecayTimeMCjpsifromB->SetBinError(75,14.28286);
478    hDecayTimeMCjpsifromB->SetBinError(76,80.58536);
479    hDecayTimeMCjpsifromB->SetBinError(77,72.44998);
480    hDecayTimeMCjpsifromB->SetBinError(78,67.30527);
481    hDecayTimeMCjpsifromB->SetBinError(79,61.53048);
482    hDecayTimeMCjpsifromB->SetBinError(80,58.10336);
483    hDecayTimeMCjpsifromB->SetBinError(81,53.77732);
484    hDecayTimeMCjpsifromB->SetBinError(82,50.8429);
485    hDecayTimeMCjpsifromB->SetBinError(83,47.54997);
486    hDecayTimeMCjpsifromB->SetBinError(84,44.68781);
487    hDecayTimeMCjpsifromB->SetBinError(85,41.37632);
488    hDecayTimeMCjpsifromB->SetBinError(86,39);
489    hDecayTimeMCjpsifromB->SetBinError(87,36.30427);
490    hDecayTimeMCjpsifromB->SetBinError(88,35.95831);
491    hDecayTimeMCjpsifromB->SetBinError(89,32.93934);
492    hDecayTimeMCjpsifromB->SetBinError(90,30.59412);
493    hDecayTimeMCjpsifromB->SetBinError(91,29.52965);
494    hDecayTimeMCjpsifromB->SetBinError(92,27.71281);
495    hDecayTimeMCjpsifromB->SetBinError(93,26.1916);
496    hDecayTimeMCjpsifromB->SetBinError(94,24.37212);
497    hDecayTimeMCjpsifromB->SetBinError(95,23.34524);
498    hDecayTimeMCjpsifromB->SetBinError(96,23.34524);
499    hDecayTimeMCjpsifromB->SetBinError(97,21.49419);
500    hDecayTimeMCjpsifromB->SetBinError(98,20.12461);
501    hDecayTimeMCjpsifromB->SetBinError(99,19.39072);
502    hDecayTimeMCjpsifromB->SetBinError(100,19.15724);
503    hDecayTimeMCjpsifromB->SetBinError(101,18.84144);
504    hDecayTimeMCjpsifromB->SetBinError(102,15.87451);
505    hDecayTimeMCjpsifromB->SetBinError(103,15.81139);
506    hDecayTimeMCjpsifromB->SetBinError(104,14.38749);
507    hDecayTimeMCjpsifromB->SetBinError(105,14.28286);
508    hDecayTimeMCjpsifromB->SetBinError(106,13.49074);
509    hDecayTimeMCjpsifromB->SetBinError(107,12.96148);
510    hDecayTimeMCjpsifromB->SetBinError(108,11.18034);
511    hDecayTimeMCjpsifromB->SetBinError(109,11.91638);
512    hDecayTimeMCjpsifromB->SetBinError(110,10.77033);
513    hDecayTimeMCjpsifromB->SetBinError(111,11.48913);
514    hDecayTimeMCjpsifromB->SetBinError(112,10);
515    hDecayTimeMCjpsifromB->SetBinError(113,10.72381);
516    hDecayTimeMCjpsifromB->SetBinError(114,9.643651);
517    hDecayTimeMCjpsifromB->SetBinError(115,9.219544);
518    hDecayTimeMCjpsifromB->SetBinError(116,9.797959);
519    hDecayTimeMCjpsifromB->SetBinError(117,8.544004);
520    hDecayTimeMCjpsifromB->SetBinError(118,8.717798);
521    hDecayTimeMCjpsifromB->SetBinError(119,7.483315);
522    hDecayTimeMCjpsifromB->SetBinError(120,7.28011);
523    hDecayTimeMCjpsifromB->SetBinError(121,6.78233);
524    hDecayTimeMCjpsifromB->SetBinError(122,7.071068);
525    hDecayTimeMCjpsifromB->SetBinError(123,6.480741);
526    hDecayTimeMCjpsifromB->SetBinError(124,6.324555);
527    hDecayTimeMCjpsifromB->SetBinError(125,6);
528    hDecayTimeMCjpsifromB->SetBinError(126,6.557439);
529    hDecayTimeMCjpsifromB->SetBinError(127,5.567764);
530    hDecayTimeMCjpsifromB->SetBinError(128,4.472136);
531    hDecayTimeMCjpsifromB->SetBinError(129,4.242641);
532    hDecayTimeMCjpsifromB->SetBinError(130,4.795832);
533    hDecayTimeMCjpsifromB->SetBinError(131,4.358899);
534    hDecayTimeMCjpsifromB->SetBinError(132,4.242641);
535    hDecayTimeMCjpsifromB->SetBinError(133,5.385165);
536    hDecayTimeMCjpsifromB->SetBinError(134,4.472136);
537    hDecayTimeMCjpsifromB->SetBinError(135,5.09902);
538    hDecayTimeMCjpsifromB->SetBinError(136,3.605551);
539    hDecayTimeMCjpsifromB->SetBinError(137,3.464102);
540    hDecayTimeMCjpsifromB->SetBinError(138,3.464102);
541    hDecayTimeMCjpsifromB->SetBinError(139,3.741657);
542    hDecayTimeMCjpsifromB->SetBinError(140,2.828427);
543    hDecayTimeMCjpsifromB->SetBinError(141,3.741657);
544    hDecayTimeMCjpsifromB->SetBinError(142,3);
545    hDecayTimeMCjpsifromB->SetBinError(143,3.741657);
546    hDecayTimeMCjpsifromB->SetBinError(144,3.605551);
547    hDecayTimeMCjpsifromB->SetBinError(145,3.605551);
548    hDecayTimeMCjpsifromB->SetBinError(146,3.316625);
549    hDecayTimeMCjpsifromB->SetBinError(147,2.236068);
550    hDecayTimeMCjpsifromB->SetBinError(148,3.316625);
551    hDecayTimeMCjpsifromB->SetBinError(149,2.645751);
552    hDecayTimeMCjpsifromB->SetBinError(150,2.44949);
553    hDecayTimeMCjpsifromB->SetBinError(151,11.6619);
554    return hDecayTimeMCjpsifromB;
555 }
556
557