BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDsToKKenu.cc
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2// Environment:
3// This software is part of the EvtGen package developed jointly
4// for the BaBar and CLEO collaborations. If you use all or part
5// of it, please give an appropriate acknowledgement.
6//
7// Copyright Information: See EvtGen/COPYRIGHT
8// Copyright (C) 1998 Caltech, UCSB
9//
10// Module: EvtDsToKKenu.cc
11//
12// Description: Routine to handle three-body decays of Ds+/Ds-
13//
14// Modification history:
15// Liaoyuan Dong Feb 12 2024 Module created
16//------------------------------------------------------------------------
20#include "EvtGenBase/EvtPDL.hh"
25#include <stdlib.h>
26
28
29void EvtDsToKKenu::getName(std::string& model_name){
30 model_name="DsToKKenu";
31}
32
36
38 checkNArg(0);
39 checkNDaug(4);
43
44 //std::cout << "EvtDsToKKenu ==> Initialization !" << std::endl;
45 nAmps = 2;
46
47 //rS = -11.57; // S-wave
48 //rS1 = 0.08;
49 //a_delta = 1.94;
50 //b_delta = -0.81;
51 //m0_1430_S = 1.425;
52 //width0_1430_S = 0.270;
53 //type[0] = 0;
54
55 mV = 2.1;
56 mA = 2.5;
57
58//for Ds to f0 e nu
59 mf0 = 0.965;
60 m2f0 = mf0*mf0;
61 g1 = 0.165*mf0;
62 g2 = 4.21*0.165*mf0;
63 rho_f0 = 7.1762;
64 phi_f0 = -3.3642;
65 //mDs = 1.96834;
66 type[0] = 0;
67
68 V_0 = 1.6456;
69 A1_0 = 1;
70 A2_0 = 0.83083;
71
72 m0 = 1.0195; // P-wave K*
73 width0 = 0.004249;
74 rBW = 3.0;
75 rho = 1;
76 phi = 0;
77 type[1] = 1;
78
79 m0_1410 = 1.414; // P-wave K*(1410)
80 width0_1410 = 0.232;
81 rho_1410 = 0.1;
82 phi_1410 = 0.;
83 type[2] = 2;
84
85 TV_0 = 1; // D-wave K*2(1430)
86 T1_0 = 1;
87 T2_0 = 1;
88 m0_1430 = 1.4324;
89 width0_1430 = 0.109;
90 rho_1430 = 15;
91 phi_1430 = 0;
92 type[3] = 3;
93
94 mPi = 0.13957;
95 mPi0 = 0.1349766;
96 mK0 = 0.497611;
97
98 mD = 1.96834;//Ds
99 mK1 = 0.49368;//k
100 mK2 = 0.49368;//k
101 Pi = atan2(0.0,-1.0);
102 root2 = sqrt(2.);
103 root2d3 = sqrt(2./3);
104 root1d2 = sqrt(0.5);
105 root3d2 = sqrt(1.5);
106}
107
109 setProbMax(450880.0);
110}
111
113/*
114 double maxprob = 0.0;
115 for(int ir=0;ir<=30000000;ir++){
116 p->initializePhaseSpace(getNDaug(),getDaugs());
117 EvtVector4R _K = p->getDaug(0)->getP4();
118 EvtVector4R _pi = p->getDaug(1)->getP4();
119 EvtVector4R _e = p->getDaug(2)->getP4();
120 EvtVector4R _nu = p->getDaug(3)->getP4();
121 int pid = EvtPDL::getStdHep(p->getDaug(0)->getId());
122 int charm;
123 if(pid == -321) charm = 1;
124 else charm = -1;
125 double m2, q2, cosV, cosL, chi;
126 KinVGen(_K, _pi, _e, _nu, charm, m2, q2, cosV, cosL, chi);
127 double _prob = calPDF(m2, q2, cosV, cosL, chi);
128 if(_prob>maxprob) {
129 maxprob=_prob;
130 std::cout << "Max PDF = " << ir << " charm= " << charm << " prob= " << _prob << std::endl;
131 }
132 }
133 std::cout << "Max!!!!!!!!!!! " << maxprob<< std::endl;
134*/
135
137 EvtVector4R K = p->getDaug(0)->getP4();
138 EvtVector4R pi = p->getDaug(1)->getP4();
139 EvtVector4R e = p->getDaug(2)->getP4();
140 EvtVector4R nu = p->getDaug(3)->getP4();
141
142 int pid = EvtPDL::getStdHep(p->getDaug(0)->getId());
143 int charm;
144 if(pid == -321) charm = 1;
145 else charm = -1;
146 double m2, q2, cosV, cosL, chi;
147 KinVGen(K, pi, e, nu, charm, m2, q2, cosV, cosL, chi);
148 double prob = calPDF(m2, q2, cosV, cosL, chi);
149 setProb(prob);
150 return;
151}
152
153void EvtDsToKKenu::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)
154{
155 EvtVector4R vp4_KPi = vp4_K + vp4_Pi;
156 EvtVector4R vp4_W = vp4_Lep + vp4_Nu;
157
158 m2 = vp4_KPi.mass2();
159 q2 = vp4_W.mass2();
160
161 EvtVector4R boost;
162 boost.set(vp4_KPi.get(0), -vp4_KPi.get(1), -vp4_KPi.get(2), -vp4_KPi.get(3));
163 EvtVector4R vp4_Kp = boostTo(vp4_K, boost);
164 cosV = vp4_Kp.dot(vp4_KPi)/(vp4_Kp.d3mag()*vp4_KPi.d3mag());
165
166 boost.set(vp4_W.get(0), -vp4_W.get(1), -vp4_W.get(2), -vp4_W.get(3));
167 EvtVector4R vp4_Lepp = boostTo(vp4_Lep, boost);
168 cosL = vp4_Lepp.dot(vp4_W)/(vp4_Lepp.d3mag()*vp4_W.d3mag());
169
170 EvtVector4R V = vp4_KPi/vp4_KPi.d3mag();
171 EvtVector4R C = vp4_Kp.cross(V);
172 C/=C.d3mag();
173 EvtVector4R D = vp4_Lepp.cross(V);
174 D/=D.d3mag();
175 double sinx = C.cross(V).dot(D);
176 double cosx = C.dot(D);
177 chi = sinx>0? acos(cosx): -acos(cosx);
178 if(charm==-1) chi=-chi;
179}
180
181double EvtDsToKKenu::calPDF(double m2, double q2, double cosV, double cosL, double chi) {
182 double delta[5] = {0};
183 double amplitude[5] = {0};
184 double m = sqrt(m2);
185 double q = sqrt(q2);
186
187 //begin to calculate form factor
188 EvtComplex F10(0.0,0.0);
189 EvtComplex F11(0.0,0.0);
190 EvtComplex F21(0.0,0.0);
191 EvtComplex F31(0.0,0.0);
192 EvtComplex F12(0.0,0.0);
193 EvtComplex F22(0.0,0.0);
194 EvtComplex F32(0.0,0.0);
195
196 EvtComplex f10(0.0,0.0);
197 EvtComplex f11(0.0,0.0);
198 EvtComplex f21(0.0,0.0);
199 EvtComplex f31(0.0,0.0);
200 EvtComplex f12(0.0,0.0);
201 EvtComplex f22(0.0,0.0);
202 EvtComplex f32(0.0,0.0);
203 EvtComplex coef(0.0,0.0);
204 double amplitude_temp, delta_temp;
205
206 for(int index=0; index<nAmps; index++)
207 {
208 switch(type[index])
209 {
210 //case 0: //calculate form factor of S wave
211 // {
212 // NRS(m, q, rS, rS1, a_delta, b_delta, mA, m0_1430_S, width0_1430_S, amplitude_temp, delta_temp, f10);
213 // amplitude[index] = amplitude_temp;
214 // delta[index] = delta_temp;
215 // F10 = F10 + f10;
216 // break;
217 // }
218 case 0://calculate Flatte
219 {
220 ResonanceSf0(m, q, mA, m0, g1, g2, amplitude_temp, delta_temp, f10);
221 amplitude[index] = amplitude_temp;
222 delta[index] = delta_temp;
223 coef = getCoef(rho_f0, phi_f0);
224 F10 = F10 + coef*f10;
225 break;
226
227 }
228 case 1: //calculate form factor of P wave (K*)
229 {
230 ResonanceP(m, q, mV, mA, V_0, A1_0, A2_0, m0, width0, rBW, amplitude_temp, delta_temp, f11, f21, f31);
231 amplitude[index] = amplitude_temp;
232 delta[index] = delta_temp;
233 coef = getCoef(rho, phi);
234 F11 = F11+ coef*f11;
235 F21 = F21+ coef*f21;
236 F31 = F31+ coef*f31;
237 break;
238 }
239 case 2: //calculate form factor of P wave (K*(1410))
240 {
241
242 ResonanceP(m, q, mV, mA, V_0, A1_0, A2_0, m0_1410, width0_1410, rBW, amplitude_temp, delta_temp, f11, f21, f31);
243 amplitude[index] = amplitude_temp;
244 delta[index] = delta_temp;
245 coef = getCoef(rho_1410, phi_1410);
246 F11 = F11+ coef*f11;
247 F21 = F21+ coef*f21;
248 F31 = F31+ coef*f31;
249 break;
250 }
251 case 3: //calculate form factor of D wave
252 {
253 ResonanceD(m, q, mV, mA, TV_0, T1_0, T2_0, m0_1430, width0_1430, rBW, amplitude_temp, delta_temp, f12, f22, f32);
254 amplitude[index] = amplitude_temp;
255 delta[index] = delta_temp;
256 coef = getCoef(rho_1430, phi_1430);
257 F12 = F12+ coef*f12;
258 F22 = F22+ coef*f22;
259 F32 = F32+ coef*f32;
260 break;
261 }
262 default:
263 {
264 std::cout<<"No such form factor type!!!"<<std::endl;
265 break;
266 }
267 }
268 }
269
270 //begin to calculate pdf value
271 double I,I1,I2,I3,I4,I5,I6,I7,I8,I9;
272
273 double cosV2 = cosV*cosV;
274 double sinV = sqrt(1.0-cosV2);
275 double sinV2 = sinV*sinV;
276
277 EvtComplex F1 = F10 + F11*cosV + F12*(1.5*cosV2-0.5);
278 EvtComplex F2 = F21*root1d2 + F22*cosV*root3d2;
279 EvtComplex F3 = F31*root1d2 + F32*cosV*root3d2;
280
281 I1 = 0.25*( abs2(F1) + 1.5*sinV2*(abs2(F2) + abs2(F3) ) );
282 I2 = -0.25*( abs2(F1) - 0.5*sinV2*(abs2(F2) + abs2(F3) ) );
283 I3 = -0.25*( abs2(F2) - abs2(F3) )*sinV2;
284 I4 = real( conj(F1)*F2 )*sinV*0.5;
285 I5 = real( conj(F1)*F3 )*sinV;
286 I6 = real( conj(F2)*F3 )*sinV2;
287 I7 = imag( conj(F2)*F1 )*sinV;
288 I8 = imag( conj(F3)*F1 )*sinV*0.5;
289 I9 = imag( conj(F3)*F2 )*sinV2*(-0.5);
290
291 double coschi = cos(chi);
292 double sinchi = sin(chi);
293 double sin2chi = 2.0*sinchi*coschi;
294 double cos2chi = 1.0 - 2.0*sinchi*sinchi;
295
296 double sinL = sqrt(1.-cosL*cosL);
297 double sinL2 = sinL*sinL;
298 double sin2L = 2.0*sinL*cosL;
299 double cos2L = 1.0 - 2.0*sinL2;
300
301 I = I1 + I2*cos2L + I3*sinL2*cos2chi + I4*sin2L*coschi + I5*sinL*coschi + I6*cosL + I7*sinL*sinchi + I8*sin2L*sinchi + I9*sinL2*sin2chi;
302 return I;
303}
304
305void EvtDsToKKenu::ResonanceP(double m, double q, double mV, double mA, double V_0, double A1_0, double A2_0, double m0, double width0, double rBW, double& amplitude, double& delta, EvtComplex& F11, EvtComplex& F21, EvtComplex& F31)
306{
307 double pKPi = getPStar(mD, m, q);
308 double mD2 = mD*mD;
309 double m2 = m*m;
310 double m02 = m0*m0;
311 double q2 = q*q;
312 double mV2 = mV*mV;
313 double mA2 = mA*mA;
314 double summDm = mD+m;
315 double V = V_0 /(1.0-q2/(mV2));
316 double A1 = A1_0/(1.0-q2/(mA2));
317 double A2 = A2_0/(1.0-q2/(mA2));
318 double A = summDm*A1;
319 double B = 2.0*mD*pKPi/summDm*V;
320
321 //construct the helicity form factor
322 double H0 = 0.5/(m*q)*((mD2-m2-q2)*summDm*A1-4.0*(mD2*pKPi*pKPi)/summDm*A2);
323 double Hp = A-B;
324 double Hm = A+B;
325
326 //calculate alpha
327 //double B_Kstar = 2./3.; //B_Kstar = Br(Kstar(892)->k pi)
328 double BF = 0.4920; //BF = Br(phi->k k)
329 double pStar0 = getPStar(m0, mK1, mK2);
330 //double alpha = sqrt(3.*Pi*B_Kstar/(pStar0*width0));
331 double alpha = sqrt(3.0*Pi*BF/(pStar0*width0));
332
333 //construct amplitudes of (non)resonance
334 double F = getF1(m, m0, mK1, mK2, rBW);
335 double width = getWidth1(m, m0, mK1, mK2, width0, rBW);
336
337 EvtComplex C(m0*width0*F,0.0);
338 double AA = m02-m2;
339 double BB = -m0*width;
340 EvtComplex amp = C/EvtComplex(AA,BB);
341 amplitude = abs(amp);
342 delta = atan2(imag(amp), real(amp));
343
344 double alpham2 = alpha*2.0;
345 F11 = amp*alpham2*q*H0*root2;
346 F21 = amp*alpham2*q*(Hp+Hm);
347 F31 = amp*alpham2*q*(Hp-Hm);
348}
349
350//void EvtDsToKKenu::NRS(double m, double q, double rS, double rS1, double a_delta, double b_delta, double mA, double m0, double width0, double& amplitude, double& delta, EvtComplex& F10)
351//{
352// static const double tmp = (mK+mK1)*(mK+mK2);
353//
354// double m2 = m*m;
355// double q2 = q*q;
356// double mA2 = mA*mA;
357// double pKPi = getPStar(mD, m, q);
358// double m_K0_1430 = m0;
359// double width_K0_1430 = width0;
360// double m2_K0_1430 = m_K0_1430*m_K0_1430;
361// double width = getWidth0(m, m_K0_1430, mK1, mK2, width_K0_1430);
362//
363// //calculate modul of the amplitude
364// double x,Pm;
365// if(m<m_K0_1430) {
366// x = sqrt(m2/tmp-1.0);
367// Pm = 1.0+rS1*x;
368// } else {
369// x = sqrt(m2_K0_1430/tmp-1.0);
370// Pm = 1.0+rS1*x;
371// Pm *= m_K0_1430*width_K0_1430/sqrt((m2_K0_1430-m2)*(m2_K0_1430-m2)+m2_K0_1430*width*width);
372// }
373//
374// //calculate phase of the amplitude
375// double pStar = getPStar(m, mK1, mK2);
376// double delta_bg = atan(2.*a_delta*pStar/(2.+a_delta*b_delta*pStar*pStar));
377// delta_bg = (delta_bg>0)? delta_bg: (delta_bg+Pi);
378//
379// double delta_K0_1430 = atan(m_K0_1430*width/(m2_K0_1430-m2));
380// delta_K0_1430 = (delta_K0_1430>0)? delta_K0_1430: (delta_K0_1430+Pi);
381// delta = delta_bg + delta_K0_1430;
382//
383// EvtComplex ci(cos(delta), sin(delta));
384// EvtComplex amp = ci*rS*Pm;
385// amplitude = rS*Pm;
386//
387// F10 = amp*pKPi*mD/(1.-q2/mA2);
388//}
389
390//LHCb
391// void EvtDsToKKenu::ResonanceSf0(double m, double q, EvtComplex& F10)
392//{
393// double pPiPi = getPStar(mD, m, q);
394// double m2 = m*m;
395// double q2 = q*q;
396// double mA2 = mA*mA;
397//
398// EvtComplex ciR(1.0, 0.0);
399// EvtComplex ciM(0.0, 1.0);
400//
401// EvtComplex rhopiPpiM = getrho(m2,mPi); //rho_pi+pi-
402// EvtComplex rhopi0pi0 = getrho(m2,mPi0); //rho_pi0pi0
403// EvtComplex rhoKPKM = getrho(m2,mK1); //rho_K+K-
404// EvtComplex rhoK0K0 = getrho(m2,mK0); //rho_K0K0
405//
406// EvtComplex rhopipi = (2.0/3.0)*rhopiPpiM + (1.0/3.0)*rhopi0pi0;
407// EvtComplex rhoKK = 0.5*rhoKPKM + 0.5*rhoK0K0;
408//
409// EvtComplex amp = ciR/(ciR*(m2f0-m2)-ciM*(g1*rhopipi+g2*rhoKK));
410// F10 = amp*pPiPi*mD/(1.0-q2/mA2);
411//}
412
413
414//BESII
415void EvtDsToKKenu::ResonanceSf0(double m, double q, double mA, double m0, double g1, double g2, double& amplitude, double& delta, EvtComplex& F10)
416{
417 double pKPi = getPStar(mD, m, q);
418 double m2 = m*m;
419 double q2 = q*q;
420 double mA2 = mA*mA;
421
422 EvtComplex ciR(1.0,0.0);
423 EvtComplex ciM(0.0,1.0);
424 EvtComplex bb3(0.0,0.0);
425
426 double aa3 = m2f0-m2;
427 bb3 = getrho(m2,mPi)*g1 + getrho(m2,mK1)*g2;
428
429 EvtComplex amp = ciR/(ciR*aa3-ciM*bb3);
430 amplitude = abs(amp);
431 delta = atan2(imag(amp), real(amp));
432
433 F10 = amp*pKPi*mD/(1.0-q2/mA2);
434
435}
436
437void EvtDsToKKenu::ResonanceD(double m, double q, double mV, double mA, double TV_0, double T1_0, double T2_0, double m0, double width0, double rBW, double& amplitude, double& delta, EvtComplex& F12, EvtComplex& F22, EvtComplex& F32)
438{
439 double pKPi = getPStar(mD, m, q);
440 double mD2 = mD*mD;
441 double m2 = m*m;
442 double m02 = m0*m0;
443 double q2 = q*q;
444 double mV2 = mV*mV;
445 double mA2 = mA*mA;
446 double summDm = mD+m;
447 double TV = TV_0 /(1.0-q2/(mV2));
448 double T1 = T1_0/(1.0-q2/(mA2));
449 double T2 = T2_0/(1.0-q2/(mA2));
450
451 //construct amplitudes of (non)resonance
452 double F = getF2(m, m0, mK1, mK2, rBW);
453 double width = getWidth2(m, m0, mK1, mK2, width0, rBW);
454 EvtComplex C(m0*width0*F,0.0);
455 double AA = m02-m2;
456 double BB = -m0*width;
457
458 EvtComplex amp = C/EvtComplex(AA,BB);
459 amplitude = abs(amp);
460 delta = atan2(imag(amp), real(amp));
461
462 F12 = amp*mD*pKPi/3.*((mD2-m2-q2)*summDm*T1-mD2*pKPi*pKPi/summDm*T2);
463 F22 = amp*root2d3*mD*m*q*pKPi*summDm*T1;
464 F32 = amp*root2d3*2.*mD2*m*q*pKPi*pKPi/summDm*TV;
465}
466
467
468EvtComplex EvtDsToKKenu::getrho(double m2, double mX)
469{
470 EvtComplex rho(0.0, 0.0);
471 EvtComplex ciR(1.0, 0.0);
472 EvtComplex ciM(0.0, 1.0);
473 if ((1.0-4.0*mX*mX/m2)>0)
474 rho = ciR*sqrt(1.0-4.0*mX*mX/m2);
475 else rho = ciM*sqrt(4.0*mX*mX/m2-1.0);
476 return rho;
477}
478
479double EvtDsToKKenu::getPStar(double m, double m1, double m2)
480{
481 double s = m*m;
482 double s1 = m1*m1;
483 double s2 = m2*m2;
484 double x = s + s1 - s2;
485 double t = 0.25*x*x/s - s1;
486 double p;
487 if (t>0.0) {
488 p = sqrt(t);
489 } else {
490 //std::cout << " Hello, pstar is less than 0.0" << std::endl;
491 p = 0.04;
492 }
493 return p;
494}
495
496double EvtDsToKKenu::getF1(double m, double m0, double m_c1, double m_c2, double rBW)
497{
498 double pStar = getPStar(m, m_c1, m_c2);
499 double pStar0 = getPStar(m0, m_c1, m_c2);
500 double rBW2 = rBW*rBW;
501 double pStar2 = pStar*pStar;
502 double pStar02 = pStar0*pStar0;
503 double B = 1./sqrt(1.+rBW2*pStar2);
504 double B0 = 1./sqrt(1.+rBW2*pStar02);
505 double F = pStar/pStar0*B/B0;
506 return F;
507}
508
509double EvtDsToKKenu::getF2(double m, double m0, double m_c1, double m_c2, double rBW)
510{
511 double pStar = getPStar(m, m_c1, m_c2);
512 double pStar0 = getPStar(m0, m_c1, m_c2);
513 double rBW2 = rBW*rBW;
514 double pStar2 = pStar*pStar;
515 double pStar02 = pStar0*pStar0;
516 double B = 1./sqrt((rBW2*pStar2-3.)*(rBW2*pStar2-3.)+9.*rBW2*pStar2);
517 double B0 = 1./sqrt((rBW2*pStar02-3.)*(rBW2*pStar02-3.)+9.*rBW2*pStar02);
518 double F = pStar2/pStar02*B/B0;
519 return F;
520}
521
522double EvtDsToKKenu::getWidth0(double m, double m0, double m_c1, double m_c2, double width0)
523{
524 double pStar = getPStar(m, m_c1, m_c2);
525 double pStar0 = getPStar(m0, m_c1, m_c2);
526 double width = width0*pStar/pStar0*m0/m;
527 return width;
528}
529
530double EvtDsToKKenu::getWidth1(double m, double m0, double m_c1, double m_c2, double width0, double rBW)
531{
532 double pStar = getPStar(m, m_c1, m_c2);
533 double pStar0 = getPStar(m0, m_c1, m_c2);
534 double F = getF1(m, m0, m_c1, m_c2, rBW);
535 double width = width0*pStar/pStar0*m0/m*F*F;
536 return width;
537}
538
539double EvtDsToKKenu::getWidth2(double m, double m0, double m_c1, double m_c2, double width0, double rBW)
540{
541 double pStar = getPStar(m, m_c1, m_c2);
542 double pStar0 = getPStar(m0, m_c1, m_c2);
543 double F = getF2(m, m0, m_c1, m_c2, rBW);
544 double width = width0*pStar/pStar0*m0/m*F*F;
545 return width;
546}
547
548EvtComplex EvtDsToKKenu::getCoef(double rho, double phi)
549{
550 EvtComplex coef ( rho*cos(phi), rho*sin(phi) );
551 return coef;
552}
double sin(const BesAngle a)
Definition BesAngle.h:210
double cos(const BesAngle a)
Definition BesAngle.h:213
const double delta
TF1 * g1
Double_t x[10]
const DifComplex I
Evt3Rank3C conj(const Evt3Rank3C &t2)
double imag(const EvtComplex &c)
double abs2(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 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)
void getName(std::string &name)
EvtDecayBase * clone()
virtual ~EvtDsToKKenu()
void decay(EvtParticle *p)
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
const float pi
Definition vector3.h:133