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