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