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