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