BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDTopipienu.cc
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2// Environment:
3// This software is part of models developed at BES collaboration
4// based on the EvtGen framework. If you use all or part
5// of it, please give an appropriate acknowledgement.
6//
7// Copyright Information: See EvtGen/BesCopyright
8// Copyright (A) 2006 Ping Rong-Gang @IHEP
9//
10// Module: EvtDTopipienu.cc
11// the necessary file: EvtDTopipienu.hh
12//
13// Description: D -> pi pi e+ nu,
14// PHYSICAL REVIEW Letter 122, 062001 (2019)
15//
16// Modification history:
17//
18// Liaoyuan Dong Jan.16, 2020 Module created
19//
20//------------------------------------------------------------------------
24#include "EvtGenBase/EvtPDL.hh"
29#include <stdlib.h>
30
32
33void EvtDTopipienu::getName(std::string& model_name){
34 model_name="DTopipienu";
35}
36
40
42 checkNArg(0);
43 checkNDaug(4);
47
48 first = 1;
49 last = 4;
50 ProbMax = 591000;
51 //cout << "Initializing EvtDTopipienu: ProbMax = " << ProbMax << endl;
52
53 type[0] = 0;
54 type[1] = 1;
55 type[2] = 2;
56 type[3] = 3;
57
58 mV = 2.01;
59 mA = 2.42;
60 V_0 = 1.6948;
61 A1_0 = 1;
62 A2_0 = 0.84489;
63
64 m0 = 0.77526;
65 width0 = 0.14910;
66 rBW = 3.0;
67 rho = 1.0;
68 phi = 0.0;
69 BF = 1.0;
70
71 m0_omega = 0.78265;
72 width0_omega = 0.00849;
73 rho_omega = 0.12902;
74 phi_omega = 2.9285;
75 BF_omega = 1.0;
76
77 m0_S = 0.953;
78 rho_S = 135.27;
79 phi_S = 3.4044;
80
81 Dp_mD = 1.86962;
82 Dp_mPi1 = 0.13957;
83 Dp_mPi2 = 0.13957;
84 D0_mD = 1.86486;
85 D0_mPi1 = 0.13957;
86 D0_mPi2 = 0.1349766;
87
88 Pi = atan2(0.0,-1.0);
89 root2 = sqrt(2.);
90 root2d3 = sqrt(2./3);
91 root1d2 = sqrt(0.5);
92 root3d2 = sqrt(1.5);
93
94 mKa = 0.493677;
95 mPi = 0.13957;
96 mEt = 0.547853;
97}
98
100 setProbMax(ProbMax);
101}
102
104/*
105 double maxprob = 0.0;
106 for(int ir=0;ir<=60000000;ir++){
107 p->initializePhaseSpace(getNDaug(),getDaugs());
108 EvtVector4R _pi1 = p->getDaug(0)->getP4();
109 EvtVector4R _pi2 = p->getDaug(1)->getP4();
110 EvtVector4R _e = p->getDaug(2)->getP4();
111 EvtVector4R _nu = p->getDaug(3)->getP4();
112
113 int pid = EvtPDL::getStdHep(p->getDaug(0)->getId());
114 int charm;
115 if(pid == -211) charm = 1;
116 else charm = -1;
117 double m2, q2, cosV, cosL, chi;
118 KinVGen(_pi1, _pi2, _e, _nu, charm, m2, q2, cosV, cosL, chi);
119 double _prob = calPDF(m2, q2, cosV, cosL, chi);
120 if(_prob>maxprob) {
121 maxprob=_prob;
122 std::cout << "Max PDF = " << ir << " charm= " << charm << " prob= " << _prob << std::endl;
123 }
124 }
125 std::cout << "Max!!!!!!!!!!! " << maxprob<< std::endl;
126*/
128 EvtVector4R pi1 = p->getDaug(0)->getP4();
129 EvtVector4R pi2 = p->getDaug(1)->getP4();
130 EvtVector4R e = p->getDaug(2)->getP4();
131 EvtVector4R nu = p->getDaug(3)->getP4();
132
133 int pid = EvtPDL::getStdHep(p->getDaug(0)->getId());
134 int charm;
135 if(pid == -211) charm = 1;
136 else charm = -1;
137 double m2, q2, cosV, cosL, chi;
138 KinVGen(pi1, pi2, e, nu, charm, m2, q2, cosV, cosL, chi);
139 double prob = calPDF(m2, q2, cosV, cosL, chi);
140 setProb(prob);
141 return;
142}
143
144void EvtDTopipienu::KinVGen(EvtVector4R vp4_K, EvtVector4R vp4_Pi, EvtVector4R vp4_Lep, EvtVector4R vp4_Nu, int charm, double& m2, double& q2, double& cosV, double& cosL, double& chi)
145{
146 EvtVector4R vp4_KPi = vp4_K + vp4_Pi;
147 EvtVector4R vp4_W = vp4_Lep + vp4_Nu;
148 m2 = vp4_KPi.mass2();
149 q2 = vp4_W.mass2();
150
151 EvtVector4R boost;
152 boost.set(vp4_KPi.get(0), -vp4_KPi.get(1), -vp4_KPi.get(2), -vp4_KPi.get(3));
153 EvtVector4R vp4_Kp = boostTo(vp4_K, boost);
154 cosV = vp4_Kp.dot(vp4_KPi)/(vp4_Kp.d3mag()*vp4_KPi.d3mag());
155
156 boost.set(vp4_W.get(0), -vp4_W.get(1), -vp4_W.get(2), -vp4_W.get(3));
157 EvtVector4R vp4_Lepp = boostTo(vp4_Lep, boost);
158 cosL = vp4_Lepp.dot(vp4_W)/(vp4_Lepp.d3mag()*vp4_W.d3mag());
159
160 EvtVector4R V = vp4_KPi/vp4_KPi.d3mag();
161 EvtVector4R C = vp4_Kp.cross(V);
162 C/=C.d3mag();
163 EvtVector4R D = vp4_Lepp.cross(V);
164 D/=D.d3mag();
165 double sinx = C.cross(V).dot(D);
166 double cosx = C.dot(D);
167 chi = sinx>0? acos(cosx): -acos(cosx);
168 if(charm==-1) chi=-chi;
169}
170
171double EvtDTopipienu::calPDF(double m2, double q2, double cosV, double cosL, double chi) {
172 double m = sqrt(m2);
173 double q = sqrt(q2);
174
175 //begin to calculate form factor
176 EvtComplex F10(0.0,0.0);
177 EvtComplex F11(0.0,0.0);
178 EvtComplex F21(0.0,0.0);
179 EvtComplex F31(0.0,0.0);
180 EvtComplex F12(0.0,0.0);
181 EvtComplex F22(0.0,0.0);
182 EvtComplex F32(0.0,0.0);
183
184 EvtComplex f10(0.0,0.0);
185 EvtComplex f11(0.0,0.0);
186 EvtComplex f21(0.0,0.0);
187 EvtComplex f31(0.0,0.0);
188
189 EvtComplex coef(0.0,0.0);
190
191 for(int index=first; index<last; index++)
192 {
193 switch(type[index])
194 {
195 case 0: //calculate form factor of P wave (rho)
196 {
197 ResonanceGS(m, q, D0_mD, D0_mPi1, D0_mPi2, f11, f21, f31);
198 coef = getCoef(rho, phi);
199 F11 = F11+ coef*f11;
200 F21 = F21+ coef*f21;
201 F31 = F31+ coef*f31;
202 break;
203 }
204 case 1: //calculate form factor of P wave (GS)
205 {
206 ResonanceGS(m, q, Dp_mD, Dp_mPi1, Dp_mPi2, f11, f21, f31);
207 coef = getCoef(rho, phi);
208 F11 = F11+ coef*f11;
209 F21 = F21+ coef*f21;
210 F31 = F31+ coef*f31;
211 break;
212 }
213 case 2: //calculate form factor of omega
214 {
215 ResonancePGScbw(m, q, f11, f21, f31);
216 coef = getCoef(rho_omega, phi_omega);
217 F11 = F11+ coef*f11;
218 F21 = F21+ coef*f21;
219 F31 = F31+ coef*f31;
220 break;
221 }
222 case 3: //calculate form factor of S wave (BW corrected by Bugg)
223 {
224 ResonanceSBugg(m, q, f10);
225 coef = getCoef(rho_S, phi_S);
226 F10 = F10 + coef*f10;
227 break;
228 }
229 default:
230 {
231 std::cout<<"No such form factor type!!!"<<std::endl;
232 break;
233 }
234 }
235 }
236
237 //begin to calculate pdf value
238 double I,I1,I2,I3,I4,I5,I6,I7,I8,I9;
239
240 double cosV2 = cosV*cosV;
241 double sinV = sqrt(1.0-cosV2);
242 double sinV2 = sinV*sinV;
243
244 EvtComplex F1 = F10 + F11*cosV + F12*(1.5*cosV2-0.5);
245 EvtComplex F2 = F21*root1d2 + F22*cosV*root3d2;
246 EvtComplex F3 = F31*root1d2 + F32*cosV*root3d2;
247
248 I1 = 0.25*( abs2(F1) + 1.5*sinV2*(abs2(F2) + abs2(F3) ) );
249 I2 = -0.25*( abs2(F1) - 0.5*sinV2*(abs2(F2) + abs2(F3) ) );
250 I3 = -0.25*( abs2(F2) - abs2(F3) )*sinV2;
251 I4 = real( conj(F1)*F2 )*sinV*0.5;
252 I5 = real( conj(F1)*F3 )*sinV;
253 I6 = real( conj(F2)*F3 )*sinV2;
254 I7 = imag( conj(F2)*F1 )*sinV;
255 I8 = imag( conj(F3)*F1 )*sinV*0.5;
256 I9 = imag( conj(F3)*F2 )*sinV2*(-0.5);
257
258 double coschi = cos(chi);
259 double sinchi = sin(chi);
260 double sin2chi = 2.0*sinchi*coschi;
261 double cos2chi = 1.0 - 2.0*sinchi*sinchi;
262
263 double sinL = sqrt(1.-cosL*cosL);
264 double sinL2 = sinL*sinL;
265 double sin2L = 2.0*sinL*cosL;
266 double cos2L = 1.0 - 2.0*sinL2;
267
268 I = I1 + I2*cos2L + I3*sinL2*cos2chi + I4*sin2L*coschi + I5*sinL*coschi + I6*cosL + I7*sinL*sinchi + I8*sin2L*sinchi + I9*sinL2*sin2chi;
269 return I;
270}
271
272void EvtDTopipienu::ResonanceGS(double m, double q, double massD, double massPi1, double massPi2, EvtComplex& F11, EvtComplex& F21, EvtComplex& F31)
273{
274 double pKPi = getPStar(massD, m, q);
275 double mD2 = massD*massD;
276 double m2 = m*m;
277 double m02 = m0*m0;
278 double q2 = q*q;
279 double mV2 = mV*mV;
280 double mA2 = mA*mA;
281 double summDm = massD+m;
282 double V = V_0 /(1.0-q2/(mV2));
283 double A1 = A1_0/(1.0-q2/(mA2));
284 double A2 = A2_0/(1.0-q2/(mA2));
285 double A = summDm*A1;
286 double B = 2.0*massD*pKPi/summDm*V;
287
288 //construct the helicity form factor
289 double H0 = 0.5/(m*q)*((mD2-m2-q2)*summDm*A1-4.0*(mD2*pKPi*pKPi)/summDm*A2);
290 double Hp = A-B;
291 double Hm = A+B;
292
293 //calculate alpha
294 //double B_Kstar = 2./3.0; //B_Kstar = Br(Kstar(892)->k pi)
295 double pStar = getPStar(m, massPi1, massPi2);
296 double pStar0 = getPStar(m0, massPi1, massPi2);
297 double alpha = sqrt(3.0*Pi*BF/(pStar0*width0));
298
299 //construct amplitudes of (non)resonance
300 double F = getF1(m, m0, massPi1, massPi2, rBW);
301 EvtComplex C(m02*(1.0+width0*getGx(m0, pStar0, massPi1, massPi2)/m0)*F, 0.0);
302 double AA = m02-m2+width0*getFx(m02, m2, pStar, pStar0, massPi1, massPi2);
303 double BB = -m0*getWidthrho(m, m0, width0, pStar, pStar0);
304 EvtComplex tmp(AA,BB);
305 EvtComplex amp = C/tmp;
306
307 double alpham2 = alpha*2.0;
308 F11 = amp*alpham2*q*H0*root2;
309 F21 = amp*alpham2*q*(Hp+Hm);
310 F31 = amp*alpham2*q*(Hp-Hm);
311}
312
313void EvtDTopipienu::ResonancePGScbw(double m, double q, EvtComplex& F11, EvtComplex& F21, EvtComplex& F31)
314{
315 double pKPi = getPStar(Dp_mD, m, q);
316 double mD2 = Dp_mD*Dp_mD;
317 double m2 = m*m;
318 double m02 = m0_omega*m0_omega;
319 double mR2 = m0*m0;
320 double q2 = q*q;
321 double mV2 = mV*mV;
322 double mA2 = mA*mA;
323 double summDm = Dp_mD+m;
324 double V = V_0 /(1.0-q2/(mV2));
325 double A1 = A1_0/(1.0-q2/(mA2));
326 double A2 = A2_0/(1.0-q2/(mA2));
327 double A = summDm*A1;
328 double B = 2.0*Dp_mD*pKPi/summDm*V;
329 //construct the helicity form factor
330 double H0 = 0.5/(m*q)*((mD2-m2-q2)*summDm*A1-4.0*(mD2*pKPi*pKPi)/summDm*A2);
331 double Hp = A-B;
332 double Hm = A+B;
333
334 //calculate alpha
335 //double B_Kstar = 2./3.0; //B_Kstar = Br(Kstar(892)->k pi)
336 double pStar = getPStar(m, Dp_mPi1, Dp_mPi2);
337 double pStar0 = getPStar(m0_omega, Dp_mPi1, Dp_mPi2);
338 double alpha = sqrt(3.0*Pi*BF_omega/(pStar0*width0_omega));
339
340 //construct amplitudes of (non)resonance
341 double F = getF1(m, m0_omega, Dp_mPi1, Dp_mPi2, rBW);
342 EvtComplex amp1(0.0, 0.0);
343 EvtComplex amp2(0.0, 0.0);
344 // CBW
345 EvtComplex C1(m0_omega*width0_omega*F,0.0);
346 double AA1 = m02-m2;
347 double BB1 = -m0_omega*width0_omega;
348 EvtComplex tmp1(AA1,BB1);
349 amp1 = C1/tmp1;
350 // GS
351 pStar0 = getPStar(m0, Dp_mPi1, Dp_mPi2);
352 EvtComplex C2(mR2*(1.0+width0*getGx(m0, pStar0, Dp_mPi1, Dp_mPi2)/m0), 0.0);
353 double AA2 = mR2-m2+width0*getFx(mR2, m2, pStar, pStar0, Dp_mPi1, Dp_mPi2);
354 double BB2 = -m0*getWidthrho(m, m0, width0, pStar, pStar0);
355 EvtComplex tmp2(AA2,BB2);
356 amp2 = C2/tmp2;
357 EvtComplex amp = amp1*amp2;
358
359 double alpham2 = alpha*2.0;
360 F11 = amp*alpham2*q*H0*root2;
361 F21 = amp*alpham2*q*(Hp+Hm);
362 F31 = amp*alpham2*q*(Hp-Hm);
363}
364
365void EvtDTopipienu::ResonanceSBugg(double m, double q, EvtComplex& F10)
366{
367 double pKPi = getPStar(Dp_mD, m, q);
368 double m2 = m*m;
369 double q2 = q*q;
370 double mA2 = mA*mA;
371
372 double sA = 0.41*mPi*mPi;
373 double mr = m0_S; //0.953;
374 double mr2 = m0_S*m0_S; //0.908209;// 0.953*0.953;
375 double alpha = 1.3;
376 double g4pi = 0.011;
377
378 EvtComplex ciR(1.0, 0.0);
379 EvtComplex ciM(0.0, 1.0);
380 EvtComplex Gamma1(0.0, 0.0);
381 EvtComplex Gamma2(0.0, 0.0);
382 EvtComplex Gamma3(0.0, 0.0);
383 EvtComplex Gamma4(0.0, 0.0);
384
385 Gamma1 = getrho(m2,mPi)*getG1(m2,mr)*(m2-sA)/(mr2-sA);
386 Gamma2 = getrho(m2,mKa)*0.6*getG1(m2,mr)*(m2/mr2)*exp(-alpha*fabs(m2-4.0*mKa*mKa));
387 Gamma3 = getrho(m2,mEt)*0.2*getG1(m2,mr)*(m2/mr2)*exp(-alpha*fabs(m2-4.0*mEt*mEt));
388 if (m > 4*mPi) Gamma4 = ciR*mr*g4pi*(1.0+exp(7.082-2.845*mr2))/(1.0+exp(7.082-2.845*m2));
389
390 double AA = mr2-m2-getG1(m2,mr)*(m2-sA)*getZ(m2, mr2)/(mr2-sA);
391 EvtComplex amp = ciR/(ciR*AA-ciM*(Gamma1+Gamma2+Gamma3+Gamma4));
392 F10 = amp*pKPi*Dp_mD/(1.0-q2/mA2);
393}
394
395double EvtDTopipienu::getPStar(double m, double m1, double m2)
396{
397 double s = m*m;
398 double s1 = m1*m1;
399 double s2 = m2*m2;
400 double x = s + s1 - s2;
401 double t = 0.25*x*x/s - s1;
402 double p;
403 if (t>0.0) {
404 p = sqrt(t);
405 } else {
406 // std::cout << " Hello, pstar is less than 0.0" << std::endl;
407 p = 0.04;
408 }
409 return p;
410}
411
412double EvtDTopipienu::getF1(double m, double m0, double m_c1, double m_c2, double rBW)
413{
414 double pStar = getPStar(m, m_c1, m_c2);
415 double pStar0 = getPStar(m0, m_c1, m_c2);
416 double rBW2 = rBW*rBW;
417 double pStar2 = pStar*pStar;
418 double pStar02 = pStar0*pStar0;
419 double B = 1./sqrt(1.+rBW2*pStar2);
420 double B0 = 1./sqrt(1.+rBW2*pStar02);
421 double F = pStar/pStar0*B/B0;
422 return F;
423}
424
425double EvtDTopipienu::getF2(double m, double m0, double m_c1, double m_c2, double rBW)
426{
427 double pStar = getPStar(m, m_c1, m_c2);
428 double pStar0 = getPStar(m0, m_c1, m_c2);
429 double rBW2 = rBW*rBW;
430 double pStar2 = pStar*pStar;
431 double pStar02 = pStar0*pStar0;
432 double B = 1./sqrt((rBW2*pStar2-3.)*(rBW2*pStar2-3.)+9.*rBW2*pStar2);
433 double B0 = 1./sqrt((rBW2*pStar02-3.)*(rBW2*pStar02-3.)+9.*rBW2*pStar02);
434 double F = pStar2/pStar02*B/B0;
435 return F;
436}
437
438double EvtDTopipienu::getWidth0(double m, double m0, double m_c1, double m_c2, double width0)
439{
440 double pStar = getPStar(m, m_c1, m_c2);
441 double pStar0 = getPStar(m0, m_c1, m_c2);
442 double width = width0*pStar/pStar0*m0/m;
443 return width;
444}
445
446double EvtDTopipienu::getWidth1(double m, double m0, double m_c1, double m_c2, double width0, double rBW)
447{
448 double pStar = getPStar(m, m_c1, m_c2);
449 double pStar0 = getPStar(m0, m_c1, m_c2);
450 double F = getF1(m, m0, m_c1, m_c2, rBW);
451 double width = width0*pStar/pStar0*m0/m*F*F;
452 return width;
453}
454
455double EvtDTopipienu::getWidth2(double m, double m0, double m_c1, double m_c2, double width0, double rBW)
456{
457 double pStar = getPStar(m, m_c1, m_c2);
458 double pStar0 = getPStar(m0, m_c1, m_c2);
459 double F = getF2(m, m0, m_c1, m_c2, rBW);
460 double width = width0*pStar/pStar0*m0/m*F*F;
461 return width;
462}
463
464EvtComplex EvtDTopipienu::getCoef(double rho, double phi)
465{
466 EvtComplex coef ( rho*cos(phi), rho*sin(phi) );
467 return coef;
468}
469
470inline double EvtDTopipienu::getGx(double m0, double p0, double m_c1, double m_c2)
471{
472 double Gg = 0;
473 double MPI = 0.5*(m_c1+m_c2);
474 Gg = 3*MPI*MPI*log((m0+2*p0)/(2*MPI))/(Pi*p0*p0) + m0/(2*Pi*p0) - MPI*MPI*m0/(Pi*p0*p0*p0);
475 return Gg;
476}
477
478inline double EvtDTopipienu::getFx(double mr2, double sx, double p, double p0, double m_c1, double m_c2)
479{
480 double Fx = 0;
481 Fx = mr2/(pow(p0,3))*(p*p*(getHx(sx,p,m_c1,m_c2)-getHx(mr2,p0,m_c1,m_c2))+(mr2-sx)*p0*p0*getdh(mr2,p0,m_c1,m_c2));
482 return Fx;
483}
484
485inline double EvtDTopipienu::getHx(double sx, double p, double m_c1, double m_c2)
486{
487 double m = sqrt(sx);
488 double Hx = 0;
489 Hx = 2*p*log((m+2*p)/(m_c1+m_c2))/(Pi*m);
490 return Hx;
491}
492
493inline double EvtDTopipienu::getdh(double mr2, double p0, double m_c1, double m_c2)
494{
495 double mass = sqrt(mr2);
496 double dh = getHx(mass,p0,m_c1,m_c2)*(1.0/(8*p0*p0)-1.0/(2*mass*mass))+1.0/(2*Pi*mass*mass);
497 return dh;
498}
499
500inline double EvtDTopipienu::getG1(double m2, double Mr)
501{
502 double b1 = 1.302;
503 double b2 = 0.340;
504 double A = 2.426;
505 double Mr2 = Mr*Mr;// 0.953*0.953;
506 double gg1 = Mr*(b1+b2*m2)*exp(-(m2-Mr2)/A);
507 return gg1;
508}
509
510inline double EvtDTopipienu::getZ(double m2, double Mr2)
511{
512 double zz = (getRho(m2,mPi)*log((1.0-getRho(m2,mPi))/(1.0+getRho(m2,mPi)))
513 -getRho(Mr2,mPi)*log((1.0-getRho(Mr2,mPi))/(1.0+getRho(Mr2,mPi))))/Pi;
514
515 return zz;
516}
517
518inline double EvtDTopipienu::getRho(double m2, double mX)
519{
520 double rho = 0.0;
521 if ((1.0-4.0*mX*mX/m2)>0)
522 rho = sqrt(1.0-4.0*mX*mX/m2);
523 else rho = 0.0;
524 return rho;
525}
526
527inline EvtComplex EvtDTopipienu::getrho(double m2, double mX)
528{
529 EvtComplex rho(0.0, 0.0);
530 EvtComplex ciR(1.0, 0.0);
531 EvtComplex ciM(0.0, 1.0);
532 if ((1.0-4.0*mX*mX/m2)>0)
533 rho = ciR*sqrt(1.0-4.0*mX*mX/m2);
534 else rho = ciM*sqrt(4.0*mX*mX/m2-1.0);
535 return rho;
536}
537inline double EvtDTopipienu::getWidthrho(double m, double m0, double width0, double p, double p0 )
538{
539 double widthRho = 0.0;
540 widthRho = width0*pow(p/p0, 3)*(m0/m);
541 return widthRho;
542}
double sin(const BesAngle a)
Definition BesAngle.h:210
double cos(const BesAngle a)
Definition BesAngle.h:213
double mass
Double_t x[10]
const DifComplex I
character *LEPTONflag integer iresonances real pi2
Evt3Rank3C conj(const Evt3Rank3C &t2)
double imag(const EvtComplex &c)
double abs2(const EvtComplex &c)
EvtComplex exp(const EvtComplex &c)
EvtDiracSpinor boostTo(const EvtDiracSpinor &sp, const EvtVector4R p4)
const double alpha
XmlRpcServer s
****INTEGER imax DOUBLE PRECISION m_pi *DOUBLE PRECISION m_amfin DOUBLE PRECISION m_Chfin DOUBLE PRECISION m_Xenph DOUBLE PRECISION m_sinw2 DOUBLE PRECISION m_GFermi DOUBLE PRECISION m_MfinMin DOUBLE PRECISION m_ta2 INTEGER m_out INTEGER m_KeyFSR INTEGER m_KeyQCD *COMMON c_Semalib $ !copy of input $ !CMS energy $ !beam mass $ !final mass $ !beam charge $ !final charge $ !smallest final mass $ !Z mass $ !Z width $ !EW mixing angle $ !Gmu Fermi $ alphaQED at q
Definition KKsem.h:33
***************************************************************************************Pseudo Class RRes *****************************************************************************************Parameters and physical constants **Maarten sept ************************************************************************DOUBLE PRECISION xsmu **************************************************************************PARTICLE DATA all others are from PDG *Only resonances with known widths into electron pairs are sept ************************************************************************C Declarations C
Definition RRes.h:29
TTree * t
Definition binning.cxx:23
void getName(std::string &name)
void decay(EvtParticle *p)
EvtDecayBase * clone()
virtual ~EvtDTopipienu()
void checkSpinDaughter(int d1, EvtSpinType::spintype sp)
void checkSpinParent(EvtSpinType::spintype sp)
void setProbMax(double prbmx)
void checkNDaug(int d1, int d2=-1)
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
void setProb(double prob)
static int getStdHep(EvtId id)
Definition EvtPDL.hh:56
EvtId getId() const
const EvtVector4R & getP4() const
EvtParticle * getDaug(int i)
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
double dot(const EvtVector4R &v2) const
double get(int i) const
EvtVector4R cross(const EvtVector4R &v2)
double d3mag() const
double mass2() const
void set(int i, double d)
double * m1
Definition qcdloop1.h:75
double double * m2
Definition qcdloop1.h:75