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