BOSS 7.0.1
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtVubNLO.cc
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2//
3// Environment:
4// This software is part of the EvtGen package developed jointly
5// for the BaBar and CLEO collaborations. If you use all or part
6// of it, please give an appropriate acknowledgement.
7//
8// Copyright Information:
9// Copyright (C) 1998 Caltech, UCSB
10//
11// Module: EvtVubNLO.cc
12//
13// Description: Routine to decay B->Xulnu according to Bosch, Lange, Neubert, and Paz hep-ph/0402094
14// Equation numbers refer to this paper
15//
16// Modification history:
17//
18// Riccardo Faccini Feb. 11, 2004
19//
20//------------------------------------------------------------------------
21//
23#include <stdlib.h>
26#include "EvtGenBase/EvtPDL.hh"
29#include <string>
36using std::cout;
37using std::endl;
38
40 delete [] _masses;
41 delete [] _weights;
42 cout <<" max pdf : "<<_gmax<<endl;
43 cout <<" efficiency : "<<(float)_ngood/(float)_ntot<<endl;
44
45}
46
47extern "C" {
48 double ddilog_(const double *);
49}
50
51void EvtVubNLO::getName(std::string& model_name){
52
53 model_name="VUB_NLO";
54
55}
56
58
59 return new EvtVubNLO;
60
61}
62
63
65
66 // max pdf
67 _gmax=0;
68 _ntot=0;
69 _ngood=0;
70 _lbar=-1000;
71 _mupi2=-1000;
72
73 // check that there are at least 6 arguments
74 int npar = 8;
75 if (getNArg()<npar) {
76
77 report(ERROR,"EvtGen") << "EvtVubNLO generator expected "
78 << " at least npar arguments but found: "
79 <<getNArg()<<endl;
80 report(ERROR,"EvtGen") << "Will terminate execution!"<<endl;
81 ::abort();
82
83 }
84 // this is the shape function parameter
85 _mb = getArg(0);
86 _b = getArg(1);
87 _lambdaSF = getArg(2);// shape function lambda is different from lambda
88 _mui = 1.5;// GeV (scale)
89 _kpar = getArg(3);// 0
90 _idSF = abs((int)getArg(4));// type of shape function 1: exponential (from Neubert)
91 _nbins = abs((int)getArg(5));
92 _masses = new double[_nbins];
93 _weights = new double[_nbins];
94
95 // Shape function normalization
96 _mB=5.28;// temporary B meson mass for normalization
97
98 std::vector<double> sCoeffs(11);
99 sCoeffs[3] = _b;
100 sCoeffs[4] = _mb;
101 sCoeffs[5] = _mB;
102 sCoeffs[6] = _idSF;
103 sCoeffs[7] = lambda_SF();
104 sCoeffs[8] = mu_h();
105 sCoeffs[9] = mu_i();
106 sCoeffs[10] = 1.;
107 _SFNorm = SFNorm(sCoeffs) ; // SF normalization;
108
109
110 cout << " pdf 0.66, 1.32 , 4.32 "<<tripleDiff(0.66, 1.32 , 4.32)<<endl;
111 cout << " pdf 0.23,0.37,3.76 "<<tripleDiff(0.23,0.37,3.76)<<endl;
112 cout << " pdf 0.97,4.32,4.42 "<<tripleDiff(0.97,4.32,4.42)<<endl;
113 cout << " pdf 0.52,1.02,2.01 "<<tripleDiff(0.52,1.02,2.01)<<endl;
114 cout << " pdf 1.35,1.39,2.73 "<<tripleDiff(1.35,1.39,2.73)<<endl;
115
116
117 if (getNArg()-npar+2 != 2*_nbins) {
118 report(ERROR,"EvtGen") << "EvtVubNLO generator expected "
119 << _nbins << " masses and weights but found: "
120 <<(getNArg()-npar)/2 <<endl;
121 report(ERROR,"EvtGen") << "Will terminate execution!"<<endl;
122 ::abort();
123 }
124 int i,j = npar-2;
125 double maxw = 0.;
126 for (i=0;i<_nbins;i++) {
127 _masses[i] = getArg(j++);
128 if (i>0 && _masses[i] <= _masses[i-1]) {
129 report(ERROR,"EvtGen") << "EvtVubNLO generator expected "
130 << " mass bins in ascending order!"
131 << "Will terminate execution!"<<endl;
132 ::abort();
133 }
134 _weights[i] = getArg(j++);
135 if (_weights[i] < 0) {
136 report(ERROR,"EvtGen") << "EvtVubNLO generator expected "
137 << " weights >= 0, but found: "
138 <<_weights[i] <<endl;
139 report(ERROR,"EvtGen") << "Will terminate execution!"<<endl;
140 ::abort();
141 }
142 if ( _weights[i] > maxw ) maxw = _weights[i];
143 }
144 if (maxw == 0) {
145 report(ERROR,"EvtGen") << "EvtVubNLO generator expected at least one "
146 << " weight > 0, but found none! "
147 << "Will terminate execution!"<<endl;
148 ::abort();
149 }
150 for (i=0;i<_nbins;i++) _weights[i]/=maxw;
151
152 // the maximum dGamma*p2 value depends on alpha_s only:
153
154
155 // _dGMax = 0.05;
156 _dGMax = 150.;
157
158 // for the Fermi Motion we need a B-Meso\n mass - but it's not critical
159 // to get an exact value; in order to stay in the phase space for
160 // B+- and B0 use the smaller mass
161
162
163 // check that there are 3 daughters
164 checkNDaug(3);
165}
166
168
169 noProbMax();
170
171}
172
174
175 int j;
176 // B+ -> u-bar specflav l+ nu
177
178 EvtParticle *xuhad, *lepton, *neutrino;
179 EvtVector4R p4;
180
181 double pp,pm,pl,ml,El(0.0),Eh(0.0),sh(0.0);
182
183
184
186
187 xuhad=p->getDaug(0);
188 lepton=p->getDaug(1);
189 neutrino=p->getDaug(2);
190
191 _mB = p->mass();
192 ml = lepton->mass();
193
194 bool tryit = true;
195
196 while (tryit) {
197 // pm=(E_H+P_H)
198 pm= EvtRandom::Flat(0.,1);
199 pm= pow(pm,1./3.)*_mB;
200 // pl=mB-2*El
201 pl = EvtRandom::Flat(0.,1);
202 pl=sqrt(pl)*pm;
203 // pp=(E_H-P_H)
204 pp = EvtRandom::Flat(0.,pl);
205
206 _ntot++;
207
208 El = (_mB-pl)/2.;
209 Eh = (pp+pm)/2;
210 sh = pp*pm;
211
212 double pdf(0.);
213 if (pp<pl && El>ml&& sh > _masses[0]*_masses[0]&& _mB*_mB + sh - 2*_mB*Eh > ml*ml) {
214 double xran = EvtRandom::Flat(0,_dGMax);
215 pdf = tripleDiff(pp,pl,pm); // triple differential distribution
216 // cout <<" P+,P-,Pl,Pdf= "<<pp <<" "<<pm<<" "<<pl<<" "<<pdf<<endl;
217 if(pdf>_dGMax){
218 report(ERROR,"EvtGen") << "EvtVubNLO pdf above maximum: " <<pdf
219 <<" P+,P-,Pl,Pdf= "<<pp <<" "<<pm<<" "<<pl<<" "<<pdf<<endl;
220 //::abort();
221
222 }
223 if ( pdf >= xran ) tryit = false;
224
225 if(pdf>_gmax)_gmax=pdf;
226 } else {
227 // cout <<" EvtVubNLO incorrect kinematics sh= "<<sh<<"EH "<<Eh<<endl;
228 }
229
230
231 // reweight the Mx distribution
232 if(!tryit && _nbins>0){
233 _ngood++;
234 double xran1 = EvtRandom::Flat();
235 double m = sqrt(sh);j=0;
236 while ( j < _nbins && m > _masses[j] ) j++;
237 double w = _weights[j-1];
238 if ( w < xran1 ) tryit = true;// through away this candidate
239 }
240 }
241
242 // cout <<" max prob "<<gmax<<" " << pp<<" "<<y<<" "<<x<<endl;
243
244 // o.k. we have the three kineamtic variables
245 // now calculate a flat cos Theta_H [-1,1] distribution of the
246 // hadron flight direction w.r.t the B flight direction
247 // because the B is a scalar and should decay isotropic.
248 // Then chose a flat Phi_H [0,2Pi] w.r.t the B flight direction
249 // and and a flat Phi_L [0,2Pi] in the W restframe w.r.t the
250 // W flight direction.
251
252 double ctH = EvtRandom::Flat(-1,1);
253 double phH = EvtRandom::Flat(0,2*M_PI);
254 double phL = EvtRandom::Flat(0,2*M_PI);
255
256 // now compute the four vectors in the B Meson restframe
257
258 double ptmp,sttmp;
259 // calculate the hadron 4 vector in the B Meson restframe
260
261 sttmp = sqrt(1-ctH*ctH);
262 ptmp = sqrt(Eh*Eh-sh);
263 double pHB[4] = {Eh,ptmp*sttmp*cos(phH),ptmp*sttmp*sin(phH),ptmp*ctH};
264 p4.set(pHB[0],pHB[1],pHB[2],pHB[3]);
265 xuhad->init( getDaug(0), p4);
266
267
268 // calculate the W 4 vector in the B Meson restrframe
269
270 double apWB = ptmp;
271 double pWB[4] = {_mB-Eh,-pHB[1],-pHB[2],-pHB[3]};
272
273 // first go in the W restframe and calculate the lepton and
274 // the neutrino in the W frame
275
276 double mW2 = _mB*_mB + sh - 2*_mB*Eh;
277 // if(mW2<0.1){
278 // cout <<" low Q2! "<<pp<<" "<<epp<<" "<<x<<" "<<y<<endl;
279 //}
280 double beta = ptmp/pWB[0];
281 double gamma = pWB[0]/sqrt(mW2);
282
283 double pLW[4];
284
285 ptmp = (mW2-ml*ml)/2/sqrt(mW2);
286 pLW[0] = sqrt(ml*ml + ptmp*ptmp);
287
288 double ctL = (El - gamma*pLW[0])/beta/gamma/ptmp;
289 if ( ctL < -1 ) ctL = -1;
290 if ( ctL > 1 ) ctL = 1;
291 sttmp = sqrt(1-ctL*ctL);
292
293 // eX' = eZ x eW
294 double xW[3] = {-pWB[2],pWB[1],0};
295 // eZ' = eW
296 double zW[3] = {pWB[1]/apWB,pWB[2]/apWB,pWB[3]/apWB};
297
298 double lx = sqrt(xW[0]*xW[0]+xW[1]*xW[1]);
299 for (j=0;j<2;j++)
300 xW[j] /= lx;
301
302 // eY' = eZ' x eX'
303 double yW[3] = {-pWB[1]*pWB[3],-pWB[2]*pWB[3],pWB[1]*pWB[1]+pWB[2]*pWB[2]};
304 double ly = sqrt(yW[0]*yW[0]+yW[1]*yW[1]+yW[2]*yW[2]);
305 for (j=0;j<3;j++)
306 yW[j] /= ly;
307
308 // p_lep = |p_lep| * ( sin(Theta) * cos(Phi) * eX'
309 // + sin(Theta) * sin(Phi) * eY'
310 // + cos(Theta) * eZ')
311 for (j=0;j<3;j++)
312 pLW[j+1] = sttmp*cos(phL)*ptmp*xW[j]
313 + sttmp*sin(phL)*ptmp*yW[j]
314 + ctL *ptmp*zW[j];
315
316 double apLW = ptmp;
317
318 // boost them back in the B Meson restframe
319
320 double appLB = beta*gamma*pLW[0] + gamma*ctL*apLW;
321
322 ptmp = sqrt(El*El-ml*ml);
323 double ctLL = appLB/ptmp;
324
325 if ( ctLL > 1 ) ctLL = 1;
326 if ( ctLL < -1 ) ctLL = -1;
327
328 double pLB[4] = {El,0,0,0};
329 double pNB[8] = {pWB[0]-El,0,0,0};
330
331 for (j=1;j<4;j++) {
332 pLB[j] = pLW[j] + (ctLL*ptmp - ctL*apLW)/apWB*pWB[j];
333 pNB[j] = pWB[j] - pLB[j];
334 }
335
336 p4.set(pLB[0],pLB[1],pLB[2],pLB[3]);
337 lepton->init( getDaug(1), p4);
338
339 p4.set(pNB[0],pNB[1],pNB[2],pNB[3]);
340 neutrino->init( getDaug(2), p4);
341
342 return ;
343}
344
345double
346EvtVubNLO::tripleDiff ( double pp, double pl, double pm){
347
348 std::vector<double> sCoeffs(11);
349 sCoeffs[0] = pp;
350 sCoeffs[1] = pl;
351 sCoeffs[2] = pm;
352 sCoeffs[3] = _b;
353 sCoeffs[4] = _mb;
354 sCoeffs[5] = _mB;
355 sCoeffs[6] = _idSF;
356 sCoeffs[7] = lambda_SF();
357 sCoeffs[8] = mu_h();
358 sCoeffs[9] = mu_i();
359 sCoeffs[10] = _SFNorm; // SF normalization;
360
361
362 double c1=(_mB+pl-pp-pm)*(pm-pl);
363 double c2=2*(pl-pp)*(pm-pl);
364 double c3=(_mB-pm)*(pm-pp);
365 double aF1=F10(sCoeffs);
366 double aF2=F20(sCoeffs);
367 double aF3=F30(sCoeffs);
368 double td0=c1*aF1+c2*aF2+c3*aF3;
369
370
371 EvtItgPtrFunction *func = new EvtItgPtrFunction(&integrand, 0., _mB, sCoeffs);
372 EvtItgAbsIntegrator *jetSF = new EvtItgSimpsonIntegrator(*func,0.01,25);
373 double smallfrac=0.000001;// stop a bit before the end to avoid problems with numerical integration
374 double tdInt = jetSF->evaluate(0,pp*(1-smallfrac));
375 delete jetSF;
376
377 double SU=U1lo(mu_h(),mu_i())*pow((pm-pp)/(_mB-pp),alo(mu_h(),mu_i()));
378 double TD=(_mB-pp)*SU*(td0+tdInt);
379
380 return TD;
381
382}
383
384double
385EvtVubNLO::integrand(double omega, const std::vector<double> &coeffs){
386 //double pp=coeffs[0];
387 double c1=(coeffs[5]+coeffs[1]-coeffs[0]-coeffs[2])*(coeffs[2]-coeffs[1]);
388 double c2=2*(coeffs[1]-coeffs[0])*(coeffs[2]-coeffs[1]);
389 double c3=(coeffs[5]-coeffs[2])*(coeffs[2]-coeffs[0]);
390
391 return c1*F1Int(omega,coeffs)+c2*F2Int(omega,coeffs)+c3*F3Int(omega,coeffs);
392}
393
394double
395EvtVubNLO::F10(const std::vector<double> &coeffs){
396 double pp=coeffs[0];
397 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
398 double mui=coeffs[9];
399 double muh=coeffs[8];
400 double z=1-y;
401 double result= U1nlo(muh,mui)/ U1lo(muh,mui);
402
403 result += anlo(muh,mui)*log(y);
404
405 result += C_F(muh)*(-4*pow(log(y*coeffs[4]/muh),2)+10*log(y*coeffs[4]/muh)-4*log(y)-2*log(y)/(1-y)-4.0*ddilog_(&z)-pow(EvtConst::pi,2)/6.-12 );
406
407 result += C_F(mui)*(2*pow(log(y*coeffs[4]*pp/pow(mui,2)),2)-3*log(y*coeffs[4]*pp/pow(mui,2))+7-pow(EvtConst::pi,2) );
408 result *=shapeFunction(pp,coeffs);
409 // changes due to SSF
410 result += (-subS(coeffs)+2*subT(coeffs)+(subU(coeffs)-subV(coeffs))*(1/y-1.))/(coeffs[5]-pp);
411 result += shapeFunction(pp,coeffs)/pow((coeffs[5]-coeffs[0]),2)*(-5*(lambda1()+3*lambda2())/6+2*(2*lambda1()/3-lambda2())/pow(y,2));
412 // result += (subS(coeffs)+subT(coeffs)+(subU(coeffs)-subV(coeffs))/y)/(coeffs[5]-pp);
413 // this part has been added after Feb '05
414
415 //result += shapeFunction(pp,coeffs)/pow((coeffs[5]-coeffs[0]),2)*((lambda1()+3*lambda2())/6+2*(2*lambda1()/3-lambda2())/pow(y,2));
416 return result;
417}
418
419double
420EvtVubNLO::F1Int(double omega,const std::vector<double> &coeffs){
421 double pp=coeffs[0];
422 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
423 // mubar == mui
424 return C_F(coeffs[9])*(
425 (shapeFunction(omega,coeffs)-shapeFunction(pp,coeffs))*(4*log(y*coeffs[4]*(pp-omega)/pow(coeffs[9],2))-3)/(pp-omega)+
426 (g1(y,(pp-omega)/(coeffs[5]-coeffs[0]))/(coeffs[5]-pp)*shapeFunction(omega,coeffs))
427 );
428}
429
430double
431EvtVubNLO::F20(const std::vector<double> &coeffs){
432 double pp=coeffs[0];
433 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
434 double result= C_F(coeffs[8])*log(y)/(1-y)*shapeFunction(pp,coeffs)-
435 1/y*(subS(coeffs)+2*subT(coeffs)-(subT(coeffs)+subV(coeffs))/y)/(coeffs[5]-pp);
436 // added after Feb '05
437 result += shapeFunction(pp,coeffs)/pow((coeffs[5]-coeffs[0])*y,2)*(2*lambda1()/3+4*lambda2()-y*(7/6*lambda1()+3*lambda2()));
438 return result;
439}
440
441double
442EvtVubNLO::F2Int(double omega,const std::vector<double> &coeffs){
443 double pp=coeffs[0];
444 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
445 return C_F(coeffs[9])*g3(y,(pp-omega)/(coeffs[5]-coeffs[0]))*shapeFunction(omega,coeffs)/(coeffs[5]-pp);
446}
447
448double
449EvtVubNLO::F30(const std::vector<double> &coeffs){
450 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
451 return shapeFunction(coeffs[0],coeffs)/pow((coeffs[5]-coeffs[0])*y,2)*(-2*lambda1()/3+lambda2());
452}
453
454double
455EvtVubNLO::F3Int(double omega,const std::vector<double> &coeffs){
456 double pp=coeffs[0];
457 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
458 return C_F(coeffs[9])*g3(y,(pp-omega)/(coeffs[5]-coeffs[0]))/2*shapeFunction(omega,coeffs)/(coeffs[2]-coeffs[0]);
459}
460
461double
462EvtVubNLO::g1(double y, double x){
463 double result=(y*(-9+10*y)+x*x*(-12+13*y)+2*x*(-8+6*y+3*y*y))/y/pow(1+x,2)/(x+y);
464 result -= 4*log((1+1/x)*y)/x;
465 result -=2*log(1+y/x)*(3*pow(x,4)*(-2+y)-2*pow(y,3)-4*pow(x,3)*(2+y)-2*x*y*y*(4+y)-x*x*y*(12+4*y+y*y))/x/pow((1+x)*y,2)/(x+y);
466 return result;
467}
468
469double
470EvtVubNLO::g2(double y, double x){
471 double result=y*(10*pow(x,4)+y*y+3*x*x*y*(10+y)+pow(x,3)*(12+19*y)+x*y*(8+4*y+y*y));
472 result -= 2*x*log(1+y/x)*(5*pow(x,4)+2*y*y+6*pow(x,3)*(1+2*y)+4*y*x*(1+2*y)+x*x*y*(18+5*y));
473 result *= 2/(pow(y*(1+x),2)*y*(x+y));
474 return result;
475}
476
477double
478EvtVubNLO::g3(double y, double x){
479 double result=(2*pow(y,3)*(-11+2*y)-10*pow(x,4)*(6-6*y+y*y)+x*y*y*(-94+29*y+2*y*y)+2*x*x*y*(-72+18*y+13*y*y)-pow(x,3)*(72+42*y-70*y*y+3*pow(y,3)))/(pow(y*(1+x),2)*y*(x+y));
480 result += 2*log(1+y/x)*(-6*x*pow(y,3)*(-5+y)+4*pow(y,4)+5*pow(x,5)*(6-6*y+y*y)-4*pow(x*y,2)*(-20+6*y+y*y)+pow(x,3)*y*(90-10*y-28*y*y+pow(y,3))+pow(x,4)*(36+36*y-50*y*y+4*pow(y,3)))/(pow((1+x)*y*y,2)*(x+y));
481 return result;
482}
483
484/* old version (before Feb 05 notebook from NNeubert
485
486double
487EvtVubNLO::F1Int(double omega,const std::vector<double> &coeffs){
488 double pp=coeffs[0];
489 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
490 // mubar == mui
491 return C_F(coeffs[9])*(
492 (shapeFunction(omega,coeffs)-shapeFunction(pp,coeffs))*(4*log(y*coeffs[4]*(pp-omega)/pow(coeffs[9],2))-3)/(pp-omega)-
493 (1./y/(coeffs[5]-pp)*shapeFunction(omega,coeffs)*(5-6*y+4*(3-y)*log((pp-omega)/y/coeffs[4])))
494 );
495}
496
497
498double
499EvtVubNLO::F2Int(double omega,const std::vector<double> &coeffs){
500 double pp=coeffs[0];
501 double y=(coeffs[2]-coeffs[0])/(coeffs[5]-coeffs[0]);
502 return C_F(coeffs[9])*shapeFunction(omega,coeffs)*(2-11/y-4/y*log((pp-omega)/y/coeffs[4]))/(coeffs[5]-pp);
503}
504
505double
506EvtVubNLO::F3(const std::vector<double> &coeffs){
507 return C_F(coeffs[9])*shapeFunction(omega,coeffs)/(coeffs[2]-coeffs[0]);
508}
509*/
510
511double EvtVubNLO::SFNorm( const std::vector<double> &coeffs){
512
513 double omega0=1.68;//normalization scale (mB-2*1.8)
514 if(_idSF==1){ // exponential SF
515 double omega0=1.68;//normalization scale (mB-2*1.8)
516 return M0(mu_i(),omega0)*pow(_b,_b)/lambda_SF()/ (Gamma(_b)-Gamma(_b,_b*omega0/lambda_SF()));
517 } else if(_idSF==2){ // Gaussian SF
518 double c=cGaus(_b);
519 return M0(mu_i(),omega0)*2/lambda_SF()/pow(c,-(1+_b)/2.)/
520 (Gamma((1+_b)/2)-Gamma((1+_b)/2,pow(omega0/lambda_SF(),2)*c));
521 } else {
522 report(ERROR,"EvtGen") << "unknown SF "<<_idSF<<endl;
523 return -1;
524 }
525}
526
527double
528EvtVubNLO::shapeFunction ( double omega, const std::vector<double> &sCoeffs){
529 if( sCoeffs[6]==1){
530 return sCoeffs[10]*expShapeFunction(omega, sCoeffs);
531 } else if( sCoeffs[6]==2) {
532 return sCoeffs[10]*gausShapeFunction(omega, sCoeffs);
533 } else {
534 report(ERROR,"EvtGen") << "EvtVubNLO : unknown shape function # "
535 <<sCoeffs[6]<<endl;
536 }
537 return -1.;
538}
539
540
541// SSF
542double
543EvtVubNLO::subS(const std::vector<double> &c){ return (lambda_bar(1.68)-c[0])*shapeFunction(c[0],c);}
544double
545EvtVubNLO::subT(const std::vector<double> &c){ return -3*lambda2()*subS(c)/mu_pi2(1.68);}
546double
547EvtVubNLO::subU(const std::vector<double> &c){ return -2*subS(c);}
548double
549EvtVubNLO::subV(const std::vector<double> &c){ return -subT(c);}
550
551
552double
553EvtVubNLO::lambda_bar(double omega0){
554 if(_lbar<0){
555 if(_idSF==1){ // exponential SF
556 double rat=omega0*_b/lambda_SF();
557 _lbar=lambda_SF()/_b*(Gamma(1+_b)-Gamma(1+_b,rat))/(Gamma(_b)-Gamma(_b,rat));
558 } else if(_idSF==2){ // Gaussian SF
559 double c=cGaus(_b);
560 _lbar=lambda_SF()*(Gamma(1+_b/2)-Gamma(1+_b/2,pow(omega0/lambda_SF(),2)*c))/(Gamma((1+_b)/2)-Gamma((1+_b)/2,pow(omega0/lambda_SF(),2)*c))/sqrt(c);
561 }
562 }
563 return _lbar;
564}
565
566
567double
568EvtVubNLO::mu_pi2(double omega0){
569 if(_mupi2<0){
570 if(_idSF==1){ // exponential SF
571 double rat=omega0*_b/lambda_SF();
572 _mupi2= 3*(pow(lambda_SF()/_b,2)*(Gamma(2+_b)-Gamma(2+_b,rat))/(Gamma(_b)-Gamma(_b,rat))-pow(lambda_bar(omega0),2));
573 } else if(_idSF==2){ // Gaussian SF
574 double c=cGaus(_b);
575 double m1=Gamma((3+_b)/2)-Gamma((3+_b)/2,pow(omega0/lambda_SF(),2)*c);
576 double m2=Gamma(1+_b/2)-Gamma(1+_b/2,pow(omega0/lambda_SF(),2)*c);
577 double m3=Gamma((1+_b)/2)-Gamma((1+_b)/2,pow(omega0/lambda_SF(),2)*c);
578 _mupi2= 3*pow(lambda_SF(),2)*(m1/m3-pow(m2/m3,2))/c;
579 }
580 }
581 return _mupi2;
582}
583
584double
585EvtVubNLO::M0(double mui,double omega0){
586 double mf=omega0-lambda_bar(omega0);
587 return 1+4*C_F(mui)*(-pow(log(mf/mui),2)-log(mf/mui)-pow(EvtConst::pi/2,2)/6.+mu_pi2(omega0)/3/pow(mf,2)*(log(mf/mui)-0.5));
588}
589
590double
591EvtVubNLO::alphas(double mu){
592 double Lambda4=0.302932;
593 double lg=2*log(mu/Lambda4);
594 return 4*EvtConst::pi/lg/beta0()*(1-beta1()*log(lg)/pow(beta0(),2)/lg+pow(beta1()/lg,2)/pow(beta0(),4)*(pow(log(lg)-0.5,2)-1.25+beta2()*beta0()/pow(beta1(),2)));
595}
596
597double
598EvtVubNLO::gausShapeFunction ( double omega, const std::vector<double> &sCoeffs){
599 double b=sCoeffs[3];
600 double l=sCoeffs[7];
601 double wL=omega/l;
602
603 return pow(wL,b)*exp(-cGaus(b)*wL*wL);
604}
605
606double
607EvtVubNLO::expShapeFunction ( double omega, const std::vector<double> &sCoeffs){
608 double b=sCoeffs[3];
609 double l=sCoeffs[7];
610 double wL=omega/l;
611
612 return pow(wL,b-1)*exp(-b*wL);
613}
614
615double
616EvtVubNLO::Gamma(double z) {
617
618 std::vector<double> gammaCoeffs(6);
619 gammaCoeffs[0]=76.18009172947146;
620 gammaCoeffs[1]=-86.50532032941677;
621 gammaCoeffs[2]=24.01409824083091;
622 gammaCoeffs[3]=-1.231739572450155;
623 gammaCoeffs[4]=0.1208650973866179e-2;
624 gammaCoeffs[5]=-0.5395239384953e-5;
625
626 //Lifted from Numerical Recipies in C
627 double x, y, tmp, ser;
628
629 int j;
630 y = z;
631 x = z;
632
633 tmp = x + 5.5;
634 tmp = tmp - (x+0.5)*log(tmp);
635 ser=1.000000000190015;
636
637 for (j=0;j<6;j++) {
638 y = y +1.0;
639 ser = ser + gammaCoeffs[j]/y;
640 }
641
642 return exp(-tmp+log(2.5066282746310005*ser/x));
643
644}
645
646
647
648double
649EvtVubNLO::Gamma(double z, double tmin) {
650 std::vector<double> c(1);
651 c[0]=z;
652 EvtItgPtrFunction *func = new EvtItgPtrFunction(&dgamma, tmin, 100., c);
653 EvtItgAbsIntegrator *jetSF = new EvtItgSimpsonIntegrator(*func,0.001);
654 return jetSF->evaluate(tmin,100.);
655}
TF1 * g1
Double_t x[10]
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252
ostream & report(Severity severity, const char *facility)
Definition: EvtReport.cc:36
@ ERROR
Definition: EvtReport.hh:49
double ddilog_(const double *)
double ddilog_(const double &sh)
double sin(const BesAngle a)
double cos(const BesAngle a)
#define M_PI
Definition: TConstant.h:4
static const double pi
Definition: EvtConst.hh:28
double getArg(int j)
void noProbMax()
void checkNDaug(int d1, int d2=-1)
EvtId * getDaugs()
Definition: EvtDecayBase.hh:65
EvtId getDaug(int i)
Definition: EvtDecayBase.hh:66
double evaluate(double lower, double upper) const
virtual void init(EvtId part_n, const EvtVector4R &p4)=0
EvtParticle * getDaug(int i)
Definition: EvtParticle.cc:85
double mass() const
Definition: EvtParticle.cc:127
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
static double Flat()
Definition: EvtRandom.cc:73
void set(int i, double d)
Definition: EvtVector4R.hh:183
void getName(std::string &name)
Definition: EvtVubNLO.cc:51
void initProbMax()
Definition: EvtVubNLO.cc:167
virtual ~EvtVubNLO()
Definition: EvtVubNLO.cc:39
EvtDecayBase * clone()
Definition: EvtVubNLO.cc:57
void init()
Definition: EvtVubNLO.cc:64
void decay(EvtParticle *p)
Definition: EvtVubNLO.cc:173