BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
D0ToKSpipipi0.cxx
Go to the documentation of this file.
1#include "D0ToKSpipipi0.h"
2#include "TLorentzVector.h"
3#include "TMath.h"
4#include <complex>
5#include <map>
6#include <math.h>
7#include <stdlib.h>
8#include <string>
9#include <vector>
10
11#include "CLHEP/Matrix/Matrix.h"
12#include "CLHEP/Matrix/SymMatrix.h"
13#include "CLHEP/Matrix/Vector.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
22#include <cassert>
23#include <iostream>
24
25using namespace std; //::endl;
26
27template <typename T> string toString( const T& t ) {
28 stringstream oss;
29 oss << t;
30 return oss.str();
31}
32
33std::vector<double> operator+( std::vector<double>& lhs, std::vector<double>& rhs ) {
34 assert( int( lhs.size() ) == int( rhs.size() ) );
35 std::vector<double> sum;
36 sum.clear();
37
38 for ( int i = 0; i < int( lhs.size() ); i++ ) { sum.push_back( lhs[i] + rhs[i] ); }
39
40 return sum;
41}
42
44
46
47 // std::cout << "D0ToKSpipipi0 ==> Initialization !" << std::endl;
48
49 _nd = 4;
50 rRes = 3.0;
51 rD = 5.0;
52 math_pi = 3.1415926;
53 mpip = 0.13957;
54 mpim = 0.13957;
55 mpion = 0.13957;
56 mkaon = 0.493677;
57}
58
59void D0ToKSpipipi0::readInputFile() {
60
61 resonance_par["kstarm_892_mass"] = 0.89166;
62 resonance_par["kstarm_892_width"] = 0.0508;
63 resonance_par["kstar0_892_mass"] = 0.89581;
64 resonance_par["kstar0_892_width"] = 0.0474;
65 resonance_par["kstarp_892_mass"] = 0.89166;
66 resonance_par["kstarp_892_width"] = 0.0508;
67 resonance_par["k1m_1270_mass"] = 1.272;
68 resonance_par["k1m_1270_width"] = 0.087;
69 resonance_par["k10_1270_mass"] = 1.272;
70 resonance_par["k10_1270_width"] = 0.087;
71 resonance_par["k1m_1400_mass"] = 1.403;
72 resonance_par["k1m_1400_width"] = 0.174;
73 resonance_par["k10_1400_mass"] = 1.403;
74 resonance_par["k10_1400_width"] = 0.174;
75 resonance_par["kstarm_1410_mass"] = 1.414;
76 resonance_par["kstarm_1410_width"] = 0.232;
77 resonance_par["kstar0_1410_mass"] = 1.414;
78 resonance_par["kstar0_1410_width"] = 0.232;
79 resonance_par["k0starm_1430_mass"] = 1.425;
80 resonance_par["k0starm_1430_width"] = 0.27;
81 resonance_par["k0star0_1430_mass"] = 1.425;
82 resonance_par["k0star0_1430_width"] = 0.27;
83 resonance_par["k2starm_1430_mass"] = 1.4256;
84 resonance_par["k2starm_1430_width"] = 0.0985;
85 resonance_par["k2star0_1430_mass"] = 1.4324;
86 resonance_par["k2star0_1430_width"] = 0.109;
87 resonance_par["km_1460_mass"] = 1.46;
88 resonance_par["km_1460_width"] = 0.26;
89 resonance_par["k0_1460_mass"] = 1.46;
90 resonance_par["k0_1460_width"] = 0.26;
91 resonance_par["k2m_1580_mass"] = 1.58;
92 resonance_par["k2m_1580_width"] = 0.11;
93 resonance_par["k20_1580_mass"] = 1.58;
94 resonance_par["k20_1580_width"] = 0.11;
95 resonance_par["km_1630_mass"] = 1.595;
96 resonance_par["km_1630_width"] = 0.016;
97 resonance_par["k0_1630_mass"] = 1.595;
98 resonance_par["k0_1630_width"] = 0.016;
99 resonance_par["k1m_1650_mass"] = 1.65;
100 resonance_par["k1m_1650_width"] = 0.15;
101 resonance_par["k10_1650_mass"] = 1.65;
102 resonance_par["k10_1650_width"] = 0.15;
103 resonance_par["kstarm_1680_mass"] = 1.717;
104 resonance_par["kstarm_1680_width"] = 0.322;
105 resonance_par["kstar0_1680_mass"] = 1.717;
106 resonance_par["kstar0_1680_width"] = 0.322;
107 resonance_par["k2m_1770_mass"] = 1.773;
108 resonance_par["k2m_1770_width"] = 0.186;
109 resonance_par["k20_1770_mass"] = 1.773;
110 resonance_par["k20_1770_width"] = 0.186;
111 resonance_par["sigma_500_mass"] = 0.9264;
112 resonance_par["sigma_500_width"] = 0.45623;
113 resonance_par["sigma_500_gf"] = 0.0024;
114 resonance_par["f0_980_mass"] = 0.965;
115 resonance_par["f0_980_width"] = 0.055;
116 resonance_par["f0_980_g1"] = 0.165;
117 resonance_par["f0_980_g2"] = 4.21;
118 resonance_par["f0_1370_mass"] = 1.37;
119 resonance_par["f0_1370_width"] = 0.26;
120 resonance_par["f2_1270_mass"] = 1.275;
121 resonance_par["f2_1270_width"] = 0.185;
122 resonance_par["rhop_770_mass"] = 0.77511;
123 resonance_par["rhop_770_width"] = 0.1491;
124 resonance_par["rhom_770_mass"] = 0.77511;
125 resonance_par["rhom_770_width"] = 0.1491;
126 resonance_par["rho0_770_mass"] = 0.77526;
127 resonance_par["rho0_770_width"] = 0.1478;
128 resonance_par["omega_782_mass"] = 0.78265;
129 resonance_par["omega_782_width"] = 0.00849;
130 resonance_par["eta_547_mass"] = 0.547862;
131 resonance_par["eta_547_width"] = 0.00131;
132 resonance_par["phi_1020_mass"] = 1.01946;
133 resonance_par["phi_1020_width"] = 0.00426;
134 resonance_par["a10_1260_mass"] = 1.366;
135 resonance_par["a10_1260_width"] = 0.542;
136 resonance_par["kstarm_phsp_mass"] = 2;
137 resonance_par["kstarm_phsp_width"] = 90;
138 resonance_par["kstar0_phsp_mass"] = 2;
139 resonance_par["kstar0_phsp_width"] = 90;
140 resonance_par["kstarp_phsp_mass"] = 2;
141 resonance_par["kstarp_phsp_width"] = 90;
142 resonance_par["rhop_phsp_mass"] = 2;
143 resonance_par["rhop_phsp_width"] = 90;
144 resonance_par["rho0_phsp_mass"] = 2;
145 resonance_par["rho0_phsp_width"] = 90;
146 resonance_par["rhom_phsp_mass"] = 2;
147 resonance_par["rhom_phsp_width"] = 90;
148 resonance_par["k1m_phsp_mass"] = 2;
149 resonance_par["k1m_phsp_width"] = 90;
150 resonance_par["k10_phsp_mass"] = 2;
151 resonance_par["k10_phsp_width"] = 90;
152 resonance_par["k1p_phsp_mass"] = 2;
153 resonance_par["k1p_phsp_width"] = 90;
154 resonance_par["a10_phsp_mass"] = 2;
155 resonance_par["a10_phsp_width"] = 90;
156
157 readInputCoeff();
158}
159
160void D0ToKSpipipi0::initPar() {
161
162 g.clear();
163 for ( int i = 0; i < 4; i++ )
164 {
165 for ( int j = 0; j < 4; j++ )
166 {
167 if ( i != j ) { g.push_back( 0.0 ); }
168 else if ( i < 3 ) { g.push_back( -1.0 ); }
169 else if ( i == 3 ) { g.push_back( 1.0 ); }
170 }
171 }
172
173 epsilon.clear();
174 for ( int i = 0; i < 4; i++ )
175 {
176 for ( int j = 0; j < 4; j++ )
177 {
178 for ( int k = 0; k < 4; k++ )
179 {
180 for ( int l = 0; l < 4; l++ )
181 {
182 if ( i == j || i == k || i == l || j == k || j == l || k == l )
183 { epsilon.push_back( 0.0 ); }
184 else
185 {
186 if ( i == 0 && j == 1 && k == 2 && l == 3 ) epsilon.push_back( 1.0 );
187 if ( i == 0 && j == 1 && k == 3 && l == 2 ) epsilon.push_back( -1.0 );
188 if ( i == 0 && j == 2 && k == 1 && l == 3 ) epsilon.push_back( -1.0 );
189 if ( i == 0 && j == 2 && k == 3 && l == 1 ) epsilon.push_back( 1.0 );
190 if ( i == 0 && j == 3 && k == 1 && l == 2 ) epsilon.push_back( 1.0 );
191 if ( i == 0 && j == 3 && k == 2 && l == 1 ) epsilon.push_back( -1.0 );
192
193 if ( i == 1 && j == 0 && k == 2 && l == 3 ) epsilon.push_back( -1.0 );
194 if ( i == 1 && j == 0 && k == 3 && l == 2 ) epsilon.push_back( 1.0 );
195 if ( i == 1 && j == 2 && k == 0 && l == 3 ) epsilon.push_back( 1.0 );
196 if ( i == 1 && j == 2 && k == 3 && l == 0 ) epsilon.push_back( -1.0 );
197 if ( i == 1 && j == 3 && k == 0 && l == 2 ) epsilon.push_back( -1.0 );
198 if ( i == 1 && j == 3 && k == 2 && l == 0 ) epsilon.push_back( 1.0 );
199
200 if ( i == 2 && j == 0 && k == 1 && l == 3 ) epsilon.push_back( 1.0 );
201 if ( i == 2 && j == 0 && k == 3 && l == 1 ) epsilon.push_back( -1.0 );
202 if ( i == 2 && j == 1 && k == 0 && l == 3 ) epsilon.push_back( -1.0 );
203 if ( i == 2 && j == 1 && k == 3 && l == 0 ) epsilon.push_back( 1.0 );
204 if ( i == 2 && j == 3 && k == 0 && l == 1 ) epsilon.push_back( 1.0 );
205 if ( i == 2 && j == 3 && k == 1 && l == 0 ) epsilon.push_back( -1.0 );
206
207 if ( i == 3 && j == 0 && k == 1 && l == 2 ) epsilon.push_back( -1.0 );
208 if ( i == 3 && j == 0 && k == 2 && l == 1 ) epsilon.push_back( 1.0 );
209 if ( i == 3 && j == 1 && k == 0 && l == 2 ) epsilon.push_back( 1.0 );
210 if ( i == 3 && j == 1 && k == 2 && l == 0 ) epsilon.push_back( -1.0 );
211 if ( i == 3 && j == 2 && k == 0 && l == 1 ) epsilon.push_back( -1.0 );
212 if ( i == 3 && j == 2 && k == 1 && l == 0 ) epsilon.push_back( 1.0 );
213 }
214 }
215 }
216 }
217 }
218
219 totalAmp = ( 0., 0. );
220
221 readInputFile();
222}
223
224//============================================================//
225// Function //
226//============================================================//
227void D0ToKSpipipi0::addPartialWave( complex<double> amp, double mag, double pha ) {
228
229 complex<double> coeff( mag * cos( pha ), mag * sin( pha ) );
230 totalAmp += ( amp * coeff );
231}
232
233void D0ToKSpipipi0::addPartialWave( complex<double> amp1, complex<double> amp2, double mag,
234 double pha ) {
235
236 complex<double> coeff( mag * cos( pha ), mag * sin( pha ) );
237 totalAmp += ( amp1 * amp2 * coeff );
238}
239
240void D0ToKSpipipi0::addPartialWave( double t1, complex<double> resX, complex<double> resY,
241 double mag, double pha ) {
242 complex<double> coeff( mag * cos( pha ), mag * sin( pha ) );
243 totalAmp += ( t1 * resX * resY * coeff );
244}
245
246// For iso-spin conjugate
247void D0ToKSpipipi0::addPartialWave( double t1, complex<double> resX1, complex<double> resY1,
248 double t2, complex<double> resX2, complex<double> resY2,
249 double fact, double mag, double pha ) {
250
251 complex<double> coeff( mag * cos( pha ), mag * sin( pha ) );
252 totalAmp += ( t1 * resX1 * resY1 + fact * t2 * resX2 * resY2 ) * coeff;
253}
254
255//============================================================//
256// Spin Factor //
257//============================================================//
258// void D0ToKSpipipi0::createSpinfactor(EvtVector4R Ks0, EvtVector4R Pip, EvtVector4R Pim,
259// EvtVector4R Pi0) {
260void D0ToKSpipipi0::createSpinfactor( vector<double> ks0, vector<double> pip,
261 vector<double> pim, vector<double> pi0 ) {
262
263 // double m2_ks0 = Ks0.mass2();
264 // double m2_pip = Pip.mass2();
265 // double m2_pim = Pim.mass2();
266 // double m2_pi0 = Pi0.mass2();
267 // double m2_ks0pip = (Ks0 + Pip).mass2();
268 // double m2_ks0pim = (Ks0 + Pim).mass2();
269 // double m2_ks0pi0 = (Ks0 + Pi0).mass2();
270 // double m2_pippim = (Pip + Pim).mass2();
271 // double m2_pippi0 = (Pip + Pi0).mass2();
272 // double m2_pimpi0 = (Pim + Pi0).mass2();
273 // double m2_ks0pippim = (Ks0 + Pip + Pim).mass2();
274 // double m2_ks0pippi0 = (Ks0 + Pip + Pi0).mass2();
275 // double m2_ks0pimpi0 = (Ks0 + Pim + Pi0).mass2();
276 // double m2_pippimpi0 = (Pip + Pim + Pi0).mass2();
277
278 // std::vector<double> ks0; ks0.clear();
279 // ks0.push_back(Ks0.get(1)); ks0.push_back(Ks0.get(2)); ks0.push_back(Ks0.get(3));
280 // ks0.push_back(Ks0.get(0)); std::vector<double> pip; pip.clear();
281 // pip.push_back(Pip.get(1)); pip.push_back(Pip.get(2)); pip.push_back(Pip.get(3));
282 // pip.push_back(Pip.get(0)); std::vector<double> pim; pim.clear();
283 // pim.push_back(Pim.get(1)); pim.push_back(Pim.get(2)); pim.push_back(Pim.get(3));
284 // pim.push_back(Pim.get(0)); std::vector<double> pi0; pi0.clear();
285 // pi0.push_back(Pi0.get(1)); pi0.push_back(Pi0.get(2)); pi0.push_back(Pi0.get(3));
286 // pi0.push_back(Pi0.get(0));
287
288 //------------------------------------------------------------//
289 // Cascade Decay //
290 //------------------------------------------------------------//
291 //---- D2PP_P2SP
292 spinfactor["D2PP_P2SP_pippimpi0_pimpi0_0"] = D2PP_P2SP();
293 spinfactor["D2PP_P2SP_pippimpi0_pippi0_0"] = D2PP_P2SP();
294 spinfactor["D2PP_P2SP_pippimpi0_pippim_0"] = D2PP_P2SP(); //
295 spinfactor["D2PP_P2SP_ks0pimpi0_ks0pim_0"] = D2PP_P2SP();
296 spinfactor["D2PP_P2SP_ks0pimpi0_ks0pi0_0"] = D2PP_P2SP();
297 spinfactor["D2PP_P2SP_ks0pimpi0_pimpi0_0"] = D2PP_P2SP(); //
298 spinfactor["D2PP_P2SP_ks0pippi0_ks0pip_0"] = D2PP_P2SP();
299 spinfactor["D2PP_P2SP_ks0pippi0_ks0pi0_0"] = D2PP_P2SP();
300 spinfactor["D2PP_P2SP_ks0pippi0_pippi0_0"] = D2PP_P2SP(); //
301 spinfactor["D2PP_P2SP_ks0pippim_ks0pip_0"] = D2PP_P2SP();
302 spinfactor["D2PP_P2SP_ks0pippim_ks0pim_0"] = D2PP_P2SP();
303 spinfactor["D2PP_P2SP_ks0pippim_pippim_0"] = D2PP_P2SP(); //
304
305 //---- D2PP_P2VP
306 spinfactor["D2PP_P2VP_pippimpi0_pimpi0_1"] =
307 D2PP_P2VP( pi0, pim, pip, ks0, 1 ); // sequence of pi0 pi-
308 spinfactor["D2PP_P2VP_pippimpi0_pippi0_1"] = D2PP_P2VP( pip, pi0, pim, ks0, 1 );
309 spinfactor["D2PP_P2VP_pippimpi0_pippim_1"] = D2PP_P2VP( pip, pim, pi0, ks0, 1 ); //
310 spinfactor["D2PP_P2VP_ks0pimpi0_ks0pim_1"] = D2PP_P2VP( ks0, pim, pi0, pip, 1 );
311 spinfactor["D2PP_P2VP_ks0pimpi0_ks0pi0_1"] = D2PP_P2VP( ks0, pi0, pim, pip, 1 );
312 spinfactor["D2PP_P2VP_ks0pimpi0_pimpi0_1"] = D2PP_P2VP( pim, pi0, ks0, pip, 1 ); //
313 spinfactor["D2PP_P2VP_ks0pippi0_ks0pip_1"] = D2PP_P2VP( ks0, pip, pi0, pim, 1 );
314 spinfactor["D2PP_P2VP_ks0pippi0_ks0pi0_1"] = D2PP_P2VP( ks0, pi0, pip, pim, 1 );
315 spinfactor["D2PP_P2VP_ks0pippi0_pippi0_1"] = D2PP_P2VP( pip, pi0, ks0, pim, 1 ); //
316 spinfactor["D2PP_P2VP_ks0pippim_ks0pip_1"] = D2PP_P2VP( ks0, pip, pim, pi0, 1 );
317 spinfactor["D2PP_P2VP_ks0pippim_ks0pim_1"] = D2PP_P2VP( ks0, pim, pip, pi0, 1 );
318 spinfactor["D2PP_P2VP_ks0pippim_pippim_1"] = D2PP_P2VP( pip, pim, ks0, pi0, 1 ); //
319
320 //---- D2VP_V2VP
321 spinfactor["D2VP_V2VP_pippimpi0_pimpi0_1"] =
322 D2VP_V2VP( pi0, pim, pip, ks0, 1 ); // sequence of pi0 pi-
323 spinfactor["D2VP_V2VP_pippimpi0_pippi0_1"] = D2VP_V2VP( pip, pi0, pim, ks0, 1 );
324 spinfactor["D2VP_V2VP_pippimpi0_pippim_1"] = D2VP_V2VP( pip, pim, pi0, ks0, 1 ); //
325 spinfactor["D2VP_V2VP_ks0pimpi0_ks0pim_1"] = D2VP_V2VP( ks0, pim, pi0, pip, 1 );
326 spinfactor["D2VP_V2VP_ks0pimpi0_ks0pi0_1"] = D2VP_V2VP( ks0, pi0, pim, pip, 1 );
327 spinfactor["D2VP_V2VP_ks0pimpi0_pimpi0_1"] = D2VP_V2VP( pim, pi0, ks0, pip, 1 ); //
328 spinfactor["D2VP_V2VP_ks0pippi0_ks0pip_1"] = D2VP_V2VP( ks0, pip, pi0, pim, 1 );
329 spinfactor["D2VP_V2VP_ks0pippi0_ks0pi0_1"] = D2VP_V2VP( ks0, pi0, pip, pim, 1 );
330 spinfactor["D2VP_V2VP_ks0pippi0_pippi0_1"] = D2VP_V2VP( pip, pi0, ks0, pim, 1 ); //
331 spinfactor["D2VP_V2VP_ks0pippim_ks0pip_1"] = D2VP_V2VP( ks0, pip, pim, pi0, 1 );
332 spinfactor["D2VP_V2VP_ks0pippim_ks0pim_1"] = D2VP_V2VP( ks0, pim, pip, pi0, 1 );
333 spinfactor["D2VP_V2VP_ks0pippim_pippim_1"] = D2VP_V2VP( pip, pim, ks0, pi0, 1 ); //
334
335 //---- D2AP_A2SP
336 spinfactor["D2AP_A2SP_pippimpi0_pimpi0_1"] =
337 D2AP_A2SP( pi0, pim, pip, ks0, 1 ); // sequence of pi0 pi-
338 spinfactor["D2AP_A2SP_pippimpi0_pippi0_1"] = D2AP_A2SP( pip, pi0, pim, ks0, 1 );
339 spinfactor["D2AP_A2SP_pippimpi0_pippim_1"] = D2AP_A2SP( pip, pim, pi0, ks0, 1 ); //
340 spinfactor["D2AP_A2SP_ks0pimpi0_ks0pim_1"] = D2AP_A2SP( ks0, pim, pi0, pip, 1 );
341 spinfactor["D2AP_A2SP_ks0pimpi0_ks0pi0_1"] = D2AP_A2SP( ks0, pi0, pim, pip, 1 );
342 spinfactor["D2AP_A2SP_ks0pimpi0_pimpi0_1"] = D2AP_A2SP( pim, pi0, ks0, pip, 1 ); //
343 spinfactor["D2AP_A2SP_ks0pippi0_ks0pip_1"] = D2AP_A2SP( ks0, pip, pi0, pim, 1 );
344 spinfactor["D2AP_A2SP_ks0pippi0_ks0pi0_1"] = D2AP_A2SP( ks0, pi0, pip, pim, 1 );
345 spinfactor["D2AP_A2SP_ks0pippi0_pippi0_1"] = D2AP_A2SP( pip, pi0, ks0, pim, 1 ); //
346 spinfactor["D2AP_A2SP_ks0pippim_ks0pip_1"] = D2AP_A2SP( ks0, pip, pim, pi0, 1 );
347 spinfactor["D2AP_A2SP_ks0pippim_ks0pim_1"] = D2AP_A2SP( ks0, pim, pip, pi0, 1 );
348 spinfactor["D2AP_A2SP_ks0pippim_pippim_1"] = D2AP_A2SP( pip, pim, ks0, pi0, 1 ); //
349
350 //---- D2AP_A2VP, [S,D]
351 spinfactor["D2AP_A2VP_pippimpi0_pimpi0_0"] =
352 D2AP_A2VP( pi0, pim, pip, ks0, 0 ); // sequence of pi0 pi-
353 spinfactor["D2AP_A2VP_pippimpi0_pippi0_0"] = D2AP_A2VP( pip, pi0, pim, ks0, 0 );
354 spinfactor["D2AP_A2VP_pippimpi0_pippim_0"] = D2AP_A2VP( pip, pim, pi0, ks0, 0 ); //
355 spinfactor["D2AP_A2VP_ks0pimpi0_ks0pim_0"] = D2AP_A2VP( ks0, pim, pi0, pip, 0 );
356 spinfactor["D2AP_A2VP_ks0pimpi0_ks0pi0_0"] = D2AP_A2VP( ks0, pi0, pim, pip, 0 );
357 spinfactor["D2AP_A2VP_ks0pimpi0_pimpi0_0"] = D2AP_A2VP( pim, pi0, ks0, pip, 0 ); //
358 spinfactor["D2AP_A2VP_ks0pippi0_ks0pip_0"] = D2AP_A2VP( ks0, pip, pi0, pim, 0 );
359 spinfactor["D2AP_A2VP_ks0pippi0_ks0pi0_0"] = D2AP_A2VP( ks0, pi0, pip, pim, 0 );
360 spinfactor["D2AP_A2VP_ks0pippi0_pippi0_0"] = D2AP_A2VP( pip, pi0, ks0, pim, 0 ); //
361 spinfactor["D2AP_A2VP_ks0pippim_ks0pip_0"] = D2AP_A2VP( ks0, pip, pim, pi0, 0 );
362 spinfactor["D2AP_A2VP_ks0pippim_ks0pim_0"] = D2AP_A2VP( ks0, pim, pip, pi0, 0 );
363 spinfactor["D2AP_A2VP_ks0pippim_pippim_0"] = D2AP_A2VP( pip, pim, ks0, pi0, 0 ); //
364 spinfactor["D2AP_A2VP_pippimpi0_pimpi0_2"] =
365 D2AP_A2VP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
366 spinfactor["D2AP_A2VP_pippimpi0_pippi0_2"] = D2AP_A2VP( pip, pi0, pim, ks0, 2 );
367 spinfactor["D2AP_A2VP_pippimpi0_pippim_2"] = D2AP_A2VP( pip, pim, pi0, ks0, 2 ); //
368 spinfactor["D2AP_A2VP_ks0pimpi0_ks0pim_2"] = D2AP_A2VP( ks0, pim, pi0, pip, 2 );
369 spinfactor["D2AP_A2VP_ks0pimpi0_ks0pi0_2"] = D2AP_A2VP( ks0, pi0, pim, pip, 2 );
370 spinfactor["D2AP_A2VP_ks0pimpi0_pimpi0_2"] = D2AP_A2VP( pim, pi0, ks0, pip, 2 ); //
371 spinfactor["D2AP_A2VP_ks0pippi0_ks0pip_2"] = D2AP_A2VP( ks0, pip, pi0, pim, 2 );
372 spinfactor["D2AP_A2VP_ks0pippi0_ks0pi0_2"] = D2AP_A2VP( ks0, pi0, pip, pim, 2 );
373 spinfactor["D2AP_A2VP_ks0pippi0_pippi0_2"] = D2AP_A2VP( pip, pi0, ks0, pim, 2 ); //
374 spinfactor["D2AP_A2VP_ks0pippim_ks0pip_2"] = D2AP_A2VP( ks0, pip, pim, pi0, 2 );
375 spinfactor["D2AP_A2VP_ks0pippim_ks0pim_2"] = D2AP_A2VP( ks0, pim, pip, pi0, 2 );
376 spinfactor["D2AP_A2VP_ks0pippim_pippim_2"] = D2AP_A2VP( pip, pim, ks0, pi0, 2 ); //
377
378 //---- D2AP_A2TP
379 spinfactor["D2AP_A2TP_pippimpi0_pimpi0_1"] =
380 D2AP_A2TP( pi0, pim, pip, ks0, 1 ); // sequence of pi0 pi-
381 spinfactor["D2AP_A2TP_pippimpi0_pippi0_1"] = D2AP_A2TP( pip, pi0, pim, ks0, 1 );
382 spinfactor["D2AP_A2TP_pippimpi0_pippim_1"] = D2AP_A2TP( pip, pim, pi0, ks0, 1 ); //
383 spinfactor["D2AP_A2TP_ks0pimpi0_ks0pim_1"] = D2AP_A2TP( ks0, pim, pi0, pip, 1 );
384 spinfactor["D2AP_A2TP_ks0pimpi0_ks0pi0_1"] = D2AP_A2TP( ks0, pi0, pim, pip, 1 );
385 spinfactor["D2AP_A2TP_ks0pimpi0_pimpi0_1"] = D2AP_A2TP( pim, pi0, ks0, pip, 1 ); //
386 spinfactor["D2AP_A2TP_ks0pippi0_ks0pip_1"] = D2AP_A2TP( ks0, pip, pi0, pim, 1 );
387 spinfactor["D2AP_A2TP_ks0pippi0_ks0pi0_1"] = D2AP_A2TP( ks0, pi0, pip, pim, 1 );
388 spinfactor["D2AP_A2TP_ks0pippi0_pippi0_1"] = D2AP_A2TP( pip, pi0, ks0, pim, 1 ); //
389 spinfactor["D2AP_A2TP_ks0pippim_ks0pip_1"] = D2AP_A2TP( ks0, pip, pim, pi0, 1 );
390 spinfactor["D2AP_A2TP_ks0pippim_ks0pim_1"] = D2AP_A2TP( ks0, pim, pip, pi0, 1 );
391 spinfactor["D2AP_A2TP_ks0pippim_pippim_1"] = D2AP_A2TP( pip, pim, ks0, pi0, 1 ); //
392
393 //---- D2TP_T2VP
394 spinfactor["D2TP_T2VP_pippimpi0_pimpi0_2"] =
395 D2TP_T2VP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
396 spinfactor["D2TP_T2VP_pippimpi0_pippi0_2"] = D2TP_T2VP( pip, pi0, pim, ks0, 2 );
397 spinfactor["D2TP_T2VP_pippimpi0_pippim_2"] = D2TP_T2VP( pip, pim, pi0, ks0, 2 ); //
398 spinfactor["D2TP_T2VP_ks0pimpi0_ks0pim_2"] = D2TP_T2VP( ks0, pim, pi0, pip, 2 );
399 spinfactor["D2TP_T2VP_ks0pimpi0_ks0pi0_2"] = D2TP_T2VP( ks0, pi0, pim, pip, 2 );
400 spinfactor["D2TP_T2VP_ks0pimpi0_pimpi0_2"] = D2TP_T2VP( pim, pi0, ks0, pip, 2 ); //
401 spinfactor["D2TP_T2VP_ks0pippi0_ks0pip_2"] = D2TP_T2VP( ks0, pip, pi0, pim, 2 );
402 spinfactor["D2TP_T2VP_ks0pippi0_ks0pi0_2"] = D2TP_T2VP( ks0, pi0, pip, pim, 2 );
403 spinfactor["D2TP_T2VP_ks0pippi0_pippi0_2"] = D2TP_T2VP( pip, pi0, ks0, pim, 2 ); //
404 spinfactor["D2TP_T2VP_ks0pippim_ks0pip_2"] = D2TP_T2VP( ks0, pip, pim, pi0, 2 );
405 spinfactor["D2TP_T2VP_ks0pippim_ks0pim_2"] = D2TP_T2VP( ks0, pim, pip, pi0, 2 );
406 spinfactor["D2TP_T2VP_ks0pippim_pippim_2"] = D2TP_T2VP( pip, pim, ks0, pi0, 2 ); //
407
408 //---- D2TP_T2TP
409 spinfactor["D2TP_T2TP_pippimpi0_pimpi0_2"] =
410 D2TP_T2TP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
411 spinfactor["D2TP_T2TP_pippimpi0_pippi0_2"] = D2TP_T2TP( pip, pi0, pim, ks0, 2 );
412 spinfactor["D2TP_T2TP_pippimpi0_pippim_2"] = D2TP_T2TP( pip, pim, pi0, ks0, 2 ); //
413 spinfactor["D2TP_T2TP_ks0pimpi0_ks0pim_2"] = D2TP_T2TP( ks0, pim, pi0, pip, 2 );
414 spinfactor["D2TP_T2TP_ks0pimpi0_ks0pi0_2"] = D2TP_T2TP( ks0, pi0, pim, pip, 2 );
415 spinfactor["D2TP_T2TP_ks0pimpi0_pimpi0_2"] = D2TP_T2TP( pim, pi0, ks0, pip, 2 ); //
416 spinfactor["D2TP_T2TP_ks0pippi0_ks0pip_2"] = D2TP_T2TP( ks0, pip, pi0, pim, 2 );
417 spinfactor["D2TP_T2TP_ks0pippi0_ks0pi0_2"] = D2TP_T2TP( ks0, pi0, pip, pim, 2 );
418 spinfactor["D2TP_T2TP_ks0pippi0_pippi0_2"] = D2TP_T2TP( pip, pi0, ks0, pim, 2 ); //
419 spinfactor["D2TP_T2TP_ks0pippim_ks0pip_2"] = D2TP_T2TP( ks0, pip, pim, pi0, 2 );
420 spinfactor["D2TP_T2TP_ks0pippim_ks0pim_2"] = D2TP_T2TP( ks0, pim, pip, pi0, 2 );
421 spinfactor["D2TP_T2TP_ks0pippim_pippim_2"] = D2TP_T2TP( pip, pim, ks0, pi0, 2 ); //
422
423 //---- D2PTP_PT2SP
424 spinfactor["D2PTP_PT2SP_pippimpi0_pimpi0_2"] =
425 D2PTP_PT2SP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
426 spinfactor["D2PTP_PT2SP_pippimpi0_pippi0_2"] = D2PTP_PT2SP( pip, pi0, pim, ks0, 2 );
427 spinfactor["D2PTP_PT2SP_pippimpi0_pippim_2"] = D2PTP_PT2SP( pip, pim, pi0, ks0, 2 ); //
428 spinfactor["D2PTP_PT2SP_ks0pimpi0_ks0pim_2"] = D2PTP_PT2SP( ks0, pim, pi0, pip, 2 );
429 spinfactor["D2PTP_PT2SP_ks0pimpi0_ks0pi0_2"] = D2PTP_PT2SP( ks0, pi0, pim, pip, 2 );
430 spinfactor["D2PTP_PT2SP_ks0pimpi0_pimpi0_2"] = D2PTP_PT2SP( pim, pi0, ks0, pip, 2 ); //
431 spinfactor["D2PTP_PT2SP_ks0pippi0_ks0pip_2"] = D2PTP_PT2SP( ks0, pip, pi0, pim, 2 );
432 spinfactor["D2PTP_PT2SP_ks0pippi0_ks0pi0_2"] = D2PTP_PT2SP( ks0, pi0, pip, pim, 2 );
433 spinfactor["D2PTP_PT2SP_ks0pippi0_pippi0_2"] = D2PTP_PT2SP( pip, pi0, ks0, pim, 2 ); //
434 spinfactor["D2PTP_PT2SP_ks0pippim_ks0pip_2"] = D2PTP_PT2SP( ks0, pip, pim, pi0, 2 );
435 spinfactor["D2PTP_PT2SP_ks0pippim_ks0pim_2"] = D2PTP_PT2SP( ks0, pim, pip, pi0, 2 );
436 spinfactor["D2PTP_PT2SP_ks0pippim_pippim_2"] = D2PTP_PT2SP( pip, pim, ks0, pi0, 2 ); //
437
438 //---- D2PTP_PT2VP
439 spinfactor["D2PTP_PT2VP_pippimpi0_pimpi0_2"] =
440 D2PTP_PT2VP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
441 spinfactor["D2PTP_PT2VP_pippimpi0_pippi0_2"] = D2PTP_PT2VP( pip, pi0, pim, ks0, 2 );
442 spinfactor["D2PTP_PT2VP_pippimpi0_pippim_2"] = D2PTP_PT2VP( pip, pim, pi0, ks0, 2 ); //
443 spinfactor["D2PTP_PT2VP_ks0pimpi0_ks0pim_2"] = D2PTP_PT2VP( ks0, pim, pi0, pip, 2 );
444 spinfactor["D2PTP_PT2VP_ks0pimpi0_ks0pi0_2"] = D2PTP_PT2VP( ks0, pi0, pim, pip, 2 );
445 spinfactor["D2PTP_PT2VP_ks0pimpi0_pimpi0_2"] = D2PTP_PT2VP( pim, pi0, ks0, pip, 2 ); //
446 spinfactor["D2PTP_PT2VP_ks0pippi0_ks0pip_2"] = D2PTP_PT2VP( ks0, pip, pi0, pim, 2 );
447 spinfactor["D2PTP_PT2VP_ks0pippi0_ks0pi0_2"] = D2PTP_PT2VP( ks0, pi0, pip, pim, 2 );
448 spinfactor["D2PTP_PT2VP_ks0pippi0_pippi0_2"] = D2PTP_PT2VP( pip, pi0, ks0, pim, 2 ); //
449 spinfactor["D2PTP_PT2VP_ks0pippim_ks0pip_2"] = D2PTP_PT2VP( ks0, pip, pim, pi0, 2 );
450 spinfactor["D2PTP_PT2VP_ks0pippim_ks0pim_2"] = D2PTP_PT2VP( ks0, pim, pip, pi0, 2 );
451 spinfactor["D2PTP_PT2VP_ks0pippim_pippim_2"] = D2PTP_PT2VP( pip, pim, ks0, pi0, 2 ); //
452
453 //---- D2PTP_PT2TP
454 spinfactor["D2PTP_PT2TP_pippimpi0_pimpi0_2"] =
455 D2PTP_PT2TP( pi0, pim, pip, ks0, 2 ); // sequence of pi0 pi-
456 spinfactor["D2PTP_PT2TP_pippimpi0_pippi0_2"] = D2PTP_PT2TP( pip, pi0, pim, ks0, 2 );
457 spinfactor["D2PTP_PT2TP_pippimpi0_pippim_2"] = D2PTP_PT2TP( pip, pim, pi0, ks0, 2 ); //
458 spinfactor["D2PTP_PT2TP_ks0pimpi0_ks0pim_2"] = D2PTP_PT2TP( ks0, pim, pi0, pip, 2 );
459 spinfactor["D2PTP_PT2TP_ks0pimpi0_ks0pi0_2"] = D2PTP_PT2TP( ks0, pi0, pim, pip, 2 );
460 spinfactor["D2PTP_PT2TP_ks0pimpi0_pimpi0_2"] = D2PTP_PT2TP( pim, pi0, ks0, pip, 2 ); //
461 spinfactor["D2PTP_PT2TP_ks0pippi0_ks0pip_2"] = D2PTP_PT2TP( ks0, pip, pi0, pim, 2 );
462 spinfactor["D2PTP_PT2TP_ks0pippi0_ks0pi0_2"] = D2PTP_PT2TP( ks0, pi0, pip, pim, 2 );
463 spinfactor["D2PTP_PT2TP_ks0pippi0_pippi0_2"] = D2PTP_PT2TP( pip, pi0, ks0, pim, 2 ); //
464 spinfactor["D2PTP_PT2TP_ks0pippim_ks0pip_2"] = D2PTP_PT2TP( ks0, pip, pim, pi0, 2 );
465 spinfactor["D2PTP_PT2TP_ks0pippim_ks0pim_2"] = D2PTP_PT2TP( ks0, pim, pip, pi0, 2 );
466 spinfactor["D2PTP_PT2TP_ks0pippim_pippim_2"] = D2PTP_PT2TP( pip, pim, ks0, pi0, 2 ); //
467
468 //------------------------------------------------------------//
469 // Quasi Two-Body Decay //
470 //------------------------------------------------------------//
471 //---- D2SS
472 spinfactor["D2SS_ks0pim_pippi0_0"] = D2SS();
473 spinfactor["D2SS_ks0pi0_pippim_0"] = D2SS();
474 spinfactor["D2SS_ks0pip_pimpi0_0"] = D2SS();
475
476 //---- D2VS
477 spinfactor["D2VS_ks0pim_pippi0_1"] = D2VS( ks0, pim, pip, pi0, 1 );
478 spinfactor["D2VS_ks0pi0_pippim_1"] = D2VS( ks0, pi0, pip, pim, 1 );
479 spinfactor["D2VS_ks0pip_pimpi0_1"] = D2VS( ks0, pip, pim, pi0, 1 );
480 spinfactor["D2VS_pimpi0_ks0pip_1"] = D2VS( pim, pi0, ks0, pip, 1 );
481 spinfactor["D2VS_pippim_ks0pi0_1"] = D2VS( pip, pim, ks0, pi0, 1 );
482 spinfactor["D2VS_pippi0_ks0pim_1"] = D2VS( pip, pi0, ks0, pim, 1 );
483
484 //---- D2VV, [S,P,D]
485 spinfactor["D2VV_ks0pim_pippi0_0"] = D2VV( ks0, pim, pip, pi0, 0 );
486 spinfactor["D2VV_ks0pi0_pippim_0"] = D2VV( ks0, pi0, pip, pim, 0 );
487 spinfactor["D2VV_ks0pip_pimpi0_0"] = D2VV( ks0, pip, pim, pi0, 0 );
488 spinfactor["D2VV_ks0pim_pippi0_1"] = D2VV( ks0, pim, pip, pi0, 1 );
489 spinfactor["D2VV_ks0pi0_pippim_1"] = D2VV( ks0, pi0, pip, pim, 1 );
490 spinfactor["D2VV_ks0pip_pimpi0_1"] = D2VV( ks0, pip, pim, pi0, 1 );
491 spinfactor["D2VV_ks0pim_pippi0_2"] = D2VV( ks0, pim, pip, pi0, 2 );
492 spinfactor["D2VV_ks0pi0_pippim_2"] = D2VV( ks0, pi0, pip, pim, 2 );
493 spinfactor["D2VV_ks0pip_pimpi0_2"] = D2VV( ks0, pip, pim, pi0, 2 );
494
495 //---- D2TS
496 spinfactor["D2TS_ks0pim_pippi0_2"] = D2TS( ks0, pim, pip, pi0, 2 );
497 spinfactor["D2TS_ks0pi0_pippim_2"] = D2TS( ks0, pi0, pip, pim, 2 );
498 spinfactor["D2TS_ks0pip_pimpi0_2"] = D2TS( ks0, pip, pim, pi0, 2 );
499 spinfactor["D2TS_pimpi0_ks0pip_2"] = D2TS( pim, pi0, ks0, pip, 2 );
500 spinfactor["D2TS_pippim_ks0pi0_2"] = D2TS( pip, pim, ks0, pi0, 2 );
501 spinfactor["D2TS_pippi0_ks0pim_2"] = D2TS( pip, pi0, ks0, pim, 2 );
502
503 //---- D2TV, [P, D]
504 spinfactor["D2TV_ks0pim_pippi0_1"] = D2TV( ks0, pim, pip, pi0, 1 );
505 spinfactor["D2TV_ks0pi0_pippim_1"] = D2TV( ks0, pi0, pip, pim, 1 );
506 spinfactor["D2TV_ks0pip_pimpi0_1"] = D2TV( ks0, pip, pim, pi0, 1 );
507 spinfactor["D2TV_pimpi0_ks0pip_1"] = D2TV( pim, pi0, ks0, pip, 1 );
508 spinfactor["D2TV_pippim_ks0pi0_1"] = D2TV( pip, pim, ks0, pi0, 1 );
509 spinfactor["D2TV_pippi0_ks0pim_1"] = D2TV( pip, pi0, ks0, pim, 1 );
510 spinfactor["D2TV_ks0pim_pippi0_2"] = D2TV( ks0, pim, pip, pi0, 2 );
511 spinfactor["D2TV_ks0pi0_pippim_2"] = D2TV( ks0, pi0, pip, pim, 2 );
512 spinfactor["D2TV_ks0pip_pimpi0_2"] = D2TV( ks0, pip, pim, pi0, 2 );
513 spinfactor["D2TV_pimpi0_ks0pip_2"] = D2TV( pim, pi0, ks0, pip, 2 );
514 spinfactor["D2TV_pippim_ks0pi0_2"] = D2TV( pip, pim, ks0, pi0, 2 );
515 spinfactor["D2TV_pippi0_ks0pim_2"] = D2TV( pip, pi0, ks0, pim, 2 );
516
517 //---- D2TT, [S,P,D]
518 spinfactor["D2TT_ks0pim_pippi0_0"] = D2TT( ks0, pim, pip, pi0, 0 );
519 spinfactor["D2TT_ks0pi0_pippim_0"] = D2TT( ks0, pi0, pip, pim, 0 );
520 spinfactor["D2TT_ks0pip_pimpi0_0"] = D2TT( ks0, pip, pim, pi0, 0 );
521 spinfactor["D2TT_ks0pim_pippi0_1"] = D2TT( ks0, pim, pip, pi0, 1 );
522 spinfactor["D2TT_ks0pi0_pippim_1"] = D2TT( ks0, pi0, pip, pim, 1 );
523 spinfactor["D2TT_ks0pip_pimpi0_1"] = D2TT( ks0, pip, pim, pi0, 1 );
524 spinfactor["D2TT_ks0pim_pippi0_2"] = D2TT( ks0, pim, pip, pi0, 2 );
525 spinfactor["D2TT_ks0pi0_pippim_2"] = D2TT( ks0, pi0, pip, pim, 2 );
526 spinfactor["D2TT_ks0pip_pimpi0_2"] = D2TT( ks0, pip, pim, pi0, 2 );
527}
528
529//============================================================//
530// Propagator //
531//============================================================//
532// void D0ToKSpipipi0::createPropagator(EvtVector4R Ks0, EvtVector4R Pip, EvtVector4R Pim,
533// EvtVector4R Pi0) {
534void D0ToKSpipipi0::createPropagator( vector<double> Ks0, vector<double> Pip,
535 vector<double> Pim, vector<double> Pi0 ) {
536
537 HepLorentzVector ks0;
538 ks0.setX( Ks0[0] );
539 ks0.setY( Ks0[1] );
540 ks0.setZ( Ks0[2] );
541 ks0.setT( Ks0[3] );
542 HepLorentzVector pip;
543 pip.setX( Pip[0] );
544 pip.setY( Pip[1] );
545 pip.setZ( Pip[2] );
546 pip.setT( Pip[3] );
547 HepLorentzVector pim;
548 pim.setX( Pim[0] );
549 pim.setY( Pim[1] );
550 pim.setZ( Pim[2] );
551 pim.setT( Pim[3] );
552 HepLorentzVector pi0;
553 pi0.setX( Pi0[0] );
554 pi0.setY( Pi0[1] );
555 pi0.setZ( Pi0[2] );
556 pi0.setT( Pi0[3] );
557
558 // double m2_ks0 = Ks0.mass2();
559 // double m2_pip = Pip.mass2();
560 // double m2_pim = Pim.mass2();
561 // double m2_pi0 = Pi0.mass2();
562 // double m2_ks0pip = (Ks0 + Pip).mass2();
563 // double m2_ks0pim = (Ks0 + Pim).mass2();
564 // double m2_ks0pi0 = (Ks0 + Pi0).mass2();
565 // double m2_pippim = (Pip + Pim).mass2();
566 // double m2_pippi0 = (Pip + Pi0).mass2();
567 // double m2_pimpi0 = (Pim + Pi0).mass2();
568 // double m2_ks0pippim = (Ks0 + Pip + Pim).mass2();
569 // double m2_ks0pippi0 = (Ks0 + Pip + Pi0).mass2();
570 // double m2_ks0pimpi0 = (Ks0 + Pim + Pi0).mass2();
571 // double m2_pippimpi0 = (Pip + Pim + Pi0).mass2();
572 double m2_ks0 = ks0.invariantMass2();
573 double m2_pip = pip.invariantMass2();
574 double m2_pim = pim.invariantMass2();
575 double m2_pi0 = pi0.invariantMass2();
576 double m2_ks0pip = ( ks0 + pip ).invariantMass2();
577 double m2_ks0pim = ( ks0 + pim ).invariantMass2();
578 double m2_ks0pi0 = ( ks0 + pi0 ).invariantMass2();
579 double m2_pippim = ( pip + pim ).invariantMass2();
580 double m2_pippi0 = ( pip + pi0 ).invariantMass2();
581 double m2_pimpi0 = ( pim + pi0 ).invariantMass2();
582 double m2_ks0pippim = ( ks0 + pip + pim ).invariantMass2();
583 double m2_ks0pippi0 = ( ks0 + pip + pi0 ).invariantMass2();
584 double m2_ks0pimpi0 = ( ks0 + pim + pi0 ).invariantMass2();
585 double m2_pippimpi0 = ( pip + pim + pi0 ).invariantMass2();
586
587 // std::vector<double> ks0; ks0.clear();
588 // ks0.push_back(Ks0.get(1)); ks0.push_back(Ks0.get(2)); ks0.push_back(Ks0.get(3));
589 // ks0.push_back(Ks0.get(0)); std::vector<double> pip; pip.clear();
590 // pip.push_back(Pip.get(1)); pip.push_back(Pip.get(2)); pip.push_back(Pip.get(3));
591 // pip.push_back(Pip.get(0)); std::vector<double> pim; pim.clear();
592 // pim.push_back(Pim.get(1)); pim.push_back(Pim.get(2)); pim.push_back(Pim.get(3));
593 // pim.push_back(Pim.get(0)); std::vector<double> pi0; pi0.clear();
594 // pi0.push_back(Pi0.get(1)); pi0.push_back(Pi0.get(2)); pi0.push_back(Pi0.get(3));
595 // pi0.push_back(Pi0.get(0));
596
597 propagator["kstarm_892"] =
598 create_RBW_propagator( "kstarm_892", m2_ks0pim, m2_ks0, m2_pim, 1 );
599 propagator["kstar0_892"] =
600 create_RBW_propagator( "kstar0_892", m2_ks0pi0, m2_ks0, m2_pi0, 1 );
601 propagator["kstarp_892"] =
602 create_RBW_propagator( "kstarp_892", m2_ks0pip, m2_ks0, m2_pip, 1 );
603 propagator["k1m_1270"] = create_BW_propagator( "k1m_1270", m2_ks0pimpi0 );
604 propagator["k10_1270"] = create_BW_propagator( "k10_1270", m2_ks0pippim );
605 propagator["k1m_1400"] = create_BW_propagator( "k1m_1400", m2_ks0pimpi0 );
606 propagator["k10_1400"] = create_BW_propagator( "k10_1400", m2_ks0pippim );
607 propagator["kstarm2_1410"] = create_BW_propagator( "kstarm_1410", m2_ks0pim );
608 propagator["kstar02_1410"] = create_BW_propagator( "kstar0_1410", m2_ks0pi0 );
609 propagator["kstarm3_1410"] = create_BW_propagator( "kstarm_1410", m2_ks0pimpi0 );
610 propagator["kstar03_1410"] = create_BW_propagator( "kstar0_1410", m2_ks0pippim );
611 propagator["k0starm_1430"] =
612 create_KPiSLASS_propagator( "k0starm_1430", m2_ks0pim, m2_ks0, m2_pim );
613 propagator["k0star0_1430"] =
614 create_KPiSLASS_propagator( "k0star0_1430", m2_ks0pi0, m2_ks0, m2_pi0 );
615 // propagator["k0starm_1430" ] = create_BW_propagator("k0starm_1430", m2_ks0pim);
616 // propagator["k0star0_1430" ] = create_BW_propagator("k0star0_1430", m2_ks0pi0);
617 propagator["k2starm2_1430"] = create_BW_propagator( "k2starm_1430", m2_ks0pim );
618 propagator["k2star02_1430"] = create_BW_propagator( "k2star0_1430", m2_ks0pi0 );
619 propagator["k2starm3_1430"] = create_BW_propagator( "k2starm_1430", m2_ks0pimpi0 );
620 propagator["k2star03_1430"] = create_BW_propagator( "k2star0_1430", m2_ks0pippim );
621 propagator["km_1460"] = create_BW_propagator( "km_1460", m2_ks0pimpi0 );
622 propagator["k0_1460"] = create_BW_propagator( "k0_1460", m2_ks0pippim );
623 propagator["k2m_1580"] = create_BW_propagator( "k2m_1580", m2_ks0pimpi0 );
624 propagator["k20_1580"] = create_BW_propagator( "k20_1580", m2_ks0pippim );
625 propagator["km_1630"] = create_BW_propagator( "km_1630", m2_ks0pimpi0 );
626 propagator["k0_1630"] = create_BW_propagator( "k0_1630", m2_ks0pippim );
627 propagator["k1m_1650"] = create_BW_propagator( "k1m_1650", m2_ks0pimpi0 );
628 propagator["k10_1650"] = create_BW_propagator( "k10_1650", m2_ks0pippim );
629 propagator["kstarm2_1680"] = create_BW_propagator( "kstarm_1680", m2_ks0pim );
630 propagator["kstar02_1680"] = create_BW_propagator( "kstar0_1680", m2_ks0pi0 );
631 propagator["kstarm3_1680"] = create_BW_propagator( "kstarm_1680", m2_ks0pimpi0 );
632 propagator["kstar03_1680"] = create_BW_propagator( "kstar0_1680", m2_ks0pippim );
633 propagator["k2m_1770"] = create_BW_propagator( "k2m_1770", m2_ks0pimpi0 );
634 propagator["k20_1770"] = create_BW_propagator( "k20_1770", m2_ks0pippim );
635
636 propagator["sigma_500"] = create_sigma_propagator( "sigma_500", m2_pippim );
637 propagator["rhop_770"] = create_RBW_propagator( "rhop_770", m2_pippi0, m2_pip, m2_pi0, 1 );
638 propagator["rhom_770"] = create_RBW_propagator( "rhom_770", m2_pimpi0, m2_pim, m2_pi0, 1 );
639 propagator["rho0_770"] = create_GS_propagator( "rho0_770", m2_pippim );
640 propagator["omega2_782"] = create_BW_propagator( "omega_782", m2_pippim );
641 propagator["f0_980"] =
642 create_Flatte2_propagator( "f0_980", m2_pippim, mpion, mpion, mkaon, mkaon );
643 propagator["f0_1370"] = create_BW_propagator( "f0_1370", m2_pippim );
644 propagator["f2_1270"] = create_BW_propagator( "f2_1270", m2_pippim );
645
646 // // The resolution should be considered in the fit, but when calculating the physical
647 // values, we could only use the pure itude
648 // //propagator["phi_1020" ] = new GPUPropagatorBreitWigner("phi_1021", m2_pippimpi0,
649 // mean, sigma);
650 propagator["phi_1020"] = create_BW_propagator( "phi_1020", m2_pippimpi0 );
651 // //propagator["omega3_782" ] = new GPUPropagatorBreitWigner("omega_782", m2_pippimpi0,
652 // mean, sigma);
653 propagator["omega3_782"] = create_BW_propagator( "omega_782", m2_pippimpi0 );
654 // //propagator["eta_547" ] = new GPUPropagatorBreitWigner("eta_547", m2_pippimpi0,
655 // mean, sigma);
656 propagator["eta_547"] = create_BW_propagator( "eta_547", m2_pippimpi0 );
657
658 // propagator["a10_1260_rhom_770_0"] = create_RBW_propagator("a10_1260", m2_pippimpi0,
659 // m2_pimpi0, m2_pip, 0); propagator["a10_1260_rhop_770_0"] =
660 // create_RBW_propagator("a10_1260", m2_pippimpi0, m2_pippi0, m2_pim, 0);
661 // propagator["a10_1260_rhom_770_2"] = create_RBW_propagator("a10_1260", m2_pippimpi0,
662 // m2_pimpi0, m2_pip, 2); propagator["a10_1260_rhop_770_2"] =
663 // create_RBW_propagator("a10_1260", m2_pippimpi0, m2_pippi0, m2_pim, 2);
664 //
665 propagator["kstarm_phsp"] = create_BW_propagator( "kstarm_phsp", m2_ks0pim );
666 propagator["kstar0_phsp"] = create_BW_propagator( "kstar0_phsp", m2_ks0pi0 );
667 propagator["kstarp_phsp"] = create_BW_propagator( "kstarp_phsp", m2_ks0pip );
668 propagator["rhop_phsp"] = create_BW_propagator( "rhop_phsp", m2_pippi0 );
669 propagator["rhom_phsp"] = create_BW_propagator( "rhom_phsp", m2_pimpi0 );
670 propagator["rho0_phsp"] = create_BW_propagator( "rho0_phsp", m2_pippim );
671 propagator["k1m_phsp"] = create_BW_propagator( "k1m_phsp", m2_ks0pimpi0 );
672 propagator["k10_phsp"] = create_BW_propagator( "k10_phsp", m2_ks0pippim );
673 propagator["k1p_phsp"] = create_BW_propagator( "k1p_phsp", m2_ks0pippi0 );
674 propagator["a10_phsp"] = create_BW_propagator( "a10_phsp", m2_pippimpi0 );
675}
676
677//============================================================//
678// Tensor Contract Function //
679//============================================================//
680// Four-momentum (PX, PY, PZ; E)
681double D0ToKSpipipi0::contract_11_0( vector<double> pa, vector<double> pb ) {
682 assert( pa.size() == pb.size() && pa.size() == 4 );
683
684 double temp = pa[3] * pb[3] - pa[0] * pb[0] - pa[1] * pb[1] - pa[2] * pb[2];
685 return temp;
686}
687
688vector<double> D0ToKSpipipi0::contract_21_1( vector<double> pa, vector<double> pb ) {
689 assert( pa.size() == 16 && pb.size() == 4 );
690
691 vector<double> temp;
692 temp.clear();
693 for ( int i = 0; i < 4; i++ )
694 {
695 double sum = 0;
696 for ( int j = 0; j < 4; j++ )
697 {
698 int idx = i * 4 + j;
699 sum += pa[idx] * pb[j] * g[4 * j + j];
700 }
701 temp.push_back( sum );
702 }
703 return temp;
704}
705
706double D0ToKSpipipi0::contract_22_0( vector<double> pa, vector<double> pb ) {
707 assert( pa.size() == 16 && pb.size() == 16 );
708
709 double temp = 0;
710 for ( int i = 0; i < 4; i++ )
711 {
712 for ( int j = 0; j < 4; j++ )
713 {
714 int idx = i * 4 + j;
715 temp += pa[idx] * pb[idx] * g[4 * i + i] * g[4 * j + j];
716 }
717 }
718 return temp;
719}
720
721vector<double> D0ToKSpipipi0::contract_31_2( vector<double> pa, vector<double> pb ) {
722 assert( pa.size() == 64 && pb.size() == 4 );
723
724 vector<double> temp;
725 temp.clear();
726 for ( int i = 0; i < 16; i++ )
727 {
728 double sum = 0;
729 for ( int j = 0; j < 4; j++ )
730 {
731 int idx = i * 4 + j;
732 sum += pa[idx] * pb[j] * g[4 * j + j];
733 }
734 temp.push_back( sum );
735 }
736 return temp;
737}
738
739vector<double> D0ToKSpipipi0::contract_41_3( vector<double> pa, vector<double> pb ) {
740 assert( pa.size() == 256 && pb.size() == 4 );
741
742 vector<double> temp;
743 temp.clear();
744 for ( int i = 0; i < 64; i++ )
745 {
746 double sum = 0;
747 for ( int j = 0; j < 4; j++ )
748 {
749 int idx = i * 4 + j;
750 sum += pa[idx] * pb[j] * g[4 * j + j];
751 }
752 temp.push_back( sum );
753 }
754 return temp;
755}
756
757vector<double> D0ToKSpipipi0::contract_42_2( vector<double> pa, vector<double> pb ) {
758 assert( pa.size() == 256 && pb.size() == 16 );
759
760 vector<double> temp;
761 temp.clear();
762 for ( int i = 0; i < 16; i++ )
763 {
764 double sum = 0;
765 for ( int j = 0; j < 4; j++ )
766 {
767 for ( int k = 0; k < 4; k++ )
768 {
769 int idxa = i * 16 + j * 4 + k;
770 int idxb = j * 4 + k;
771 sum += pa[idxa] * pb[idxb] * g[4 * j + j] * g[4 * k + k];
772 }
773 }
774 temp.push_back( sum );
775 }
776 return temp;
777}
778
779vector<double> D0ToKSpipipi0::contract_22_2( vector<double> pa, vector<double> pb ) {
780 assert( pa.size() == 16 && pb.size() == 16 );
781
782 vector<double> temp;
783 temp.clear();
784 for ( int i = 0; i < 4; i++ )
785 {
786 for ( int j = 0; j < 4; j++ )
787 {
788 double sum = 0;
789 for ( int k = 0; k < 4; k++ )
790 {
791 int idxa = i * 4 + k;
792 int idxb = j * 4 + k;
793 sum += pa[idxa] * pb[idxb] * g[4 * k + k];
794 }
795 temp.push_back( sum );
796 }
797 }
798 return temp;
799}
800
801//------------------------------------------------------------//
802// Spin Factor //
803//------------------------------------------------------------//
804vector<double> D0ToKSpipipi0::Proj( vector<double> pa, int rank ) {
805 if ( pa.size() != 4 )
806 {
807 cout << "Error: pa" << endl;
808 exit( 0 );
809 }
810 if ( rank < 0 )
811 {
812 cout << "Error: L<0 !!!" << endl;
813 exit( 0 );
814 }
815
816 double msa = contract_11_0( pa, pa );
817
818 // Projection Operator 2-rank
819 vector<double> proj_uv;
820 proj_uv.clear();
821 for ( int i = 0; i < 4; i++ )
822 {
823 for ( int j = 0; j < 4; j++ )
824 {
825 int idx = i * 4 + j;
826 double temp = -g[idx] + pa[i] * pa[j] / msa;
827 proj_uv.push_back( temp );
828 }
829 }
830
831 vector<double> proj;
832 proj.clear();
833 // Orbital Tensors
834 if ( rank == 0 )
835 {
836 vector<double> t;
837 t.clear();
838 t.push_back( 1.0 );
839 proj = t;
840 }
841 else if ( rank == 1 ) { proj = proj_uv; }
842 else if ( rank == 2 )
843 {
844 vector<double> proj_uvmn;
845 proj_uvmn.clear();
846 for ( int i = 0; i < 4; i++ )
847 {
848 for ( int j = 0; j < 4; j++ )
849 {
850 for ( int k = 0; k < 4; k++ )
851 {
852 for ( int l = 0; l < 4; l++ )
853 {
854
855 int idx1_1 = 4 * i + k;
856 int idx1_2 = 4 * i + l;
857 int idx1_3 = 4 * i + j;
858
859 int idx2_1 = 4 * j + l;
860 int idx2_2 = 4 * j + k;
861 int idx2_3 = 4 * k + l;
862
863 double temp = ( 1.0 / 2.0 ) * ( proj_uv[idx1_1] * proj_uv[idx2_1] +
864 proj_uv[idx1_2] * proj_uv[idx2_2] ) -
865 ( 1.0 / 3.0 ) * proj_uv[idx1_3] * proj_uv[idx2_3];
866 proj_uvmn.push_back( temp );
867 }
868 }
869 }
870 }
871 proj = proj_uvmn;
872 }
873 else
874 {
875 cout << "rank>2: please add it by yourself!!!" << endl;
876 exit( 0 );
877 }
878
879 return proj;
880}
881
882// Orbital Tensor
883vector<double> D0ToKSpipipi0::OrbitalTensors( vector<double> pa, vector<double> pb,
884 vector<double> pc, double r, int rank ) {
885 assert( pa.size() == 4 && pb.size() == 4 && pc.size() == 4 );
886 if ( rank < 0 )
887 {
888 cout << "Error: L<0 !!!" << endl;
889 exit( 0 );
890 }
891
892 // relative momentum
893 vector<double> mr;
894 mr.clear();
895
896 for ( int i = 0; i < 4; i++ )
897 {
898 double temp = pb[i] - pc[i];
899 mr.push_back( temp );
900 }
901
902 // "Masses" of particles
903 double msa = contract_11_0( pa, pa );
904 double msb = contract_11_0( pb, pb );
905 double msc = contract_11_0( pc, pc );
906
907 // Q^2_abc
908 double top = msa + msb - msc;
909 double Q2abc = top * top / ( 4.0 * msa ) - msb;
910
911 // Barrier factors
912 // double Q_0 = 0.197321f/r; //mRadius with unit of fm;
913 double Q_0 = 1.0 / r; // mRadius with unit of GeV^-1;
914 double Q_02 = Q_0 * Q_0;
915 double Q_04 = Q_02 * Q_02;
916 // double Q_06 = Q_04*Q_02;
917 // double Q_08 = Q_04*Q_04;
918
919 double Q4abc = Q2abc * Q2abc;
920 // double Q6abc = Q4abc*Q2abc;
921 // double Q8abc = Q4abc*Q4abc;
922
923 double mB1 = sqrt( 2.0 / ( Q2abc + Q_02 ) );
924 double mB2 = sqrt( 13.0 / ( Q4abc + ( 3.0 * Q_02 ) * Q2abc + 9.0 * Q_04 ) );
925 // mB3 = &sqrt(277.0f/(Q6abc + (6.0f*Q_02)*Q4abc + (45.0f*Q_04)*Q2abc + 225.0f*Q_06));
926 // mB4 = &sqrt(12746.0f/(Q8abc + (10.0f*Q_02)*Q6abc + (135.0f*Q_04)*Q4abc +
927 // (1575.0f*Q_06)*Q2abc + 11025.0f*Q_08));
928
929 // Projection Operator 2-rank
930 vector<double> proj_uv;
931 proj_uv.clear();
932 for ( int i = 0; i < 4; i++ )
933 {
934 for ( int j = 0; j < 4; j++ )
935 {
936 int idx = i * 4 + j;
937 double temp = -g[idx] + pa[i] * pa[j] / msa;
938 proj_uv.push_back( temp );
939 }
940 }
941
942 vector<double> Bt;
943 Bt.clear();
944 // Orbital Tensors
945 if ( rank == 0 )
946 {
947 vector<double> t;
948 t.clear();
949 t.push_back( 1.0 );
950 Bt = t;
951 }
952 else if ( rank < 3 )
953 {
954 vector<double> t_u;
955 t_u.clear();
956 vector<double> Bt_u;
957 Bt_u.clear();
958 for ( int i = 0; i < 4; i++ )
959 {
960 double temp = 0;
961 for ( int j = 0; j < 4; j++ )
962 {
963 int idx = i * 4 + j;
964 temp += -proj_uv[idx] * mr[j] * g[j * 4 + j];
965 }
966 t_u.push_back( temp );
967 Bt_u.push_back( temp * mB1 );
968 }
969 if ( rank == 1 ) Bt = Bt_u;
970
971 double t_u2 = contract_11_0( t_u, t_u );
972
973 vector<double> Bt_uv;
974 Bt_uv.clear();
975 for ( int i = 0; i < 4; i++ )
976 {
977 for ( int j = 0; j < 4; j++ )
978 {
979 int idx = 4 * i + j;
980 double temp = t_u[i] * t_u[j] + ( 1.0 / 3.0 ) * proj_uv[idx] * t_u2;
981 Bt_uv.push_back( temp * mB2 );
982 }
983 }
984 if ( rank == 2 ) Bt = Bt_uv;
985 }
986 else
987 {
988 cout << "rank>2: please add it by yorself!!!" << endl;
989 exit( 0 );
990 }
991
992 return Bt;
993}
994
995// The sequence is 0+, 0-, 1-, 1+, 2+, 2-
996//------------------------------------------------------------//
997// Cascade Decay //
998//------------------------------------------------------------//
999//---- D -> P P1 [S], P -> S P2 [S]
1000double D0ToKSpipipi0::D2PP_P2SP() {
1001
1002 double SF_D2PP_P2SP = 1.0;
1003
1004 return SF_D2PP_P2SP;
1005}
1006
1007//----D -> P P1 [S], P -> V P2 [P], V -> P3 P4 [P]
1008double D0ToKSpipipi0::D2PP_P2VP( vector<double> p1, vector<double> p2, vector<double> p3,
1009 vector<double> p4, int l ) {
1010 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1011 assert( l >= 0 );
1012
1013 vector<double> V = p1 + p2;
1014 vector<double> P2 = p3;
1015 vector<double> P = V + P2;
1016 vector<double> P1 = p4;
1017
1018 vector<double> P_orbitals_Spin1OrbitalTensor = OrbitalTensors( P, V, P2, rRes, 1 );
1019 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p1, p2, rRes, 1 );
1020
1021 double SF_D2PP_P2VP = 0.;
1022 SF_D2PP_P2VP = contract_11_0( P_orbitals_Spin1OrbitalTensor, V_orbitals_Spin1OrbitalTensor );
1023 return SF_D2PP_P2VP;
1024}
1025
1026//----D -> V1 P1 [P], V1 -> V2 P2 [P], V2 -> P3 P4 [P]
1027double D0ToKSpipipi0::D2VP_V2VP( vector<double> p1, vector<double> p2, vector<double> p3,
1028 vector<double> p4, int l ) {
1029 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1030 assert( l >= 0 );
1031
1032 vector<double> V2 = p1 + p2;
1033 vector<double> P2 = p3;
1034 vector<double> V1 = V2 + P2;
1035 vector<double> P1 = p4;
1036 vector<double> D = V1 + P1;
1037
1038 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, V1, P1, rD, 1 );
1039 vector<double> V1_orbitals_Spin1OrbitalTensor = OrbitalTensors( V1, V2, P2, rRes, 1 );
1040 vector<double> V2_orbitals_Spin1OrbitalTensor = OrbitalTensors( V2, p1, p2, rRes, 1 );
1041 vector<double> Proj1_V1 = Proj( V1, 1 );
1042
1043 double SF_D2VP_V2VP = 0.;
1044 SF_D2VP_V2VP = contract_11_0( contract_21_1( contract_31_2( contract_41_3( epsilon, V1 ),
1045 V2_orbitals_Spin1OrbitalTensor ),
1046 V1_orbitals_Spin1OrbitalTensor ),
1047 contract_21_1( Proj1_V1, D_orbitals_Spin1OrbitalTensor ) );
1048 // SF_D2VP_V2VP = contract_11_0( contract_21_1(contract_31_2(contract_41_3(epsilon, V1),
1049 // V1_orbitals_Spin1OrbitalTensor) | contract_21_1(Proj1_V1, V2_orbitals_Spin1OrbitalTensor))
1050 // | contract_21_1(Proj1_V1, D_orbitals_Spin1OrbitalTensor) );//alternative
1051 return SF_D2VP_V2VP;
1052}
1053
1054//----D -> A P1 [P], A -> S P2 [P], S -> P3 P4 [S]
1055double D0ToKSpipipi0::D2AP_A2SP( vector<double> p1, vector<double> p2, vector<double> p3,
1056 vector<double> p4, int l ) {
1057 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1058 assert( l >= 0 );
1059
1060 vector<double> S = p1 + p2;
1061 vector<double> P2 = p3;
1062 vector<double> A = S + P2;
1063 vector<double> P1 = p4;
1064 vector<double> D = A + P1;
1065
1066 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, A, P1, rD, 1 );
1067 vector<double> A_orbitals_Spin1OrbitalTensor = OrbitalTensors( A, S, P2, rRes, 1 );
1068
1069 double SF_D2AP_A2SP = 0.;
1070 SF_D2AP_A2SP = contract_11_0( D_orbitals_Spin1OrbitalTensor, A_orbitals_Spin1OrbitalTensor );
1071 return SF_D2AP_A2SP;
1072}
1073
1074//----D -> A P1 [P], A -> V P2 [S, D], V -> P3 P4 [P]
1075double D0ToKSpipipi0::D2AP_A2VP( vector<double> p1, vector<double> p2, vector<double> p3,
1076 vector<double> p4, int l ) {
1077 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1078 assert( l >= 0 );
1079
1080 vector<double> V = p1 + p2;
1081 vector<double> P2 = p3;
1082 vector<double> A = V + P2;
1083 vector<double> P1 = p4;
1084 vector<double> D = A + P1;
1085
1086 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, A, P1, rD, 1 );
1087 vector<double> A_orbitals_Spin2OrbitalTensor = OrbitalTensors( A, V, P2, rRes, 2 );
1088 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p1, p2, rRes, 1 );
1089 vector<double> Proj1_A = Proj( A, 1 );
1090
1091 double SF_D2AP_A2VP = 0.;
1092 if ( l == 0 )
1093 {
1094 SF_D2AP_A2VP = contract_11_0( D_orbitals_Spin1OrbitalTensor,
1095 contract_21_1( Proj1_A, V_orbitals_Spin1OrbitalTensor ) );
1096 }
1097 if ( l == 2 )
1098 {
1099 SF_D2AP_A2VP = contract_11_0(
1100 D_orbitals_Spin1OrbitalTensor,
1101 contract_21_1( A_orbitals_Spin2OrbitalTensor, V_orbitals_Spin1OrbitalTensor ) );
1102 }
1103 return SF_D2AP_A2VP;
1104}
1105
1106//----D -> A P1 [P], A -> T P2 [P], T -> P3 P4 [D]
1107double D0ToKSpipipi0::D2AP_A2TP( vector<double> p1, vector<double> p2, vector<double> p3,
1108 vector<double> p4, int l ) {
1109 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1110 assert( l >= 0 );
1111
1112 vector<double> T = p1 + p2;
1113 vector<double> P2 = p3;
1114 vector<double> A = T + P2;
1115 vector<double> P1 = p4;
1116 vector<double> D = A + P1;
1117
1118 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, A, P1, rD, 1 );
1119 vector<double> A_orbitals_Spin1OrbitalTensor = OrbitalTensors( A, T, P2, rRes, 1 );
1120 vector<double> T_orbitals_Spin2OrbitalTensor = OrbitalTensors( T, p1, p2, rRes, 2 );
1121 vector<double> Proj2_A = Proj( A, 2 );
1122
1123 double SF_D2AP_A2TP = 0.;
1124 // SF_D2AP_A2TP = &(D_orbitals.Spin1OrbitalTensor() | (T_orbitals.Spin2OrbitalTensor() |
1125 // A_orbitals.Spin1OrbitalTensor()));
1126 SF_D2AP_A2TP =
1127 contract_11_0( D_orbitals_Spin1OrbitalTensor,
1128 contract_21_1( contract_42_2( Proj2_A, T_orbitals_Spin2OrbitalTensor ),
1129 A_orbitals_Spin1OrbitalTensor ) );
1130 return SF_D2AP_A2TP;
1131}
1132
1133//----D -> T P1 [D], T -> V P2 [D], V -> P3 P4 [P]
1134double D0ToKSpipipi0::D2TP_T2VP( vector<double> p1, vector<double> p2, vector<double> p3,
1135 vector<double> p4, int l ) {
1136 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1137 assert( l >= 0 );
1138
1139 vector<double> V = p1 + p2;
1140 vector<double> P2 = p3;
1141 vector<double> T = V + P2;
1142 vector<double> P1 = p4;
1143 vector<double> D = T + P1;
1144
1145 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, T, P1, rD, 2 );
1146 vector<double> T_orbitals_Spin2OrbitalTensor = OrbitalTensors( T, V, P2, rRes, 2 );
1147 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p1, p2, rRes, 1 );
1148 vector<double> Proj1_T = Proj( T, 1 );
1149 vector<double> Proj2_T = Proj( T, 2 );
1150
1151 double SF_D2TP_T2VP = 0.;
1152 SF_D2TP_T2VP = contract_22_0(
1153 contract_42_2( Proj2_T, D_orbitals_Spin2OrbitalTensor ),
1154 contract_22_2( contract_31_2( contract_41_3( epsilon, T ),
1155 contract_21_1( Proj1_T, V_orbitals_Spin1OrbitalTensor ) ),
1156 T_orbitals_Spin2OrbitalTensor ) );
1157 return SF_D2TP_T2VP;
1158}
1159
1160//----D -> T1 P1 [D], T1 -> T2 P2 [P], T2 -> P3 P4 [D]
1161double D0ToKSpipipi0::D2TP_T2TP( vector<double> p1, vector<double> p2, vector<double> p3,
1162 vector<double> p4, int l ) {
1163 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1164 assert( l >= 0 );
1165
1166 vector<double> T2 = p1 + p2;
1167 vector<double> P2 = p3;
1168 vector<double> T1 = T2 + P2;
1169 vector<double> P1 = p4;
1170 vector<double> D = T1 + P1;
1171
1172 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, T1, P1, rD, 2 );
1173 vector<double> T1_orbitals_Spin1OrbitalTensor = OrbitalTensors( T1, T2, P2, rRes, 1 );
1174 vector<double> T2_orbitals_Spin2OrbitalTensor = OrbitalTensors( T2, p1, p2, rRes, 2 );
1175 vector<double> Proj2_T1 = Proj( T1, 2 );
1176
1177 double SF_D2TP_T2TP = 0.;
1178 SF_D2TP_T2TP = contract_22_0( contract_22_2( contract_31_2( contract_41_3( epsilon, T1 ),
1179 T1_orbitals_Spin1OrbitalTensor ),
1180 T2_orbitals_Spin2OrbitalTensor ),
1181 contract_42_2( Proj2_T1, D_orbitals_Spin2OrbitalTensor ) );
1182 // SF_D2TP_T2TP = &( (((epsilon | T1) | T1_orbitals.Spin1OrbitalTensor()) || (Proj2_T1 |
1183 // D_orbitals.Spin2OrbitalTensor())) | T2_orbitals.Spin2OrbitalTensor() );//alternative
1184 return SF_D2TP_T2TP;
1185}
1186
1187//----D -> PT P1 [D], PT -> S P2 [D], S -> P3 P4 [S]
1188double D0ToKSpipipi0::D2PTP_PT2SP( vector<double> p1, vector<double> p2, vector<double> p3,
1189 vector<double> p4, int l ) {
1190 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1191 assert( l >= 0 );
1192
1193 vector<double> S = p1 + p2;
1194 vector<double> P2 = p3;
1195 vector<double> PT = S + P2;
1196 vector<double> P1 = p4;
1197 vector<double> D = PT + P1;
1198
1199 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, PT, P1, rD, 2 );
1200 vector<double> PT_orbitals_Spin2OrbitalTensor = OrbitalTensors( PT, S, P2, rRes, 2 );
1201
1202 double SF_D2PTP_PT2SP = 0.;
1203 SF_D2PTP_PT2SP =
1204 contract_22_0( D_orbitals_Spin2OrbitalTensor, PT_orbitals_Spin2OrbitalTensor );
1205 return SF_D2PTP_PT2SP;
1206}
1207
1208//----D -> PT P1 [D], PT -> V P2 [P], V -> P3 P4 [P]
1209double D0ToKSpipipi0::D2PTP_PT2VP( vector<double> p1, vector<double> p2, vector<double> p3,
1210 vector<double> p4, int l ) {
1211 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1212 assert( l >= 0 );
1213
1214 vector<double> V = p1 + p2;
1215 vector<double> P2 = p3;
1216 vector<double> PT = V + P2;
1217 vector<double> P1 = p4;
1218 vector<double> D = PT + P1;
1219
1220 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, PT, P1, rD, 2 );
1221 vector<double> PT_orbitals_Spin1OrbitalTensor = OrbitalTensors( PT, V, P2, rRes, 1 );
1222 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p1, p2, rRes, 1 );
1223 vector<double> Proj2_PT = Proj( PT, 2 );
1224
1225 double SF_D2PTP_PT2VP = 0.;
1226 SF_D2PTP_PT2VP =
1227 contract_22_0( D_orbitals_Spin2OrbitalTensor,
1228 contract_31_2( contract_41_3( Proj2_PT, V_orbitals_Spin1OrbitalTensor ),
1229 PT_orbitals_Spin1OrbitalTensor ) );
1230 return SF_D2PTP_PT2VP;
1231}
1232
1233//----D -> PT P1 [D], PT -> T P2 [S]
1234double D0ToKSpipipi0::D2PTP_PT2TP( vector<double> p1, vector<double> p2, vector<double> p3,
1235 vector<double> p4, int l ) {
1236 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1237 assert( l >= 0 );
1238
1239 vector<double> T = p1 + p2;
1240 vector<double> P2 = p3;
1241 vector<double> PT = T + P2;
1242 vector<double> P1 = p4;
1243 vector<double> D = PT + P1;
1244
1245 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, PT, P1, rD, 2 );
1246 vector<double> T_orbitals_Spin2OrbitalTensor = OrbitalTensors( T, p1, p2, rRes, 2 );
1247 vector<double> Proj2_PT = Proj( PT, 2 );
1248
1249 double SF_D2PTP_PT2TP = 0.;
1250 SF_D2PTP_PT2TP = contract_22_0( D_orbitals_Spin2OrbitalTensor,
1251 contract_42_2( Proj2_PT, T_orbitals_Spin2OrbitalTensor ) );
1252 return SF_D2PTP_PT2TP;
1253}
1254
1255//------------------------------------------------------------//
1256// Quasi Two-Body Decay //
1257//------------------------------------------------------------//
1258//----D -> S1 S2 [S]
1259double D0ToKSpipipi0::D2SS() {
1260
1261 double SF_D2SS = 1.0;
1262 return SF_D2SS;
1263}
1264
1265//----D->V S [P]
1266double D0ToKSpipipi0::D2VS( vector<double> p1, vector<double> p2, vector<double> p3,
1267 vector<double> p4, int l ) {
1268 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1269 assert( l >= 0 );
1270
1271 vector<double> V = p1 + p2;
1272 vector<double> S = p3 + p4;
1273 vector<double> D = V + S;
1274
1275 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, V, S, rD, 1 );
1276 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p1, p2, rRes, 1 );
1277
1278 double SF_D2VS = 0.;
1279 SF_D2VS = contract_11_0( D_orbitals_Spin1OrbitalTensor, V_orbitals_Spin1OrbitalTensor );
1280 return SF_D2VS;
1281}
1282
1283//----D -> V1 V2 [S, P, D], V1 -> P1 P2 [P], V2 -> P3 P4 [P]
1284double D0ToKSpipipi0::D2VV( vector<double> p1, vector<double> p2, vector<double> p3,
1285 vector<double> p4, int l ) {
1286 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1287 assert( l >= 0 );
1288
1289 vector<double> V1 = p1 + p2;
1290 vector<double> V2 = p3 + p4;
1291 vector<double> D = V1 + V2;
1292
1293 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, V1, V2, rD, 1 );
1294 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, V1, V2, rD, 2 );
1295 vector<double> V1_orbitals_Spin1OrbitalTensor = OrbitalTensors( V1, p1, p2, rRes, 1 );
1296 vector<double> V2_orbitals_Spin1OrbitalTensor = OrbitalTensors( V2, p3, p4, rRes, 1 );
1297
1298 double SF_D2VV = 0;
1299 if ( l == 0 )
1300 {
1301 SF_D2VV = contract_11_0( V1_orbitals_Spin1OrbitalTensor, V2_orbitals_Spin1OrbitalTensor );
1302 }
1303 if ( l == 1 )
1304 {
1305 SF_D2VV = contract_11_0( contract_21_1( contract_31_2( contract_41_3( epsilon, D ),
1306 V2_orbitals_Spin1OrbitalTensor ),
1307 V1_orbitals_Spin1OrbitalTensor ),
1308 D_orbitals_Spin1OrbitalTensor );
1309 }
1310 if ( l == 2 )
1311 {
1312 SF_D2VV = contract_11_0(
1313 contract_21_1( D_orbitals_Spin2OrbitalTensor, V2_orbitals_Spin1OrbitalTensor ),
1314 V1_orbitals_Spin1OrbitalTensor );
1315 }
1316
1317 return SF_D2VV;
1318}
1319
1320//----D -> T S [D], T -> P1 P2 [D], S -> P3 P4 [S]
1321double D0ToKSpipipi0::D2TS( vector<double> p1, vector<double> p2, vector<double> p3,
1322 vector<double> p4, int l ) {
1323 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1324 assert( l >= 0 );
1325
1326 vector<double> T = p1 + p2;
1327 vector<double> S = p3 + p4;
1328 vector<double> D = T + S;
1329
1330 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, T, S, rD, 2 );
1331 vector<double> T_orbitals_Spin2OrbitalTensor = OrbitalTensors( T, p1, p2, rRes, 2 );
1332
1333 double SF_D2TS = 0.;
1334 SF_D2TS = contract_22_0( D_orbitals_Spin2OrbitalTensor, T_orbitals_Spin2OrbitalTensor );
1335 return SF_D2TS;
1336}
1337
1338//----D -> T V [P, D], T -> P1 P2 [D], V -> P3 P4 [P]
1339double D0ToKSpipipi0::D2TV( vector<double> p1, vector<double> p2, vector<double> p3,
1340 vector<double> p4, int l ) {
1341 assert( p1.size() == 4 && p2.size() == 4 && p3.size() == 4 && p4.size() == 4 );
1342 assert( l >= 0 );
1343
1344 vector<double> T = p1 + p2;
1345 vector<double> V = p3 + p4;
1346 vector<double> D = T + V;
1347
1348 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, T, V, rD, 1 );
1349 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, T, V, rD, 2 );
1350 vector<double> T_orbitals_Spin2OrbitalTensor = OrbitalTensors( T, p1, p2, rRes, 2 );
1351 vector<double> V_orbitals_Spin1OrbitalTensor = OrbitalTensors( V, p3, p4, rRes, 1 );
1352
1353 double SF_D2TV = 0.;
1354 if ( l == 1 )
1355 {
1356 SF_D2TV = contract_11_0(
1357 D_orbitals_Spin1OrbitalTensor,
1358 contract_21_1( T_orbitals_Spin2OrbitalTensor, V_orbitals_Spin1OrbitalTensor ) );
1359 }
1360 if ( l == 2 )
1361 {
1362 SF_D2TV = contract_22_0(
1363 contract_31_2( contract_41_3( epsilon, D ), V_orbitals_Spin1OrbitalTensor ),
1364 contract_22_2( D_orbitals_Spin2OrbitalTensor, T_orbitals_Spin2OrbitalTensor ) );
1365 }
1366 return SF_D2TV;
1367}
1368
1369//----D -> T1 T2 [S, P, D], T1 -> P1 P2 [D], T2 -> P3 P4 [D]
1370// ok
1371double D0ToKSpipipi0::D2TT( vector<double> p1, vector<double> p2, vector<double> p3,
1372 vector<double> p4, int l ) {
1373
1374 vector<double> T1 = p1 + p2;
1375 vector<double> T2 = p3 + p4;
1376 vector<double> D = T1 + T2;
1377
1378 vector<double> D_orbitals_Spin1OrbitalTensor = OrbitalTensors( D, T1, T2, rD, 1 );
1379 vector<double> D_orbitals_Spin2OrbitalTensor = OrbitalTensors( D, T1, T2, rD, 2 );
1380 vector<double> T1_orbitals_Spin2OrbitalTensor = OrbitalTensors( T1, p1, p2, rRes, 2 );
1381 vector<double> T2_orbitals_Spin2OrbitalTensor = OrbitalTensors( T2, p3, p4, rRes, 2 );
1382
1383 double SF_D2TT = 0.;
1384 if ( l == 0 )
1385 {
1386 SF_D2TT = contract_22_0( T1_orbitals_Spin2OrbitalTensor, T2_orbitals_Spin2OrbitalTensor );
1387 }
1388 if ( l == 1 )
1389 {
1390 SF_D2TT = contract_22_0(
1391 contract_31_2( contract_41_3( epsilon, D ), D_orbitals_Spin1OrbitalTensor ),
1392 contract_22_2( T1_orbitals_Spin2OrbitalTensor, T2_orbitals_Spin2OrbitalTensor ) );
1393 }
1394 if ( l == 2 )
1395 {
1396 SF_D2TT = contract_22_0(
1397 D_orbitals_Spin2OrbitalTensor,
1398 contract_22_2( T1_orbitals_Spin2OrbitalTensor, T2_orbitals_Spin2OrbitalTensor ) );
1399 }
1400 return SF_D2TT;
1401}
1402
1403//------------------------------------------------------------//
1404// Propagators //
1405//------------------------------------------------------------//
1406double D0ToKSpipipi0::fundecaymomentum( double mr2, double m1_2, double m2_2 ) {
1407 double mr = sqrt( mr2 );
1408 double poly = mr2 * mr2 + m1_2 * m1_2 + m2_2 * m2_2 - 2.0 * m1_2 * mr2 - 2.0 * m2_2 * mr2 -
1409 2.0 * m1_2 * m2_2;
1410 double ret = sqrt( poly ) / ( 2.0 * mr );
1411 if ( poly < 0 ) ret = 0.0;
1412 return ret;
1413}
1414
1415double D0ToKSpipipi0::fundecaymomentum2( double mr2, double m1_2, double m2_2 ) {
1416 double poly = mr2 * mr2 + m1_2 * m1_2 + m2_2 * m2_2 - 2.0 * m1_2 * mr2 - 2.0 * m2_2 * mr2 -
1417 2.0 * m1_2 * m2_2;
1418 double ret = poly / ( 4.0 * mr2 );
1419 if ( poly < 0 ) ret = 0.0;
1420 return ret;
1421}
1422
1423double D0ToKSpipipi0::wid( double mass, double sa, double sb, double sc, double r, int l ) {
1424
1425 double widm = 1.0;
1426 double sa0 = mass * mass;
1427 double m = sqrt( sa );
1428 double q = fundecaymomentum2( sa, sb, sc );
1429 double q0 = fundecaymomentum2( sa0, sb, sc );
1430 double z = q * r * r;
1431 double z0 = q0 * r * r;
1432 double F = 0.0;
1433 if ( l == 0 ) F = 1.0;
1434 if ( l == 1 ) F = sqrt( ( 1.0 + z0 ) / ( 1.0 + z ) );
1435 if ( l == 2 ) F = sqrt( ( 9.0 + 3.0 * z0 + z0 * z0 ) / ( 9.0 + 3.0 * z + z * z ) );
1436 if ( l == 3 )
1437 F = sqrt( ( 225.0 + 45.0 * z0 + 6.0 * z0 * z0 + z0 * z0 * z0 ) /
1438 ( 225.0 + 45.0 * z + 6.0 * z * z + z * z * z ) );
1439 if ( l == 4 )
1440 F = sqrt(
1441 ( 11025.0 + 1575.0 * z0 + 135.0 * z0 * z0 + 10.0 * z0 * z0 * z0 + z0 * z0 * z0 * z0 ) /
1442 ( 11025.0 + 1575.0 * z + 135.0 * z * z + 10.0 * z * z * z + z * z * z * z ) );
1443
1444 double t = sqrt( q / q0 );
1445 for ( uint i = 0; i < ( 2 * l + 1 ); i++ ) { widm *= t; }
1446 widm *= ( mass / m * F * F );
1447 return widm;
1448}
1449
1450complex<double> D0ToKSpipipi0::BW( double mx2, double mr, double wr ) {
1451
1452 double mr2 = mr * mr;
1453 double denom_real = mr2 - mx2;
1454 double denom_imag = mr * wr;
1455 double denom = denom_real * denom_real + denom_imag * denom_imag;
1456
1457 double output_x = 0., output_y = 0.;
1458 if ( wr < 0 ) { output_x = 1 / sqrt( denom ); }
1459 else if ( wr < 10 )
1460 {
1461 output_x = denom_real / denom;
1462 output_y = denom_imag / denom;
1463 }
1464 else
1465 { /* phase space */
1466 output_x = 1.;
1467 output_y = 0.;
1468 }
1469
1470 complex<double> output( output_x, output_y );
1471
1472 return output;
1473}
1474
1475complex<double> D0ToKSpipipi0::RBW( double mx2, double mr, double wr, double m1_2, double m2_2,
1476 double r, int l ) {
1477
1478 double mr2 = mr * mr;
1479 double denom_real = mr2 - mx2;
1480 double denom_imag = mr * wr * wid( mr, mx2, m1_2, m2_2, r, l ); // real-i*imag;
1481
1482 double denom = denom_real * denom_real + denom_imag * denom_imag;
1483 double output_x = denom_real / denom;
1484 double output_y = denom_imag / denom;
1485
1486 complex<double> output( output_x, output_y );
1487
1488 return output;
1489}
1490
1491/* KPi S-wave LASS Model */
1492complex<double> D0ToKSpipipi0::LASS( double mx2, double m1_2, double m2_2 ) {
1493 double sa = mx2;
1494 double sb = m1_2;
1495 double sc = m2_2;
1496 double m1430 = 1.463;
1497 double w1430 = 0.233;
1498 double sa0 = m1430 * m1430;
1499 double q0 = ( sa0 + sb - sc ) * ( sa0 + sb - sc ) / ( 4.0 * sa0 ) - sb;
1500 if ( q0 < 0 ) q0 = 1e-16;
1501 double qs = ( sa + sb - sc ) * ( sa + sb - sc ) / ( 4.0 * sa ) - sb;
1502 double q = sqrt( qs );
1503 double width = w1430 * q * m1430 / sqrt( sa * q0 );
1504 double temp_R = atan( m1430 * width / ( sa0 - sa ) );
1505 if ( temp_R < 0 ) temp_R += math_pi;
1506 double deltaR = -5.31 + temp_R;
1507 double temp_F = atan( 2.0 * 1.07 * q / ( 2.0 + ( -1.8 ) * 1.07 * qs ) );
1508 if ( temp_F < 0 ) temp_F += math_pi;
1509 double deltaF = 2.33 + temp_F;
1510 double output_x =
1511 0.8 * sin( deltaF ) * cos( deltaF ) + sin( deltaR ) * cos( deltaR + 2.0 * deltaF );
1512 double output_y =
1513 0.8 * sin( deltaF ) * sin( deltaF ) + sin( deltaR ) * sin( deltaR + 2.0 * deltaF );
1514
1515 complex<double> output( output_x, output_y );
1516 return output;
1517}
1518
1519/* GS propagator */
1520const double mass_Pion = 0.13957;
1521double D0ToKSpipipi0::h( double m, double q ) {
1522 double h = 2.0 / math_pi * q / m * log( ( m + 2.0 * q ) / ( 2.0 * mass_Pion ) );
1523 return h;
1524}
1525
1526double D0ToKSpipipi0::dh( double m0, double q0 ) {
1527 double dh = h( m0, q0 ) * ( 1.0 / ( 8.0 * q0 * q0 ) - 1.0 / ( 2.0 * m0 * m0 ) ) +
1528 1.0 / ( 2.0 * math_pi * m0 * m0 );
1529 return dh;
1530}
1531
1532double D0ToKSpipipi0::f( double m0, double sx, double q0, double q ) {
1533 double m = sqrt( sx );
1534 double f =
1535 m0 * m0 / ( q0 * q0 * q0 ) *
1536 ( q * q * ( h( m, q ) - h( m0, q0 ) ) + ( m0 * m0 - sx ) * q0 * q0 * dh( m0, q0 ) );
1537 return f;
1538}
1539
1540double D0ToKSpipipi0::d( double m0, double q0 ) {
1541 double d = 3.0 / math_pi * mass_Pion * mass_Pion / ( q0 * q0 ) *
1542 log( ( m0 + 2.0 * q0 ) / ( 2.0 * mass_Pion ) ) +
1543 m0 / ( 2.0 * math_pi * q0 ) -
1544 ( mass_Pion * mass_Pion * m0 ) / ( math_pi * q0 * q0 * q0 );
1545 return d;
1546}
1547
1548complex<double> D0ToKSpipipi0::GS( double mx2, double mr, double wr, double m1_2, double m2_2,
1549 double r, int l ) {
1550
1551 double mr2 = mr * mr;
1552 double q = fundecaymomentum( mx2, m1_2, m2_2 );
1553 double q0 = fundecaymomentum( mr2, m1_2, m2_2 );
1554 double numer = 1.0 + d( mr, q0 ) * wr / mr;
1555 double denom_real = mr2 - mx2 + wr * f( mr, mx2, q0, q );
1556 double denom_imag = mr * wr * wid( mr, mx2, m1_2, m2_2, r, l ); // real-i*imag;
1557
1558 double denom = denom_real * denom_real + denom_imag * denom_imag;
1559 double output_x = denom_real * numer / denom;
1560 double output_y = denom_imag * numer / denom;
1561
1562 complex<double> output( output_x, output_y );
1563 return output;
1564}
1565
1566complex<double> D0ToKSpipipi0::create_RBW_propagator( string name, double mx2, double m1_2,
1567 double m2_2, int l ) {
1568 double mr = resonance_par[name + "_mass"];
1569 double wr = resonance_par[name + "_width"];
1570 return RBW( mx2, mr, wr, m1_2, m2_2, rRes, l );
1571}
1572
1573complex<double> D0ToKSpipipi0::create_BW_propagator( string name, double mx2 ) {
1574 double mr = resonance_par[name + "_mass"];
1575 double wr = resonance_par[name + "_width"];
1576 return BW( mx2, mr, wr );
1577}
1578
1579complex<double> D0ToKSpipipi0::create_GS_propagator( string name, double mx2 ) { // for rho0
1580 double mr = resonance_par[name + "_mass"];
1581 double wr = resonance_par[name + "_width"];
1582 return GS( mx2, mr, wr, mpip * mpip, mpim * mpim, rRes, 1 );
1583}
1584
1585complex<double> D0ToKSpipipi0::create_KPiSLASS_propagator( string name, double mx2,
1586 double m1_2, double m2_2 ) {
1587 return LASS( mx2, m1_2, m2_2 );
1588}
1589
1590/* sigma propagator for sigma(500)*/
1591double D0ToKSpipipi0::rho4pi( double s ) {
1592 // double mpi=0.13957018;
1593 double mpi = 0.13957;
1594 double tmp = 1.0 - 16.0 * mpi * mpi / s;
1595 double ret = 0.;
1596 if ( tmp <= 0 ) ret = 0.0;
1597 else ret = sqrt( tmp ) / ( 1.0 + exp( 9.8f - 3.5 * s ) );
1598
1599 return ret;
1600}
1601
1602double D0ToKSpipipi0::rho2pi( double s ) {
1603 // double mpi=0.13957018;
1604 double mpi = 0.13957;
1605 double ret = 0.;
1606 if ( s <= 4.0 * mpi * mpi ) ret = 0.0;
1607 else ret = sqrt( 1.0 - 4.0 * mpi * mpi / s );
1608
1609 return ret;
1610}
1611
1612complex<double> D0ToKSpipipi0::sigma( double mx2, double mr, double gf ) {
1613
1614 // double mr = 0.9264f;
1615 double mr2 = mr * mr;
1616 double diff = mr2 - mx2; /* m*m-s */
1617 // double mpi=0.13957018;
1618 double mpi = 0.13957;
1619 double g1 = ( 0.5843 + 1.6663 * mx2 ) * ( mx2 - mpi * mpi / 2.0 ) /
1620 ( mr2 - mpi * mpi / 2.0 ) * exp( -( mx2 - mr2 ) / 1.082 );
1621 // double g2=0.0024f;
1622 double g2 = gf;
1623 double w1 = g1 * rho2pi( mx2 ) / rho2pi( mr2 );
1624 double w2 = g2 * rho4pi( mx2 ) / rho4pi( mr2 );
1625 double ws = mr * ( w1 + w2 );
1626 double denom = diff * diff + ws * ws;
1627
1628 double output_x = diff / denom;
1629 double output_y = ws / denom;
1630
1631 complex<double> output( output_x, output_y );
1632
1633 return output;
1634}
1635
1636complex<double> D0ToKSpipipi0::create_sigma_propagator( string name, double mx2 ) {
1637
1638 double mr = resonance_par[name + "_mass"];
1639 double wr = resonance_par[name + "_width"];
1640
1641 double gf = resonance_par[name + "_gf"];
1642
1643 return sigma( mx2, mr, gf );
1644}
1645
1646/* Flatte propagator for f0(980) and a0(980)*/
1647complex<double> D0ToKSpipipi0::irho( double mr2, double m1_2, double m2_2 ) {
1648 double mr = sqrt( mr2 );
1649 double poly = mr2 * mr2 + m1_2 * m1_2 + m2_2 * m2_2 - 2 * m1_2 * mr2 - 2 * m2_2 * mr2 -
1650 2 * m1_2 * m2_2;
1651
1652 double output_x = 0., output_y = 0.;
1653 if ( poly >= 0 )
1654 {
1655 output_x = 2.0 * sqrt( poly ) / ( 2.0 * mr2 );
1656 output_y = 0.;
1657 }
1658 if ( poly < 0 )
1659 {
1660 output_x = 0.;
1661 output_y = 2.0 * sqrt( -1.0 * poly ) / ( 2.0 * mr2 );
1662 }
1663
1664 complex<double> output( output_x, output_y );
1665
1666 return output;
1667}
1668
1669complex<double> D0ToKSpipipi0::Flatte2( double mx2, double mr2, double g1, double m1a,
1670 double m1b, double g2, double m2a, double m2b ) {
1671
1672 double diff = mr2 - mx2; /* m*m-s */
1673 g2 = g2 * g1;
1674 complex<double> ps1 = irho( mx2, m1a * m1a, m1b * m1b );
1675 complex<double> ps2 = irho( mx2, m2a * m2a, m2b * m2b );
1676 complex<double> iws = g1 * ps1 + g2 * ps2; /*mass dependent width */
1677
1678 diff += imag( iws );
1679 double ws = real( iws );
1680 double denom = diff * diff + ws * ws; /* common denominator*/
1681
1682 double output_x = diff / denom; /* real part*/
1683 double output_y = ws / denom; /* imaginary part*/
1684
1685 complex<double> output( output_x, output_y );
1686
1687 return output;
1688}
1689
1690complex<double> D0ToKSpipipi0::create_Flatte2_propagator( string name, double mx2, double m1a,
1691 double m1b, double m2a,
1692 double m2b ) {
1693 double mr = resonance_par[name + "_mass"];
1694 double wr = resonance_par[name + "_width"];
1695 double mr2 = mr * mr;
1696
1697 double g1 = resonance_par[name + "_g1"];
1698 double g2 = resonance_par[name + "_g2"];
1699
1700 return Flatte2( mx2, mr2, g1, m1a, m1b, g2, m2a, m2b );
1701}
1702
1703void D0ToKSpipipi0::readInputCoeff() {
1704 coefficient["phasespace"] = 0;
1705 coefficient["<kstarm_892|rhop_770>_0_mag"] = 25;
1706 coefficient["<kstarm_892|rhop_770>_0_phase"] = 0;
1707 coefficient["<kstarm_892|rhop_770>_1_mag"] = 11.2665;
1708 coefficient["<kstarm_892|rhop_770>_1_phase"] = 0.548994;
1709 coefficient["<kstarm_892|rhop_770>_2_mag"] = 6.73954;
1710 coefficient["<kstarm_892|rhop_770>_2_phase"] = 1.96845;
1711 coefficient["<kstar0_892|rho0_770>_0_mag"] = 9.62259;
1712 coefficient["<kstar0_892|rho0_770>_0_phase"] = 4.39676;
1713 coefficient["<kstar0_892|rho0_770>_1_mag"] = 3.49377;
1714 coefficient["<kstar0_892|rho0_770>_1_phase"] = -2.56438;
1715 coefficient["<kstar0_892|rho0_770>_2_mag"] = 3.30796;
1716 coefficient["<kstar0_892|rho0_770>_2_phase"] = 5.17417;
1717 coefficient["<omega3_782|rho_770>_1_mag"] = 14.1209;
1718 coefficient["<omega3_782|rho_770>_1_phase"] = 2.06989;
1719 coefficient["<phi_1020|rho_770>_1_mag"] = 0.835627;
1720 coefficient["<phi_1020|rho_770>_1_phase"] = -2.92636;
1721 coefficient["<eta_547|rhom_phsp>_1_mag"] = 81.3879;
1722 coefficient["<eta_547|rhom_phsp>_1_phase"] = -1.51303;
1723 coefficient["<eta_547|rhom_phsp>_0_mag"] = 72.5434;
1724 coefficient["<eta_547|rhom_phsp>_0_phase"] = 3.89337;
1725 coefficient["<k1m_1270|rhom_770>_0_mag"] = 48.1147;
1726 coefficient["<k1m_1270|rhom_770>_0_phase"] = -0.0202615;
1727 coefficient["<k1m_1270|kstar_892>_0_mag"] = 5.73667;
1728 coefficient["<k1m_1270|kstar_892>_0_phase"] = 3.92441;
1729 coefficient["<k1m_1270|k0star_1430>_1_mag"] = 72.8988;
1730 coefficient["<k1m_1270|k0star_1430>_1_phase"] = -1.12332;
1731 coefficient["<k10_1270|rho0_770>_0_mag"] = 14.6425;
1732 coefficient["<k10_1270|rho0_770>_0_phase"] = 0.882522;
1733 coefficient["<k10_1270|kstarm_892>_0_mag"] = 4.40263;
1734 coefficient["<k10_1270|kstarm_892>_0_phase"] = 0.299025;
1735 coefficient["<k10_1270|k0starm_1430>_1_mag"] = 57.4271;
1736 coefficient["<k10_1270|k0starm_1430>_1_phase"] = 1.25342;
1737 coefficient["<k1m_1400|kstar_892>_0_mag"] = 7.41837;
1738 coefficient["<k1m_1400|kstar_892>_0_phase"] = 2.60505;
1739 coefficient["<k10_1400|kstarm_892>_0_mag"] = -38.728;
1740 coefficient["<k10_1400|kstarm_892>_0_phase"] = 2.64084;
1741 coefficient["<kstar02_1410|rho0_770>_1_mag"] = 32.6332;
1742 coefficient["<kstar02_1410|rho0_770>_1_phase"] = 0.266067;
1743 coefficient["<kstarm2_1410|rhop_770>_0_mag"] = 84.8711;
1744 coefficient["<kstarm2_1410|rhop_770>_0_phase"] = -1.22924;
1745 coefficient["<kstarm2_1410|rhop_770>_1_mag"] = 30.6599;
1746 coefficient["<kstarm2_1410|rhop_770>_1_phase"] = -3.02164;
1747 coefficient["<kstarm2_1410|rhop_770>_2_mag"] = 23.5133;
1748 coefficient["<kstarm2_1410|rhop_770>_2_phase"] = 4.90416;
1749 coefficient["<kstarm3_1410|kstar_892>_1_mag"] = -5.00895;
1750 coefficient["<kstarm3_1410|kstar_892>_1_phase"] = 0.75297;
1751 coefficient["<kstarm3_1410|rhom_770>_1_mag"] = 21.0146;
1752 coefficient["<kstarm3_1410|rhom_770>_1_phase"] = 2.11827;
1753 coefficient["<km_1460|kstar_892>_1_mag"] = -18.5953;
1754 coefficient["<km_1460|kstar_892>_1_phase"] = 1.72652;
1755 coefficient["<km_1460|k0star_1430>_0_mag"] = 500;
1756 coefficient["<km_1460|k0star_1430>_0_phase"] = 3.22602;
1757 coefficient["<k0_1460|kstarm_892>_1_mag"] = 14.7543;
1758 coefficient["<k0_1460|kstarm_892>_1_phase"] = 0.510956;
1759 coefficient["<k1m_1650|kstar_892>_0_mag"] = 11.7635;
1760 coefficient["<k1m_1650|kstar_892>_0_phase"] = 2.59972;
1761 coefficient["<k10_1650|kstarm_892>_0_mag"] = 27.9014;
1762 coefficient["<k10_1650|kstarm_892>_0_phase"] = -1.07364;
1763 coefficient["<k1m_1650|rhom_770>_0_mag"] = 23.2323;
1764 coefficient["<k1m_1650|rhom_770>_0_phase"] = -1.33583;
1765 coefficient["<kstar03_1680|kstarm_892>_1_mag"] = 23.4824;
1766 coefficient["<kstar03_1680|kstarm_892>_1_phase"] = 0.139691;
1767 coefficient["<kstarm3_1680|kstar_892>_1_mag"] = 11.8553;
1768 coefficient["<kstarm3_1680|kstar_892>_1_phase"] = -3.47852;
1769 coefficient["<k0star0_1430|rho0_phsp>_0_mag"] = 500;
1770 coefficient["<k0star0_1430|rho0_phsp>_0_phase"] = -5.1306;
1771 coefficient["<k1p_phsp|rhop_770>_1_mag"] = 197.192;
1772 coefficient["<k1p_phsp|rhop_770>_1_phase"] = 0.35067;
1773 coefficient["<k1m_phsp|rhom_770>_1_mag"] = 151.864;
1774 coefficient["<k1m_phsp|rhom_770>_1_phase"] = -5.11224;
1775 coefficient["<k1m_phsp|rhom_770>_2_mag"] = 21.6383;
1776 coefficient["<k1m_phsp|rhom_770>_2_phase"] = -0.843621;
1777 coefficient["<k0starm_1430|rhop_phsp>_0_mag"] = -500;
1778 coefficient["<k0starm_1430|rhop_phsp>_0_phase"] = 4.00277;
1779 coefficient["<k0_1460|sigma_500>_0_mag"] = 500;
1780 coefficient["<k0_1460|sigma_500>_0_phase"] = -5.70929;
1781 coefficient["<k10_1270|f0_980>_1_mag"] = 106.293;
1782 coefficient["<k10_1270|f0_980>_1_phase"] = 0.947063;
1783 coefficient["<k10_1400|f0_1370>_1_mag"] = 184.717;
1784 coefficient["<k10_1400|f0_1370>_1_phase"] = -0.0197254;
1785 coefficient["<rho0_770|k0star0_1430>_1_mag"] = 159.521;
1786 coefficient["<rho0_770|k0star0_1430>_1_phase"] = 0.164038;
1787 coefficient["<k2starm3_1430|kstar_892>_2_mag"] = 0.400134;
1788 coefficient["<k2starm3_1430|kstar_892>_2_phase"] = 0.533053;
1789 coefficient["<k2star02_1430|f2_1270>_1_mag"] = 45.2478;
1790 coefficient["<k2star02_1430|f2_1270>_1_phase"] = 1.81556;
1791 coefficient["<kstar0_892|f0_980>_1_mag"] = 24.2878;
1792 coefficient["<kstar0_892|f0_980>_1_phase"] = 1.93217;
1793 coefficient["<kstarm_phsp|rhop_phsp>_0_mag"] = 500;
1794 coefficient["<kstarm_phsp|rhop_phsp>_0_phase"] = 1.48706;
1795 coefficient["<rho0_770|kstar0_phsp>_1_mag"] = 58.8657;
1796 coefficient["<rho0_770|kstar0_phsp>_1_phase"] = 3.54134;
1797}
1798
1799complex<double> D0ToKSpipipi0::Amp( vector<double> ks0, vector<double> pip, vector<double> pim,
1800 vector<double> pi0 ) {
1801
1802 initPar();
1803
1804 // (E; PX, PY, PZ)
1805 // EvtVector4R ks0 = _pd[0];
1806 // EvtVector4R pip = _pd[1];
1807 // EvtVector4R pim = _pd[2];
1808 // EvtVector4R pi0 = _pd[3];
1809
1810 createPropagator( ks0, pip, pim, pi0 );
1811 createSpinfactor( ks0, pip, pim, pi0 );
1812
1813 //----------------------------------------------------------------------//
1814 // Add the partial waves here //
1815 //----------------------------------------------------------------------//
1816
1817 double f_isoConj = -1. * sqrt( 2.0 );
1818 double one = 1.0, minus_one = -1.0;
1819 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_0"], propagator["kstarm_892"],
1820 propagator["rhop_770"], coefficient["<kstarm_892|rhop_770>_0_mag"],
1821 coefficient["<kstarm_892|rhop_770>_0_phase"] );
1822 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_1"], propagator["kstarm_892"],
1823 propagator["rhop_770"], coefficient["<kstarm_892|rhop_770>_1_mag"],
1824 coefficient["<kstarm_892|rhop_770>_1_phase"] );
1825 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_2"], propagator["kstarm_892"],
1826 propagator["rhop_770"], coefficient["<kstarm_892|rhop_770>_2_mag"],
1827 coefficient["<kstarm_892|rhop_770>_2_phase"] );
1828 addPartialWave( spinfactor["D2VV_ks0pi0_pippim_0"], propagator["kstar0_892"],
1829 propagator["rho0_770"], coefficient["<kstar0_892|rho0_770>_0_mag"],
1830 coefficient["<kstar0_892|rho0_770>_0_phase"] );
1831 addPartialWave( spinfactor["D2VV_ks0pi0_pippim_1"], propagator["kstar0_892"],
1832 propagator["rho0_770"], coefficient["<kstar0_892|rho0_770>_1_mag"],
1833 coefficient["<kstar0_892|rho0_770>_1_phase"] );
1834 addPartialWave( spinfactor["D2VV_ks0pi0_pippim_2"], propagator["kstar0_892"],
1835 propagator["rho0_770"], coefficient["<kstar0_892|rho0_770>_2_mag"],
1836 coefficient["<kstar0_892|rho0_770>_2_phase"] );
1837 complex<double> omega_ks0_propagator =
1838 propagator["rhop_770"] * spinfactor["D2VP_V2VP_pippimpi0_pippi0_1"] +
1839 minus_one * propagator["rho0_770"] * spinfactor["D2VP_V2VP_pippimpi0_pippim_1"] +
1840 propagator["rhom_770"] * spinfactor["D2VP_V2VP_pippimpi0_pimpi0_1"];
1841 addPartialWave( omega_ks0_propagator, propagator["omega3_782"],
1842 coefficient["<omega3_782|rho_770>_1_mag"],
1843 coefficient["<omega3_782|rho_770>_1_phase"] );
1844 complex<double> phi_ks0_propagator =
1845 propagator["rhop_770"] * spinfactor["D2VP_V2VP_pippimpi0_pippi0_1"] +
1846 minus_one * propagator["rho0_770"] * spinfactor["D2VP_V2VP_pippimpi0_pippim_1"] +
1847 propagator["rhom_770"] * spinfactor["D2VP_V2VP_pippimpi0_pimpi0_1"];
1848 addPartialWave( phi_ks0_propagator, propagator["phi_1020"],
1849 coefficient["<phi_1020|rho_770>_1_mag"],
1850 coefficient["<phi_1020|rho_770>_1_phase"] );
1851 addPartialWave( spinfactor["D2PP_P2VP_pippimpi0_pimpi0_1"], propagator["eta_547"],
1852 propagator["rhom_phsp"], coefficient["<eta_547|rhom_phsp>_1_mag"],
1853 coefficient["<eta_547|rhom_phsp>_1_phase"] );
1854 addPartialWave( spinfactor["D2PP_P2SP_pippimpi0_pimpi0_0"], propagator["eta_547"],
1855 propagator["rhom_phsp"], coefficient["<eta_547|rhom_phsp>_0_mag"],
1856 coefficient["<eta_547|rhom_phsp>_0_phase"] );
1857 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_pimpi0_0"], propagator["k1m_1270"],
1858 propagator["rhom_770"], coefficient["<k1m_1270|rhom_770>_0_mag"],
1859 coefficient["<k1m_1270|rhom_770>_0_phase"] );
1860 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_ks0pim_0"], propagator["k1m_1270"],
1861 propagator["kstarm_892"], spinfactor["D2AP_A2VP_ks0pimpi0_ks0pi0_0"],
1862 propagator["k1m_1270"], propagator["kstar0_892"], f_isoConj,
1863 coefficient["<k1m_1270|kstar_892>_0_mag"],
1864 coefficient["<k1m_1270|kstar_892>_0_phase"] );
1865 addPartialWave( spinfactor["D2AP_A2SP_ks0pimpi0_ks0pim_1"], propagator["k1m_1270"],
1866 propagator["k0starm_1430"], spinfactor["D2AP_A2SP_ks0pimpi0_ks0pi0_1"],
1867 propagator["k1m_1270"], propagator["k0star0_1430"], f_isoConj,
1868 coefficient["<k1m_1270|k0star_1430>_1_mag"],
1869 coefficient["<k1m_1270|k0star_1430>_1_phase"] );
1870 addPartialWave( spinfactor["D2AP_A2VP_ks0pippim_pippim_0"], propagator["k10_1270"],
1871 propagator["rho0_770"], coefficient["<k10_1270|rho0_770>_0_mag"],
1872 coefficient["<k10_1270|rho0_770>_0_phase"] );
1873 addPartialWave( spinfactor["D2AP_A2VP_ks0pippim_ks0pim_0"], propagator["k10_1270"],
1874 propagator["kstarm_892"], coefficient["<k10_1270|kstarm_892>_0_mag"],
1875 coefficient["<k10_1270|kstarm_892>_0_phase"] );
1876 addPartialWave( spinfactor["D2AP_A2SP_ks0pippim_ks0pim_1"], propagator["k10_1270"],
1877 propagator["k0starm_1430"], coefficient["<k10_1270|k0starm_1430>_1_mag"],
1878 coefficient["<k10_1270|k0starm_1430>_1_phase"] );
1879 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_ks0pim_0"], propagator["k1m_1400"],
1880 propagator["kstarm_892"], spinfactor["D2AP_A2VP_ks0pimpi0_ks0pi0_0"],
1881 propagator["k1m_1400"], propagator["kstar0_892"], f_isoConj,
1882 coefficient["<k1m_1400|kstar_892>_0_mag"],
1883 coefficient["<k1m_1400|kstar_892>_0_phase"] );
1884 addPartialWave( spinfactor["D2AP_A2VP_ks0pippim_ks0pim_0"], propagator["k10_1400"],
1885 propagator["kstarm_892"], coefficient["<k10_1400|kstarm_892>_0_mag"],
1886 coefficient["<k10_1400|kstarm_892>_0_phase"] );
1887 addPartialWave( spinfactor["D2VV_ks0pi0_pippim_1"], propagator["kstar02_1410"],
1888 propagator["rho0_770"], coefficient["<kstar02_1410|rho0_770>_1_mag"],
1889 coefficient["<kstar02_1410|rho0_770>_1_phase"] );
1890 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_0"], propagator["kstarm2_1410"],
1891 propagator["rhop_770"], coefficient["<kstarm2_1410|rhop_770>_0_mag"],
1892 coefficient["<kstarm2_1410|rhop_770>_0_phase"] );
1893 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_1"], propagator["kstarm2_1410"],
1894 propagator["rhop_770"], coefficient["<kstarm2_1410|rhop_770>_1_mag"],
1895 coefficient["<kstarm2_1410|rhop_770>_1_phase"] );
1896 addPartialWave( spinfactor["D2VV_ks0pim_pippi0_2"], propagator["kstarm2_1410"],
1897 propagator["rhop_770"], coefficient["<kstarm2_1410|rhop_770>_2_mag"],
1898 coefficient["<kstarm2_1410|rhop_770>_2_phase"] );
1899 addPartialWave( spinfactor["D2VP_V2VP_ks0pimpi0_ks0pim_1"], propagator["kstarm3_1410"],
1900 propagator["kstarm_892"], spinfactor["D2VP_V2VP_ks0pimpi0_ks0pi0_1"],
1901 propagator["kstarm3_1410"], propagator["kstar0_892"], f_isoConj,
1902 coefficient["<kstarm3_1410|kstar_892>_1_mag"],
1903 coefficient["<kstarm3_1410|kstar_892>_1_phase"] );
1904 addPartialWave( spinfactor["D2VP_V2VP_ks0pimpi0_pimpi0_1"], propagator["kstarm3_1410"],
1905 propagator["rhom_770"], coefficient["<kstarm3_1410|rhom_770>_1_mag"],
1906 coefficient["<kstarm3_1410|rhom_770>_1_phase"] );
1907 addPartialWave( spinfactor["D2PP_P2VP_ks0pimpi0_ks0pim_1"], propagator["km_1460"],
1908 propagator["kstarm_892"], spinfactor["D2PP_P2VP_ks0pimpi0_ks0pi0_1"],
1909 propagator["km_1460"], propagator["kstar0_892"], f_isoConj,
1910 coefficient["<km_1460|kstar_892>_1_mag"],
1911 coefficient["<km_1460|kstar_892>_1_phase"] );
1912 addPartialWave( spinfactor["D2PP_P2SP_ks0pimpi0_ks0pim_0"], propagator["km_1460"],
1913 propagator["k0starm_1430"], spinfactor["D2PP_P2SP_ks0pimpi0_ks0pi0_0"],
1914 propagator["km_1460"], propagator["k0star0_1430"], f_isoConj,
1915 coefficient["<km_1460|k0star_1430>_0_mag"],
1916 coefficient["<km_1460|k0star_1430>_0_phase"] );
1917 addPartialWave( spinfactor["D2PP_P2VP_ks0pippim_ks0pim_1"], propagator["k0_1460"],
1918 propagator["kstarm_892"], coefficient["<k0_1460|kstarm_892>_1_mag"],
1919 coefficient["<k0_1460|kstarm_892>_1_phase"] );
1920 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_ks0pim_0"], propagator["k1m_1650"],
1921 propagator["kstarm_892"], spinfactor["D2AP_A2VP_ks0pimpi0_ks0pi0_0"],
1922 propagator["k1m_1650"], propagator["kstar0_892"], f_isoConj,
1923 coefficient["<k1m_1650|kstar_892>_0_mag"],
1924 coefficient["<k1m_1650|kstar_892>_0_phase"] );
1925 addPartialWave( spinfactor["D2AP_A2VP_ks0pippim_ks0pim_0"], propagator["k10_1650"],
1926 propagator["kstarm_892"], coefficient["<k10_1650|kstarm_892>_0_mag"],
1927 coefficient["<k10_1650|kstarm_892>_0_phase"] );
1928 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_pimpi0_0"], propagator["k1m_1650"],
1929 propagator["rhom_770"], coefficient["<k1m_1650|rhom_770>_0_mag"],
1930 coefficient["<k1m_1650|rhom_770>_0_phase"] );
1931 addPartialWave( spinfactor["D2VP_V2VP_ks0pippim_ks0pim_1"], propagator["kstar03_1680"],
1932 propagator["kstarm_892"], coefficient["<kstar03_1680|kstarm_892>_1_mag"],
1933 coefficient["<kstar03_1680|kstarm_892>_1_phase"] );
1934 addPartialWave( spinfactor["D2VP_V2VP_ks0pimpi0_ks0pim_1"], propagator["kstarm3_1680"],
1935 propagator["kstarm_892"], spinfactor["D2VP_V2VP_ks0pimpi0_ks0pi0_1"],
1936 propagator["kstarm3_1680"], propagator["kstar0_892"], f_isoConj,
1937 coefficient["<kstarm3_1680|kstar_892>_1_mag"],
1938 coefficient["<kstarm3_1680|kstar_892>_1_phase"] );
1939 addPartialWave( spinfactor["D2SS_ks0pi0_pippim_0"], propagator["k0star0_1430"],
1940 propagator["rho0_phsp"], coefficient["<k0star0_1430|rho0_phsp>_0_mag"],
1941 coefficient["<k0star0_1430|rho0_phsp>_0_phase"] );
1942 addPartialWave( spinfactor["D2PP_P2VP_ks0pippi0_pippi0_1"], propagator["k1p_phsp"],
1943 propagator["rhop_770"], coefficient["<k1p_phsp|rhop_770>_1_mag"],
1944 coefficient["<k1p_phsp|rhop_770>_1_phase"] );
1945 addPartialWave( spinfactor["D2PP_P2VP_ks0pimpi0_pimpi0_1"], propagator["k1m_phsp"],
1946 propagator["rhom_770"], coefficient["<k1m_phsp|rhom_770>_1_mag"],
1947 coefficient["<k1m_phsp|rhom_770>_1_phase"] );
1948 addPartialWave( spinfactor["D2AP_A2VP_ks0pimpi0_pimpi0_2"], propagator["k1m_phsp"],
1949 propagator["rhom_770"], coefficient["<k1m_phsp|rhom_770>_2_mag"],
1950 coefficient["<k1m_phsp|rhom_770>_2_phase"] );
1951 addPartialWave( spinfactor["D2SS_ks0pim_pippi0_0"], propagator["k0starm_1430"],
1952 propagator["rhop_phsp"], coefficient["<k0starm_1430|rhop_phsp>_0_mag"],
1953 coefficient["<k0starm_1430|rhop_phsp>_0_phase"] );
1954 addPartialWave( spinfactor["D2PP_P2SP_ks0pippim_pippim_0"], propagator["k0_1460"],
1955 propagator["sigma_500"], coefficient["<k0_1460|sigma_500>_0_mag"],
1956 coefficient["<k0_1460|sigma_500>_0_phase"] );
1957 addPartialWave( spinfactor["D2AP_A2SP_ks0pippim_pippim_1"], propagator["k10_1270"],
1958 propagator["f0_980"], coefficient["<k10_1270|f0_980>_1_mag"],
1959 coefficient["<k10_1270|f0_980>_1_phase"] );
1960 addPartialWave( spinfactor["D2AP_A2SP_ks0pippim_pippim_1"], propagator["k10_1400"],
1961 propagator["f0_1370"], coefficient["<k10_1400|f0_1370>_1_mag"],
1962 coefficient["<k10_1400|f0_1370>_1_phase"] );
1963 addPartialWave( spinfactor["D2VS_pippim_ks0pi0_1"], propagator["rho0_770"],
1964 propagator["k0star0_1430"], coefficient["<rho0_770|k0star0_1430>_1_mag"],
1965 coefficient["<rho0_770|k0star0_1430>_1_phase"] );
1966 addPartialWave( spinfactor["D2TP_T2VP_ks0pimpi0_ks0pim_2"], propagator["k2starm3_1430"],
1967 propagator["kstarm_892"], spinfactor["D2TP_T2VP_ks0pimpi0_ks0pi0_2"],
1968 propagator["k2starm3_1430"], propagator["kstar0_892"], f_isoConj,
1969 coefficient["<k2starm3_1430|kstar_892>_2_mag"],
1970 coefficient["<k2starm3_1430|kstar_892>_2_phase"] );
1971 addPartialWave( spinfactor["D2TT_ks0pi0_pippim_1"], propagator["k2star02_1430"],
1972 propagator["f2_1270"], coefficient["<k2star02_1430|f2_1270>_1_mag"],
1973 coefficient["<k2star02_1430|f2_1270>_1_phase"] );
1974 addPartialWave( spinfactor["D2VS_ks0pi0_pippim_1"], propagator["kstar0_892"],
1975 propagator["f0_980"], coefficient["<kstar0_892|f0_980>_1_mag"],
1976 coefficient["<kstar0_892|f0_980>_1_phase"] );
1977 addPartialWave( spinfactor["D2SS_ks0pim_pippi0_0"], propagator["kstarm_phsp"],
1978 propagator["rhop_phsp"], coefficient["<kstarm_phsp|rhop_phsp>_0_mag"],
1979 coefficient["<kstarm_phsp|rhop_phsp>_0_phase"] );
1980 addPartialWave( spinfactor["D2VS_pippim_ks0pi0_1"], propagator["rho0_770"],
1981 propagator["kstar0_phsp"], coefficient["<rho0_770|kstar0_phsp>_1_mag"],
1982 coefficient["<rho0_770|kstar0_phsp>_1_phase"] );
1983
1984 // double Amps = abs2(totalAmp);
1985 // return Amps;
1986 return totalAmp;
1987}
double p2[4]
double p1[4]
double P(RecMdcKalTrack *trk)
double mass
std::vector< double > operator+(std::vector< double > &lhs, std::vector< double > &rhs)
string toString(const T &t)
TF1 * g1
double imag(const EvtComplex &c)
EvtComplex exp(const EvtComplex &c)
const double mass_Pion
double mpi
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per saves r n generator level $ !Flag for chat level in output
Definition FoamA.h:89
XmlRpcServer s
const DifNumber one
****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
complex< double > Amp(vector< double > ks0, vector< double > pip, vector< double > pim, vector< double > pi0)
virtual ~D0ToKSpipipi0()
int t()
Definition t.c:1