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