BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
D0ToKLpipi.cxx
Go to the documentation of this file.
1#include "D0ToKLpipi.h"
2#include "TMath.h"
3#include <complex>
4#include <iostream>
5#include <math.h>
6#include <stdlib.h>
7#include <string>
8#include <vector>
9
10#include "CLHEP/Matrix/Matrix.h"
11#include "CLHEP/Matrix/SymMatrix.h"
12#include "CLHEP/Matrix/Vector.h"
13#include "CLHEP/Random/RandFlat.h"
14#include "CLHEP/Vector/LorentzVector.h"
15#include "CLHEP/Vector/ThreeVector.h"
16#include "CLHEP/Vector/TwoVector.h"
17using CLHEP::Hep2Vector;
18using CLHEP::Hep3Vector;
19using CLHEP::HepLorentzVector;
20using CLHEP::HepVector;
21
22using namespace std; //::endl;
23
25
27 // std::cout << "D0ToKLpipi ==> Initialization !" << std::endl;
28
29 _nd = 3;
30 tan2thetaC = ( 0.22650 * 0.22650 ) /
31 ( 1. - ( 0.22650 * 0.22650 ) ); // sin(theta_C) = 0.22650 +/- 0.00048
32 pi180inv = 1.0 * 3.1415926 / 180;
33
34 // spin[0] = 1 ; //Rho770
35 // spin[1] = 1 ; //Omega782
36 // spin[2] = 2 ; //ftwo1270
37 // spin[3] = 1 ; //Rho1450
38 // spin[4] = 1 ; //Kstar892minus
39 // spin[5] = 2 ; //K2star1430minus
40 // spin[6] = 1 ; //Kstar1680minus
41 // spin[7] = 1 ; //Kstar1410minus
42 // spin[8] = 1 ; //Kstar892plus
43 // spin[9] = 2 ; //K2star1430plus
44 // spin[10]= 1 ; //Kstar1410plus
45 mass_R[0] = 0.77526;
46 width_R[0] = 0.14740;
47 spin_R[0] = 1;
48 ar[0] = 1;
49 phir[0] = 0;
50 mass_R[1] = 0.78266;
51 width_R[1] = 0.00868;
52 spin_R[1] = 1;
53 ar[1] = 0.03928;
54 phir[1] = 106.1;
55 mass_R[2] = 1.27550;
56 width_R[2] = 0.18670;
57 spin_R[2] = 2;
58 ar[2] = 1.302;
59 phir[2] = -39.97;
60 mass_R[3] = 1.46500;
61 width_R[3] = 0.40000;
62 spin_R[3] = 1;
63 ar[3] = 1.656;
64 phir[3] = 113;
65 mass_R[4] = 0.89167;
66 width_R[4] = 0.0514;
67 spin_R[4] = 1;
68 ar[4] = 1.843;
69 phir[4] = 138;
70 mass_R[5] = 1.42730;
71 width_R[5] = 0.10000;
72 spin_R[5] = 2;
73 ar[5] = 1.537;
74 phir[5] = -49.45;
75 mass_R[6] = 1.71800;
76 width_R[6] = 0.3220;
77 spin_R[6] = 1;
78 ar[6] = 2.251;
79 phir[6] = -196.3;
80 mass_R[7] = 1.41400;
81 width_R[7] = 0.2320;
82 spin_R[7] = 1;
83 ar[7] = 0.4733;
84 phir[7] = 1659;
85 mass_R[8] = 0.89167;
86 width_R[8] = 0.0514;
87 spin_R[8] = 1;
88 ar[8] = 0.1617;
89 phir[8] = -36.25;
90 mass_R[9] = 1.42730;
91 width_R[9] = 0.1000;
92 spin_R[9] = 2;
93 ar[9] = 0.2267;
94 phir[9] = -84.64;
95 mass_R[10] = 1.41400;
96 width_R[10] = 0.2320;
97 spin_R[10] = 1;
98 ar[10] = 0.2185;
99 phir[10] = 73.12;
100 mass_R[11] = 1.42500;
101 width_R[11] = 0.2700;
102 spin_R[11] = 1;
103 ar[11] = 2.396;
104 phir[11] = 98.14;
105 // beta[kb] = EvtComplex(mag*cos(phase*pi180inv), mag*sin(phase*pi180inv)) ;
106 // fprod[kf] = EvtComplex(mag*cos(phase*pi180inv), mag*sin(phase*pi180inv)) ;
107 beta[0] = complex<double>( 8.5 * cos( 68.5 * pi180inv ), 8.5 * sin( 68.5 * pi180inv ) );
108 beta[1] = complex<double>( 12.2 * cos( 24 * pi180inv ), 12.2 * sin( 24 * pi180inv ) );
109 beta[2] = complex<double>( 29.2 * cos( -0.1 * pi180inv ), 29.2 * sin( -0.1 * pi180inv ) );
110 beta[3] = complex<double>( 10.8 * cos( -51.9 * pi180inv ), 10.8 * sin( -51.9 * pi180inv ) );
111 beta[4] = complex<double>( 0., 0. );
112
113 fprod[0] = complex<double>( 8 * cos( -126 * pi180inv ), 8 * sin( -126 * pi180inv ) );
114 fprod[1] =
115 complex<double>( 26.3 * cos( -152.3 * pi180inv ), 26.3 * sin( -152.3 * pi180inv ) );
116 fprod[2] = complex<double>( 33 * cos( -93.2 * pi180inv ), 33 * sin( -93.2 * pi180inv ) );
117 fprod[3] =
118 complex<double>( 26.2 * cos( -121.4 * pi180inv ), 26.2 * sin( -121.4 * pi180inv ) );
119 fprod[4] = complex<double>( 0., 0. );
120
121 // CP_mult[w] = (EvtComplex(1., 0.) - 2.*tan2thetaC*EvtComplex(r*cos(delta), r*sin(delta))) ;
122 CP_mult[0] =
123 complex<double>( 1., 0. ) - 2. * tan2thetaC *
124 complex<double>( 1.851 * cos( -94.07 * pi180inv ),
125 1.851 * sin( -94.07 * pi180inv ) );
126 CP_mult[1] =
127 complex<double>( 1., 0. ) -
128 2. * tan2thetaC *
129 complex<double>( 6.332 * cos( 2.103 * pi180inv ), 6.332 * sin( 2.103 * pi180inv ) );
130 CP_mult[2] =
131 complex<double>( 1., 0. ) - 2. * tan2thetaC *
132 complex<double>( 3.229 * cos( -60.05 * pi180inv ),
133 3.229 * sin( -60.05 * pi180inv ) );
134 CP_mult[3] =
135 complex<double>( 1., 0. ) -
136 2. * tan2thetaC *
137 complex<double>( 12.75 * cos( 73.62 * pi180inv ), 12.75 * sin( 73.62 * pi180inv ) );
138 CP_mult[4] =
139 complex<double>( 1., 0. ) - 2. * tan2thetaC *
140 complex<double>( 0.7116 * cos( -177.149 * pi180inv ),
141 0.7116 * sin( -177.149 * pi180inv ) );
142
143 ma[0] = 0.651;
144 g[0][0] = 0.22889;
145 g[0][1] = -0.55377;
146 g[0][2] = 0;
147 g[0][3] = -0.39899;
148 g[0][4] = -0.34639;
149 ma[1] = 1.20360;
150 g[1][0] = 0.94128;
151 g[1][1] = 0.55095;
152 g[1][2] = 0;
153 g[1][3] = 0.39065;
154 g[1][4] = 0.31503;
155 ma[2] = 1.55817;
156 g[2][0] = 0.36856;
157 g[2][1] = 0.23888;
158 g[2][2] = 0.55639;
159 g[2][3] = 0.18340;
160 g[2][4] = 0.18681;
161 ma[3] = 1.21000;
162 g[3][0] = 0.33650;
163 g[3][1] = 0.40907;
164 g[3][2] = 0.85679;
165 g[3][3] = 0.19906;
166 g[3][4] = -0.00984;
167 ma[4] = 1.82206;
168 g[4][0] = 0.18171;
169 g[4][1] = -0.17558;
170 g[4][2] = -0.79658;
171 g[4][3] = -0.00355;
172 g[4][4] = 0.22358;
173
174 // Hadronic parameters for tag modes: 0=no-specified, 1=Kpi, 2=Kpipi0, 3=K3pi
175 rd[0] = 0.0;
176 rd[1] = 0.0586;
177 rd[2] = 0.0440;
178 rd[3] = 0.0546;
179 deltad[0] = 0.0;
180 deltad[1] = 194.7 * pi180inv;
181 deltad[2] = 196.0 * pi180inv;
182 deltad[3] = 167.0 * pi180inv;
183 Rf[0] = 0.0;
184 Rf[1] = 1.0;
185 Rf[2] = 0.78;
186 Rf[3] = 0.52;
187
188 return;
189}
190
191complex<double> D0ToKLpipi::Amp_PFT( vector<double> k0l, vector<double> pip,
192 vector<double> pim ) {
193
194 vector<double> pD;
195 pD.clear();
196 if ( k0l.size() != 4 || pip.size() != 4 || pim.size() != 4 )
197 cout << "ERROR in KSPIPI daughter 4 momentum" << endl;
198 for ( int i = 0; i < k0l.size(); i++ ) { pD.push_back( k0l[i] + pip[i] + pim[i] ); }
199
200 // EvtResonance2 DK2piRes0(pD, pip, pim, ar[0], phir[0], width_R[0], mass_R[0], spin[0]);
201 // //ar, phir, width, mass, spin //Rho770 EvtResonance2 DK2piRes1(pD, pip, pim, ar[1],
202 // phir[1], width_R[1], mass_R[1], spin[1]); //ar, phir, width, mass, spin //Omega782
203 // EvtResonance2 DK2piRes2(pD, pip, pim, ar[2], phir[2], width_R[2], mass_R[2], spin[2]);
204 // //ar, phir, width, mass, spin //ftwo1270 EvtResonance2 DK2piRes3(pD, pip, pim, ar[3],
205 // phir[3], width_R[3], mass_R[3], spin[3]); //ar, phir, width, mass, spin //Rho1450
206 // EvtResonance2 DK2piRes4(pD, k0l, pim, ar[4], phir[4], width_R[4], mass_R[4], spin[4]);
207 // //ar, phir, width, mass, spin //Kstar892minus EvtResonance2 DK2piRes5(pD, k0l, pim,
208 // ar[5], phir[5], width_R[5], mass_R[5], spin[5]); //ar, phir, width, mass, spin
209 // //K2star1430minus EvtResonance2 DK2piRes6(pD, k0l, pim, ar[6], phir[6], width_R[6],
210 // mass_R[6], spin[6]); //ar, phir, width, mass, spin //Kstar1680minus EvtResonance2
211 // DK2piRes7(pD, k0l, pim, ar[7], phir[7], width_R[7], mass_R[7], spin[7]); //ar, phir,
212 // width, mass, spin //Kstar1410minus EvtResonance2 DK2piRes8(pD, k0l, pip, ar[8], phir[8],
213 // width_R[8], mass_R[8], spin[8]); //ar, phir, width, mass, spin //Kstar892plus
214 // EvtResonance2 DK2piRes9(pD, k0l, pip, ar[9], phir[9], width_R[9], mass_R[9], spin[9]);
215 // //ar, phir, width, mass, spin //K2star1430plus EvtResonance2 DK2piRes10(pD, k0l, pip,
216 // ar[10], phir[10], width_R[10], mass_R[10], spin[10]); //ar, phir, width, mass, spin
217 // //Kstar1410plus
218 complex<double> DK2piRes0 = Resonance2( pD, pip, pim, ar[0], phir[0], width_R[0], mass_R[0],
219 spin_R[0] ); // ar, phir, width, mass, spin Rho770
220 complex<double> DK2piRes1 = Resonance2( pD, pip, pim, ar[1], phir[1], width_R[1], mass_R[1],
221 spin_R[1] ); // ar, phir, width, mass, spin Omega782
222 complex<double> DK2piRes2 = Resonance2( pD, pip, pim, ar[2], phir[2], width_R[2], mass_R[2],
223 spin_R[2] ); // ar, phir, width, mass, spin ftwo1270
224 complex<double> DK2piRes3 = Resonance2( pD, pip, pim, ar[3], phir[3], width_R[3], mass_R[3],
225 spin_R[3] ); // ar, phir, width, mass, spin Rho1450
226 complex<double> DK2piRes4 = Resonance2( pD, k0l, pim, ar[4], phir[4], width_R[4], mass_R[4],
227 spin_R[4] ); // ar, phir, width, mass, spin Kstar892-
228 complex<double> DK2piRes5 =
229 Resonance2( pD, k0l, pim, ar[5], phir[5], width_R[5], mass_R[5],
230 spin_R[5] ); // ar, phir, width, mass, spin K2star1430-
231 complex<double> DK2piRes6 =
232 Resonance2( pD, k0l, pim, ar[6], phir[6], width_R[6], mass_R[6],
233 spin_R[6] ); // ar, phir, width, mass, spin Kstar1680-
234 complex<double> DK2piRes7 =
235 Resonance2( pD, k0l, pim, ar[7], phir[7], width_R[7], mass_R[7],
236 spin_R[7] ); // ar, phir, width, mass, spin Kstar1410-
237 complex<double> DK2piRes8 = Resonance2( pD, k0l, pip, ar[8], phir[8], width_R[8], mass_R[8],
238 spin_R[8] ); // ar, phir, width, mass, spin Kstar892+
239 complex<double> DK2piRes9 =
240 Resonance2( pD, k0l, pip, ar[9], phir[9], width_R[9], mass_R[9],
241 spin_R[9] ); // ar, phir, width, mass, spin K2star1430+
242 complex<double> DK2piRes10 =
243 Resonance2( pD, k0l, pip, ar[10], phir[10], width_R[10], mass_R[10],
244 spin_R[10] ); // ar, phir, width, mass, spin Kstar1410+
245
246 complex<double> pipi_s_wave = K_matrix( pip, pim );
247 if ( pipi_s_wave == complex<double>( 9999., 9999. ) ) return 1e-20;
248
249 complex<double> kpi_s_wave =
250 amplitude_LASS( k0l, pip, pim, "k0lpim", ar[11], phir[11] * pi180inv );
251 // complex<double> kpi_s_wave_dcs = amplitude_LASS(k0l, pip, pim, "k0spip", ar[12],
252 // phir[12]*pi180inv) ; should be there but not observed yet GUESS
253
254 complex<double> TOT_PFT_AMP = DK2piRes0 * CP_mult[0] + DK2piRes1 * CP_mult[1] +
255 DK2piRes2 * CP_mult[2] + DK2piRes3 * CP_mult[3] + DK2piRes4 +
256 DK2piRes5 + DK2piRes6 + DK2piRes7 + DK2piRes8 * ( -1. ) +
257 DK2piRes9 * ( -1. ) + DK2piRes10 * ( -1. ) +
258 pipi_s_wave * CP_mult[4] + kpi_s_wave;
259
260 return TOT_PFT_AMP;
261}
262
263complex<double> D0ToKLpipi::Resonance2( vector<double> p4_p, vector<double> p4_d1,
264 vector<double> p4_d2, double mag, double theta,
265 double gamma, double bwm, int spin ) {
266
267 complex<double> ampl;
268
269 // EvtVector4R p4_d3 = p4_p - p4_d1 - p4_d2;
270 HepLorentzVector _p4_p;
271 _p4_p.setX( p4_p[0] );
272 _p4_p.setY( p4_p[1] );
273 _p4_p.setZ( p4_p[2] );
274 _p4_p.setT( p4_p[3] );
275 HepLorentzVector _p4_d1;
276 _p4_d1.setX( p4_d1[0] );
277 _p4_d1.setY( p4_d1[1] );
278 _p4_d1.setZ( p4_d1[2] );
279 _p4_d1.setT( p4_d1[3] );
280 HepLorentzVector _p4_d2;
281 _p4_d2.setX( p4_d2[0] );
282 _p4_d2.setY( p4_d2[1] );
283 _p4_d2.setZ( p4_d2[2] );
284 _p4_d2.setT( p4_d2[3] );
285 HepLorentzVector _p4_d3 = _p4_p - _p4_d1 - _p4_d2;
286
287 double mAB = ( _p4_d1 + _p4_d2 ).invariantMass();
288 double mBC = ( _p4_d2 + _p4_d3 ).invariantMass();
289 double mAC = ( _p4_d1 + _p4_d3 ).invariantMass();
290 double mA = _p4_d1.invariantMass();
291 double mB = _p4_d2.invariantMass();
292 double mD = _p4_p.invariantMass();
293 double mC = _p4_d3.invariantMass();
294
295 double mR = bwm;
296 double gammaR = gamma;
297 double pAB =
298 sqrt( ( ( ( mAB * mAB - mA * mA - mB * mB ) * ( mAB * mAB - mA * mA - mB * mB ) / 4.0 ) -
299 mA * mA * mB * mB ) /
300 ( mAB * mAB ) );
301 double pR =
302 sqrt( ( ( ( mR * mR - mA * mA - mB * mB ) * ( mR * mR - mA * mA - mB * mB ) / 4.0 ) -
303 mA * mA * mB * mB ) /
304 ( mR * mR ) );
305
306 double pD = ( ( ( mD * mD - mR * mR - mC * mC ) * ( mD * mD - mR * mR - mC * mC ) / 4.0 ) -
307 mR * mR * mC * mC ) /
308 ( mD * mD );
309 if ( pD > 0 ) { pD = sqrt( pD ); }
310 else { pD = 0; }
311 double pDAB =
312 sqrt( ( ( ( mD * mD - mAB * mAB - mC * mC ) * ( mD * mD - mAB * mAB - mC * mC ) / 4.0 ) -
313 mAB * mAB * mC * mC ) /
314 ( mD * mD ) );
315 double fR = 1;
316 double fD = 1;
317 int power = 0;
318 switch ( spin )
319 {
320 case 0:
321 fR = 1.0;
322 fD = 1.0;
323 power = 1;
324 break;
325 case 1:
326 fR = sqrt( 1.0 + 1.5 * 1.5 * pR * pR ) / sqrt( 1.0 + 1.5 * 1.5 * pAB * pAB );
327 fD = sqrt( 1.0 + 5.0 * 5.0 * pD * pD ) / sqrt( 1.0 + 5.0 * 5.0 * pDAB * pDAB );
328 power = 3;
329 break;
330 case 2:
331 fR = sqrt( ( 9 + 3 * pow( ( 1.5 * pR ), 2 ) + pow( ( 1.5 * pR ), 4 ) ) /
332 ( 9 + 3 * pow( ( 1.5 * pAB ), 2 ) + pow( ( 1.5 * pAB ), 4 ) ) );
333 fD = sqrt( ( 9 + 3 * pow( ( 5.0 * pD ), 2 ) + pow( ( 5.0 * pD ), 4 ) ) /
334 ( 9 + 3 * pow( ( 5.0 * pDAB ), 2 ) + pow( ( 5.0 * pDAB ), 4 ) ) );
335 power = 5;
336 break;
337 default: cout << "Incorrect spin in D0ToKLpipi::EvtResonance2.cc\n" << endl;
338 }
339
340 double gammaAB = gammaR * pow( pAB / pR, power ) * ( mR / mAB ) * fR * fR;
341 switch ( spin )
342 {
343 case 0:
344 ampl = mag * complex<double>( cos( theta * pi180inv ), sin( theta * pi180inv ) ) * fR *
345 fD / ( mR * mR - mAB * mAB - complex<double>( 0.0, mR * gammaAB ) );
346 break;
347 case 1:
348 ampl = mag * complex<double>( cos( theta * pi180inv ), sin( theta * pi180inv ) ) *
349 ( fR * fD *
350 ( mAC * mAC - mBC * mBC +
351 ( ( mD * mD - mC * mC ) * ( mB * mB - mA * mA ) / ( mAB * mAB ) ) ) /
352 ( mR * mR - mAB * mAB - complex<double>( 0.0, mR * gammaAB ) ) );
353 break;
354 case 2:
355 ampl = mag * complex<double>( cos( theta * pi180inv ), sin( theta * pi180inv ) ) *
356 ( fR * fD / ( mR * mR - mAB * mAB - complex<double>( 0.0, mR * gammaAB ) ) ) *
357 ( pow( ( mBC * mBC - mAC * mAC +
358 ( mD * mD - mC * mC ) * ( mA * mA - mB * mB ) / ( mAB * mAB ) ),
359 2 ) -
360 ( 1.0 / 3.0 ) *
361 ( mAB * mAB - 2 * mD * mD - 2 * mC * mC +
362 pow( ( mD * mD - mC * mC ) / mAB, 2 ) ) *
363 ( mAB * mAB - 2 * mA * mA - 2 * mB * mB +
364 pow( ( mA * mA - mB * mB ) / mAB, 2 ) ) );
365 break;
366 default: cout << "Incorrect spin in D0ToKSpipi::Resonance2.cc\n" << endl;
367 }
368
369 return ampl;
370}
371
372complex<double> D0ToKLpipi::K_matrix( vector<double> p_pip, vector<double> p_pim ) {
373
374 // double pi180inv = 1.0/EvtConst::radToDegrees;
375 bool reject = false;
376
377 const double mD0 = 1.86483;
378 const double mKl = 0.49761;
379 const double mPi = 0.13957;
380
381 HepLorentzVector _p_pip( p_pip[0], p_pip[1], p_pip[2], p_pip[3] );
382 HepLorentzVector _p_pim( p_pim[0], p_pim[1], p_pim[2], p_pim[3] );
383
384 double mAB = ( _p_pip + _p_pim ).m();
385
386 double s = mAB * mAB;
387
388 // Define the complex coupling constants
389 // Double_t g[5][5]; // g[Physical pole]Decay channel]
390
391 // pi+pi- channel
392 // g[0][0]=0.22889;
393 // g[1][0]=0.94128;
394 // g[2][0]=0.36856;
395 // g[3][0]=0.33650;
396 // g[4][0]=0.18171;
397
398 //// K+K- channel
399 // g[0][1]=-0.55377;
400 // g[1][1]=0.55095;
401 // g[2][1]=0.23888;
402 // g[3][1]=0.40907;
403 // g[4][1]=-0.17558;
404
405 //// 4pi channel
406 // g[0][2]=0;
407 // g[1][2]=0;
408 // g[2][2]=0.55639;
409 // g[3][2]=0.85679;
410 // g[4][2]=-0.79658;
411
412 //// eta eta channel
413 // g[0][3]=-0.39899;
414 // g[1][3]=0.39065;
415 // g[2][3]=0.18340;
416 // g[3][3]=0.19906;
417 // g[4][3]=-0.00355;
418
419 ////eta eta' channel
420 // g[0][4]=-0.34639;
421 // g[1][4]=0.31503;
422 // g[2][4]=0.18681;
423 // g[3][4]=-0.00984;
424 // g[4][4]=0.22358;
425
426 // Define masses of the physical poles (in GeV)
427 // Double_t ma[5];
428
429 // ma[0]=0.651;
430 // ma[1]=1.20360;
431 // ma[2]=1.55817;
432 // ma[3]=1.21000;
433 // ma[4]=1.82206;
434
435 // Define variables
436 complex<double> n11, n12, n13, n14, n15, n21, n22, n23, n24, n25, n31, n32, n33, n34, n35,
437 n41, n42, n43, n44, n45, n51, n52, n53, n54, n55;
438 double rho1sq, rho2sq, rho4sq, rho5sq;
439 complex<double> rho1, rho2, rho3, rho4, rho5;
440 complex<double> rho[5];
441 complex<double> pole, SVT, Adler;
442 complex<double> det;
443 complex<double> i[5][5];
444 double f[5][5];
445
446 // pi+, K+, eta, and eta' PDG masses
447 double mpi = 0.13957;
448 double mK = 0.493677;
449 double meta = 0.54775;
450 double metap = 0.95778;
451
452 // Init matrices and vectors with zeros
453 complex<double> K[5][5];
454 for ( Int_t k = 0; k < 5; k++ )
455 {
456 for ( Int_t l = 0; l < 5; l++ )
457 {
458 i[k][l] = complex<double>( 0., 0. );
459 K[k][l] = complex<double>( 0., 0. );
460 f[k][l] = 0.;
461 }
462 rho[k] = 0.;
463 }
464
465 // Fill scattering data values
466 Double_t s_scatt = -3.92637;
467 Double_t sa = 1.0;
468 Double_t sa_0 = -0.15;
469
470 // f_scattering
471 f[0][0] = 0.23399;
472 f[0][1] = 0.15044;
473 f[0][2] = -0.20545;
474 f[0][3] = 0.32825;
475 f[0][4] = 0.35412;
476
477 f[1][0] = f[0][1];
478 f[2][0] = f[0][2];
479 f[3][0] = f[0][3];
480 f[4][0] = f[0][4];
481
482 // Compute phase space factors
483 rho1sq = ( 1.0 - ( pow( ( mpi + mpi ), 2 ) / s ) );
484 if ( rho1sq >= 0. ) { rho1 = complex<double>( sqrt( rho1sq ), 0. ); }
485 else { rho1 = complex<double>( 0., sqrt( -rho1sq ) ); }
486 rho[0] = rho1;
487
488 rho2sq = ( 1.0 - ( pow( ( mK + mK ), 2 ) / s ) );
489 if ( rho2sq >= 0. ) { rho2 = complex<double>( sqrt( rho2sq ), 0. ); }
490 else { rho2 = complex<double>( 0., sqrt( -rho2sq ) ); }
491
492 rho[1] = rho2;
493
494 rho3 = complex<double>( 0., 0. );
495
496 if ( s <= 1 )
497 {
498 Double_t real = 1.2274 + 0.00370909 / ( s * s ) - ( 0.111203 ) / (s)-6.39017 * s +
499 16.8358 * s * s - 21.8845 * s * s * s + 11.3153 * s * s * s * s;
500 Double_t cont32 = sqrt( 1.0 - ( 16.0 * mpi * mpi ) );
501 rho3 = complex<double>( cont32 * real, 0. );
502 }
503 else { rho3 = complex<double>( sqrt( 1.0 - ( 16.0 * mpi * mpi / s ) ), 0. ); }
504 rho[2] = rho3;
505
506 rho4sq = ( 1.0 - ( pow( ( meta + meta ), 2 ) / s ) );
507 if ( rho4sq >= 0. ) { rho4 = complex<double>( sqrt( rho4sq ), 0. ); }
508 else { rho4 = complex<double>( 0., sqrt( -rho4sq ) ); }
509 rho[3] = rho4;
510
511 rho5sq = ( 1.0 - ( pow( ( meta + metap ), 2 ) / s ) );
512 if ( rho5sq >= 0. ) { rho5 = complex<double>( sqrt( rho5sq ), 0. ); }
513 else { rho5 = complex<double>( 0., sqrt( -rho5sq ) ); }
514 rho[4] = rho5;
515
516 // Sum over the poles
517 for ( Int_t k = 0; k < 5; k++ )
518 {
519 for ( Int_t l = 0; l < 5; l++ )
520 {
521 for ( Int_t pole_index = 0; pole_index < 5; pole_index++ )
522 {
523 Double_t A = g[pole_index][k] * g[pole_index][l];
524 Double_t B = ma[pole_index] * ma[pole_index] - s;
525 K[k][l] = K[k][l] + complex<double>( A / B, 0. );
526 }
527 }
528 }
529
530 for ( Int_t k = 0; k < 5; k++ )
531 {
532 for ( Int_t l = 0; l < 5; l++ )
533 {
534 Double_t C = f[k][l] * ( 1.0 - s_scatt );
535 Double_t D = ( s - s_scatt );
536 K[k][l] = K[k][l] + complex<double>( C / D, 0. );
537 }
538 }
539
540 for ( Int_t k = 0; k < 5; k++ )
541 {
542 for ( Int_t l = 0; l < 5; l++ )
543 {
544 Double_t E = ( s - ( sa * mpi * mpi * 0.5 ) ) * ( 1.0 - sa_0 );
545 Double_t F = ( s - sa_0 );
546 K[k][l] = K[k][l] * complex<double>( E / F, 0. );
547 }
548 }
549
550 n11 = complex<double>( 1., 0. ) - complex<double>( 0., 1. ) * K[0][0] * rho[0];
551 n12 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[0][1] * rho[1];
552 n13 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[0][2] * rho[2];
553 n14 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[0][3] * rho[3];
554 n15 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[0][4] * rho[4];
555
556 n21 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[1][0] * rho[0];
557 n22 = complex<double>( 1., 0. ) - complex<double>( 0., 1. ) * K[1][1] * rho[1];
558 n23 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[1][2] * rho[2];
559 n24 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[1][3] * rho[3];
560 n25 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[1][4] * rho[4];
561
562 n31 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[2][0] * rho[0];
563 n32 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[2][1] * rho[1];
564 n33 = complex<double>( 1., 0. ) - complex<double>( 0., 1. ) * K[2][2] * rho[2];
565 n34 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[2][3] * rho[3];
566 n35 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[2][4] * rho[4];
567
568 n41 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[3][0] * rho[0];
569 n42 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[3][1] * rho[1];
570 n43 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[3][2] * rho[2];
571 n44 = complex<double>( 1., 0. ) - complex<double>( 0., 1. ) * K[3][3] * rho[3];
572 n45 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[3][4] * rho[4];
573
574 n51 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[4][0] * rho[0];
575 n52 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[4][1] * rho[1];
576 n53 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[4][2] * rho[2];
577 n54 = complex<double>( 0., 0. ) - complex<double>( 0., 1. ) * K[4][3] * rho[3];
578 n55 = complex<double>( 1., 0. ) - complex<double>( 0., 1. ) * K[4][4] * rho[4];
579
580 // Compute the determinant
581 det = ( n15 * n24 * n33 * n42 * n51 - n14 * n25 * n33 * n42 * n51 -
582 n15 * n23 * n34 * n42 * n51 + n13 * n25 * n34 * n42 * n51 +
583 n14 * n23 * n35 * n42 * n51 - n13 * n24 * n35 * n42 * n51 -
584 n15 * n24 * n32 * n43 * n51 + n14 * n25 * n32 * n43 * n51 +
585 n15 * n22 * n34 * n43 * n51 - n12 * n25 * n34 * n43 * n51 -
586 n14 * n22 * n35 * n43 * n51 + n12 * n24 * n35 * n43 * n51 +
587 n15 * n23 * n32 * n44 * n51 - n13 * n25 * n32 * n44 * n51 -
588 n15 * n22 * n33 * n44 * n51 + n12 * n25 * n33 * n44 * n51 +
589 n13 * n22 * n35 * n44 * n51 - n12 * n23 * n35 * n44 * n51 -
590 n14 * n23 * n32 * n45 * n51 + n13 * n24 * n32 * n45 * n51 +
591 n14 * n22 * n33 * n45 * n51 - n12 * n24 * n33 * n45 * n51 -
592 n13 * n22 * n34 * n45 * n51 + n12 * n23 * n34 * n45 * n51 -
593 n15 * n24 * n33 * n41 * n52 + n14 * n25 * n33 * n41 * n52 +
594 n15 * n23 * n34 * n41 * n52 - n13 * n25 * n34 * n41 * n52 -
595 n14 * n23 * n35 * n41 * n52 + n13 * n24 * n35 * n41 * n52 +
596 n15 * n24 * n31 * n43 * n52 - n14 * n25 * n31 * n43 * n52 -
597 n15 * n21 * n34 * n43 * n52 + n11 * n25 * n34 * n43 * n52 +
598 n14 * n21 * n35 * n43 * n52 - n11 * n24 * n35 * n43 * n52 -
599 n15 * n23 * n31 * n44 * n52 + n13 * n25 * n31 * n44 * n52 +
600 n15 * n21 * n33 * n44 * n52 - n11 * n25 * n33 * n44 * n52 -
601 n13 * n21 * n35 * n44 * n52 + n11 * n23 * n35 * n44 * n52 +
602 n14 * n23 * n31 * n45 * n52 - n13 * n24 * n31 * n45 * n52 -
603 n14 * n21 * n33 * n45 * n52 + n11 * n24 * n33 * n45 * n52 +
604 n13 * n21 * n34 * n45 * n52 - n11 * n23 * n34 * n45 * n52 +
605 n15 * n24 * n32 * n41 * n53 - n14 * n25 * n32 * n41 * n53 -
606 n15 * n22 * n34 * n41 * n53 + n12 * n25 * n34 * n41 * n53 +
607 n14 * n22 * n35 * n41 * n53 - n12 * n24 * n35 * n41 * n53 -
608 n15 * n24 * n31 * n42 * n53 + n14 * n25 * n31 * n42 * n53 +
609 n15 * n21 * n34 * n42 * n53 - n11 * n25 * n34 * n42 * n53 -
610 n14 * n21 * n35 * n42 * n53 + n11 * n24 * n35 * n42 * n53 +
611 n15 * n22 * n31 * n44 * n53 - n12 * n25 * n31 * n44 * n53 -
612 n15 * n21 * n32 * n44 * n53 + n11 * n25 * n32 * n44 * n53 +
613 n12 * n21 * n35 * n44 * n53 - n11 * n22 * n35 * n44 * n53 -
614 n14 * n22 * n31 * n45 * n53 + n12 * n24 * n31 * n45 * n53 +
615 n14 * n21 * n32 * n45 * n53 - n11 * n24 * n32 * n45 * n53 -
616 n12 * n21 * n34 * n45 * n53 + n11 * n22 * n34 * n45 * n53 -
617 n15 * n23 * n32 * n41 * n54 + n13 * n25 * n32 * n41 * n54 +
618 n15 * n22 * n33 * n41 * n54 - n12 * n25 * n33 * n41 * n54 -
619 n13 * n22 * n35 * n41 * n54 + n12 * n23 * n35 * n41 * n54 +
620 n15 * n23 * n31 * n42 * n54 - n13 * n25 * n31 * n42 * n54 -
621 n15 * n21 * n33 * n42 * n54 + n11 * n25 * n33 * n42 * n54 +
622 n13 * n21 * n35 * n42 * n54 - n11 * n23 * n35 * n42 * n54 -
623 n15 * n22 * n31 * n43 * n54 + n12 * n25 * n31 * n43 * n54 +
624 n15 * n21 * n32 * n43 * n54 - n11 * n25 * n32 * n43 * n54 -
625 n12 * n21 * n35 * n43 * n54 + n11 * n22 * n35 * n43 * n54 +
626 n13 * n22 * n31 * n45 * n54 - n12 * n23 * n31 * n45 * n54 -
627 n13 * n21 * n32 * n45 * n54 + n11 * n23 * n32 * n45 * n54 +
628 n12 * n21 * n33 * n45 * n54 - n11 * n22 * n33 * n45 * n54 +
629 n14 * n23 * n32 * n41 * n55 - n13 * n24 * n32 * n41 * n55 -
630 n14 * n22 * n33 * n41 * n55 + n12 * n24 * n33 * n41 * n55 +
631 n13 * n22 * n34 * n41 * n55 - n12 * n23 * n34 * n41 * n55 -
632 n14 * n23 * n31 * n42 * n55 + n13 * n24 * n31 * n42 * n55 +
633 n14 * n21 * n33 * n42 * n55 - n11 * n24 * n33 * n42 * n55 -
634 n13 * n21 * n34 * n42 * n55 + n11 * n23 * n34 * n42 * n55 +
635 n14 * n22 * n31 * n43 * n55 - n12 * n24 * n31 * n43 * n55 -
636 n14 * n21 * n32 * n43 * n55 + n11 * n24 * n32 * n43 * n55 +
637 n12 * n21 * n34 * n43 * n55 - n11 * n22 * n34 * n43 * n55 -
638 n13 * n22 * n31 * n44 * n55 + n12 * n23 * n31 * n44 * n55 +
639 n13 * n21 * n32 * n44 * n55 - n11 * n23 * n32 * n44 * n55 -
640 n12 * n21 * n33 * n44 * n55 + n11 * n22 * n33 * n44 * n55 );
641
642 if ( det == complex<double>( 0., 0. ) ) reject = true;
643
644 // The 1st row of the inverse matrix {(I-iKp)^-1}_0j
645 i[0][0] = ( n25 * n34 * n43 * n52 - n24 * n35 * n43 * n52 - n25 * n33 * n44 * n52 +
646 n23 * n35 * n44 * n52 + n24 * n33 * n45 * n52 - n23 * n34 * n45 * n52 -
647 n25 * n34 * n42 * n53 + n24 * n35 * n42 * n53 + n25 * n32 * n44 * n53 -
648 n22 * n35 * n44 * n53 - n24 * n32 * n45 * n53 + n22 * n34 * n45 * n53 +
649 n25 * n33 * n42 * n54 - n23 * n35 * n42 * n54 - n25 * n32 * n43 * n54 +
650 n22 * n35 * n43 * n54 + n23 * n32 * n45 * n54 - n22 * n33 * n45 * n54 -
651 n24 * n33 * n42 * n55 + n23 * n34 * n42 * n55 + n24 * n32 * n43 * n55 -
652 n22 * n34 * n43 * n55 - n23 * n32 * n44 * n55 + n22 * n33 * n44 * n55 ) /
653 det;
654
655 i[0][1] = ( -n15 * n34 * n43 * n52 + n14 * n35 * n43 * n52 + n15 * n33 * n44 * n52 -
656 n13 * n35 * n44 * n52 - n14 * n33 * n45 * n52 + n13 * n34 * n45 * n52 +
657 n15 * n34 * n42 * n53 - n14 * n35 * n42 * n53 - n15 * n32 * n44 * n53 +
658 n12 * n35 * n44 * n53 + n14 * n32 * n45 * n53 - n12 * n34 * n45 * n53 -
659 n15 * n33 * n42 * n54 + n13 * n35 * n42 * n54 + n15 * n32 * n43 * n54 -
660 n12 * n35 * n43 * n54 - n13 * n32 * n45 * n54 + n12 * n33 * n45 * n54 +
661 n14 * n33 * n42 * n55 - n13 * n34 * n42 * n55 - n14 * n32 * n43 * n55 +
662 n12 * n34 * n43 * n55 + n13 * n32 * n44 * n55 - n12 * n33 * n44 * n55 ) /
663 det;
664
665 i[0][2] = ( n15 * n24 * n43 * n52 - n14 * n25 * n43 * n52 - n15 * n23 * n44 * n52 +
666 n13 * n25 * n44 * n52 + n14 * n23 * n45 * n52 - n13 * n24 * n45 * n52 -
667 n15 * n24 * n42 * n53 + n14 * n25 * n42 * n53 + n15 * n22 * n44 * n53 -
668 n12 * n25 * n44 * n53 - n14 * n22 * n45 * n53 + n12 * n24 * n45 * n53 +
669 n15 * n23 * n42 * n54 - n13 * n25 * n42 * n54 - n15 * n22 * n43 * n54 +
670 n12 * n25 * n43 * n54 + n13 * n22 * n45 * n54 - n12 * n23 * n45 * n54 -
671 n14 * n23 * n42 * n55 + n13 * n24 * n42 * n55 + n14 * n22 * n43 * n55 -
672 n12 * n24 * n43 * n55 - n13 * n22 * n44 * n55 + n12 * n23 * n44 * n55 ) /
673 det;
674
675 i[0][3] = ( -n15 * n24 * n33 * n52 + n14 * n25 * n33 * n52 + n15 * n23 * n34 * n52 -
676 n13 * n25 * n34 * n52 - n14 * n23 * n35 * n52 + n13 * n24 * n35 * n52 +
677 n15 * n24 * n32 * n53 - n14 * n25 * n32 * n53 - n15 * n22 * n34 * n53 +
678 n12 * n25 * n34 * n53 + n14 * n22 * n35 * n53 - n12 * n24 * n35 * n53 -
679 n15 * n23 * n32 * n54 + n13 * n25 * n32 * n54 + n15 * n22 * n33 * n54 -
680 n12 * n25 * n33 * n54 - n13 * n22 * n35 * n54 + n12 * n23 * n35 * n54 +
681 n14 * n23 * n32 * n55 - n13 * n24 * n32 * n55 - n14 * n22 * n33 * n55 +
682 n12 * n24 * n33 * n55 + n13 * n22 * n34 * n55 - n12 * n23 * n34 * n55 ) /
683 det;
684
685 i[0][4] = ( n15 * n24 * n33 * n42 - n14 * n25 * n33 * n42 - n15 * n23 * n34 * n42 +
686 n13 * n25 * n34 * n42 + n14 * n23 * n35 * n42 - n13 * n24 * n35 * n42 -
687 n15 * n24 * n32 * n43 + n14 * n25 * n32 * n43 + n15 * n22 * n34 * n43 -
688 n12 * n25 * n34 * n43 - n14 * n22 * n35 * n43 + n12 * n24 * n35 * n43 +
689 n15 * n23 * n32 * n44 - n13 * n25 * n32 * n44 - n15 * n22 * n33 * n44 +
690 n12 * n25 * n33 * n44 + n13 * n22 * n35 * n44 - n12 * n23 * n35 * n44 -
691 n14 * n23 * n32 * n45 + n13 * n24 * n32 * n45 + n14 * n22 * n33 * n45 -
692 n12 * n24 * n33 * n45 - n13 * n22 * n34 * n45 + n12 * n23 * n34 * n45 ) /
693 det;
694
695 //--------------------------------------------------------------------------------------------------------------------------------
696 double s0_prod = -0.07;
697
698 double u1j_re_limit[5][2], u1j_im_limit[5][2];
699 u1j_re_limit[0][0] = 0.;
700 u1j_re_limit[0][1] = 1.;
701 u1j_re_limit[1][0] = -0.29;
702 u1j_re_limit[1][1] = 0.12;
703 u1j_re_limit[2][0] = -0.17;
704 u1j_re_limit[2][1] = 0.065;
705 u1j_re_limit[3][0] = -0.66;
706 u1j_re_limit[3][1] = 0.1;
707 u1j_re_limit[4][0] = -1.36;
708 u1j_re_limit[4][1] = 0.18;
709
710 u1j_im_limit[0][0] = -0.58;
711 u1j_im_limit[0][1] = 0.58;
712 u1j_im_limit[1][0] = 0.00;
713 u1j_im_limit[1][1] = 0.28;
714 u1j_im_limit[2][0] = -0.135;
715 u1j_im_limit[2][1] = 0.10;
716 u1j_im_limit[3][0] = -0.13;
717 u1j_im_limit[3][1] = 0.40;
718 u1j_im_limit[4][0] = -0.36;
719 u1j_im_limit[4][1] = 0.80;
720
721 // for(int kk=0 ; kk<5 ; kk++)
722 //{
723 // cout<<"beta["<<kk<<"]: "<<abs(beta[kk])<<",
724 // \t"<<arg(beta[kk])*EvtConst::radToDegrees<<endl; cout<<"fprod["<<kk<<"]:
725 // "<<abs(fprod[kk])<<", \t"<<arg(fprod[kk])*EvtConst::radToDegrees<<endl;
726 // }
727
728 complex<double> value0( 0., 0. );
729 complex<double> value1( 0., 0. );
730
731 for ( int l = 0; l < 5; l++ )
732 {
733 double u1j_re = real( i[0][l] );
734 double u1j_im = imag( i[0][l] );
735 if ( u1j_re == 0. || u1j_im == 0. ) reject = true;
736 if ( u1j_re < u1j_re_limit[l][0] || u1j_re > u1j_re_limit[l][1] ||
737 u1j_im < u1j_im_limit[l][0] || u1j_im > u1j_im_limit[l][1] )
738 reject = true;
739
740 for ( int pole_index = 0; pole_index < 5; pole_index++ )
741 {
742 complex<double> A = beta[pole_index] * g[pole_index][l];
743 value0 += ( i[0][l] * A ) / ( ma[pole_index] * ma[pole_index] - s );
744 // cout<<"value0["<<l<<"]["<<pole_index<<"] = "<<value0<<endl;
745 }
746 }
747
748 value1 += i[0][0] * fprod[0];
749 value1 += i[0][1] * fprod[1];
750 value1 += i[0][2] * fprod[2];
751 value1 += i[0][3] * fprod[3];
752 value1 += i[0][4] * fprod[4];
753 value1 *= ( 1. - s0_prod ) / ( s - s0_prod );
754 // cout<<"value1 = "<<value1<<endl;
755
756 if ( reject == true ) return complex<double>( 9999., 9999. );
757 else return ( value0 + value1 );
758 // return (value0+value1) ;
759}
760
761complex<double> D0ToKLpipi::amplitude_LASS( vector<double> p_k0l, vector<double> p_pip,
762 vector<double> p_pim, string reso, double A_r,
763 double Phi_r ) {
764
765 double mR = 1.425;
766 double gammaR = 0.27;
767 double mab2 = 0.0;
768 HepLorentzVector _p_k0l( p_k0l[0], p_k0l[1], p_k0l[2], p_k0l[3] );
769 HepLorentzVector _p_pip( p_pip[0], p_pip[1], p_pip[2], p_pip[3] );
770 HepLorentzVector _p_pim( p_pim[0], p_pim[1], p_pim[2], p_pim[3] );
771 if ( reso == "k0lpim" ) mab2 = pow( ( _p_k0l + _p_pim ).m(), 2 );
772 else if ( reso == "k0lpip" ) mab2 = pow( ( _p_k0l + _p_pip ).m(), 2 );
773
774 double s = mab2;
775
776 const double mD0 = 1.86483;
777 const double mKl = 0.49761;
778 const double mPi = 0.13957;
779
780 double _a = 0.113;
781 double _r = -33.8;
782 double _R = 1.0;
783 double _F = 0.96;
784 double _phiR = -1.9146;
785 double _phiF = 0.0017;
786
787 double fR = 1.0; // K*0(1430) has spin zero
788 int power = 1; // Power is 1 for spin zero
789
790 double mAB = sqrt( mab2 ); // (_p4_d1+_p4_d2).mass();
791
792 double mA = mKl; // _p4_d1.mass();
793 double mB = mPi; // _p4_d2.mass();
794 double mC = mPi;
795 double mD = mD0;
796
797 double pAB =
798 sqrt( ( ( ( mAB * mAB - mA * mA - mB * mB ) * ( mAB * mAB - mA * mA - mB * mB ) / 4.0 ) -
799 mA * mA * mB * mB ) /
800 ( mAB * mAB ) );
801 double q = pAB;
802
803 double pR =
804 sqrt( ( ( ( mR * mR - mA * mA - mB * mB ) * ( mR * mR - mA * mA - mB * mB ) / 4.0 ) -
805 mA * mA * mB * mB ) /
806 ( mR * mR ) );
807 double q0 = pR;
808
809 // Running width.
810 double g = gammaR * pow( q / q0, power ) * ( mR / mAB ) * fR * fR;
811
812 complex<double> propagator_relativistic_BreitWigner =
813 1. / ( mR * mR - mAB * mAB - complex<double>( 0., mR * g ) );
814
815 // Non-resonant phase shift
816 double cot_deltaF = 1.0 / ( _a * q ) + 0.5 * _r * q;
817 double qcot_deltaF = 1.0 / _a + 0.5 * _r * q * q;
818
819 // Compute resonant part
820 complex<double> expi2deltaF =
821 complex<double>( qcot_deltaF, q ) / complex<double>( qcot_deltaF, -q );
822
823 complex<double> resonant_term_T =
824 _R * complex<double>( cos( _phiR + 2 * _phiF ), sin( _phiR + 2 * _phiF ) ) *
825 propagator_relativistic_BreitWigner * mR * gammaR * mR / q0 * expi2deltaF;
826
827 // Compute non-resonant part
828 complex<double> non_resonant_term_F = _F * complex<double>( cos( _phiF ), sin( _phiF ) ) *
829 ( cos( _phiF ) + cot_deltaF * sin( _phiF ) ) *
830 sqrt( s ) / complex<double>( qcot_deltaF, -q );
831
832 // Add non-resonant and resonant terms
833 complex<double> LASS_contribution = non_resonant_term_F + resonant_term_T;
834
835 return complex<double>( A_r * cos( Phi_r ), A_r * sin( Phi_r ) ) * LASS_contribution;
836}
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
double imag(const EvtComplex &c)
double meta
double mpi
double mPi
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
const double mD0
Definition MyConst.h:5
***************************************************************************************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
complex< double > Amp_PFT(vector< double > k0l, vector< double > pip, vector< double > pim)
virtual ~D0ToKLpipi()
void init()