BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDsToKSKmPipPip.cc
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2// Environment:
3// This software is part of the EvtGen package developed jointly
4// for the BaBar and CLEO collaborations. If you use all or part
5// of it, please give an appropriate acknowledgement.
6//
7// Copyright Information: See EvtGen/COPYRIGHT
8// Copyright (C) 1998 Caltech, UCSB
9//
10// Module: EvtDsToKSKmPipPip.cc
11//
12// Description: Routine to handle three-body decays of Ds+/Ds-
13//
14// Modification history:
15// Liaoyuan Dong Sun Jun 02 01:23:25 2024 Module created
16//------------------------------------------------------------------------
18#include <fstream>
19#include <stdlib.h>
22#include "EvtGenBase/EvtPDL.hh"
27#include <string>
31using std::endl;
32
34
35void EvtDsToKSKmPipPip::getName(std::string& model_name){
36 model_name="DsToKSKmPipPip";
37}
38
42
44 // check that there are 0 arguments
45 checkNArg(0);
46 checkNDaug(4);
47
53
54 int mode = 0;
55 phi[1] = -1.5871e+00;
56 phi[2] = -6.4308e+00;
57 phi[3] = -4.4006e+00;
58 phi[4] = 4.7540e+00;
59 phi[5] = 4.2785e+00;
60 phi[6] = 6.3523e+00;
61 phi[7] = 6.3523e+00;
62 phi[8] = 5.4902e+00;
63 phi[9] = 5.2024e+00;
64 phi[10]= -3.9776e+00;
65
66 rho[1] = 3.6108e-01;
67 rho[2] = 4.5990e-01;
68 rho[3] = 1.7722e+00;
69 rho[4] = 2.1434e+00;
70 rho[5] = 1.2504e+00;
71 rho[6] = 1.0069e+00;
72 rho[7] = 1.0069e+00;
73 rho[8] = 1.7351e+00;
74 rho[9] = 4.1005e+00;
75 rho[10]= 3.7821e+00;
76
77 modetype[0] = 1;
78 modetype[1] = 2;
79 modetype[2] = 3;
80 modetype[3] = 21;
81 modetype[4] = 22;
82 modetype[5] = 8;
83 modetype[6] = 4;
84 modetype[7] = 14;
85 modetype[8] = 17;
86 modetype[9] = 5;
87 modetype[10]= 20;
88
89 phi[0] = 0;
90 rho[0] = 1;
91
92 //cout << "DsToKSKmPipPip :" << endl;
93 //for (int i=0; i<11; i++) {
94 // cout << i << " rho= " << rho[i] << " phi= " << phi[i] << endl;
95 //}
96
97 G1510 = 1.0000e+00;
98 GKst0 = 4.7300e-02;
99 GKstp = 5.0300e-02;
100 Ga0_980 = 5.0000e-02;
101 Geta1475 = 9.0000e-02;
102 Gf1285 = 2.2700e-02;
103 m1510 = 1.6200e+00;
104 mKst0 = 8.9555e-01;
105 mKstp = 8.9176e-01;
106 ma0_980 = 9.9000e-01;
107 meta1475 = 1.4750e+00;
108 mf1285 = 1.2819e+00;
109
110 mD = 1.86486;
111 metap = 0.95778;
112 mk0 = 0.497614;
113 mass_Kaon = 0.49368;
114 mass_Eta = 0.547862;
115 mass_Pion = 0.13957;
116 math_pi = 3.1415926;
117 mass_Pion2 = 0.0194797849;
118 mass_2Pion = 0.27914;
119 math_2pi = 6.2831852;
120 rD2 = 25.0; // 5*5
121 rRes2 = 9.0; // 3*3
122 g1 = 0.5468;
123 g2 = 0.23;
124 GS1 = 0.636619783;
125 GS2 = 0.01860182466;
126 GS3 = 0.1591549458; // 1/(2*math_2pi)
127 GS4 = 0.00620060822; // mass_Pion2/math_pi
128
129 int GG[4][4] = { {1,0,0,0}, {0,-1,0,0}, {0,0,-1,0}, {0,0,0,-1} };
130 int EE[4][4][4][4] =
131 { { {{0,0,0,0}, {0,0,0,0}, {0,0,0,0}, {0,0,0,0} },
132 {{0,0,0,0}, {0,0,0,0}, {0,0,0,1}, {0,0,-1,0}},
133 {{0,0,0,0}, {0,0,0,-1}, {0,0,0,0}, {0,1,0,0} },
134 {{0,0,0,0}, {0,0,1,0}, {0,-1,0,0}, {0,0,0,0} } },
135 { {{0,0,0,0}, {0,0,0,0}, {0,0,0,-1}, {0,0,1,0} },
136 {{0,0,0,0}, {0,0,0,0}, {0,0,0,0}, {0,0,0,0} },
137 {{0,0,0,1}, {0,0,0,0}, {0,0,0,0}, {-1,0,0,0}},
138 {{0,0,-1,0}, {0,0,0,0}, {1,0,0,0}, {0,0,0,0} } },
139 { {{0,0,0,0}, {0,0,0,1}, {0,0,0,0}, {0,-1,0,0}},
140 {{0,0,0,-1}, {0,0,0,0}, {0,0,0,0}, {1,0,0,0} },
141 {{0,0,0,0}, {0,0,0,0}, {0,0,0,0}, {0,0,0,0} },
142 {{0,1,0,0}, {-1,0,0,0}, {0,0,0,0}, {0,0,0,0} } },
143 { {{0,0,0,0}, {0,0,-1,0}, {0,1,0,0}, {0,0,0,0} },
144 {{0,0,1,0}, {0,0,0,0}, {-1,0,0,0}, {0,0,0,0} },
145 {{0,-1,0,0}, {1,0,0,0}, {0,0,0,0}, {0,0,0,0} },
146 {{0,0,0,0}, {0,0,0,0}, {0,0,0,0}, {0,0,0,0} } } };
147 for (int i=0; i<4; i++) {
148 for (int j=0; j<4; j++) {
149 G[i][j] = GG[i][j];
150 for (int k=0; k<4; k++) {
151 for (int l=0; l<4; l++) {
152 E[i][j][k][l] = EE[i][j][k][l];
153 }
154 }
155 }
156 }
157}
158
159
160
161
163 //setProbMax(75000);
164 setProbMax(95050);
165 //setProbMax(1.0);
166}
167
169/*
170 double maxprob=0;
171 for(int ir=0;ir<=60000000;ir++){
172 p->initializePhaseSpace(getNDaug(),getDaugs());
173 EvtVector4R _ks = p->getDaug(0)->getP4();
174 EvtVector4R _km = p->getDaug(1)->getP4();
175 EvtVector4R _pip = p->getDaug(2)->getP4();
176 EvtVector4R _pip2 = p->getDaug(3)->getP4();
177
178 double _Pip[4],_Pip2[4],_Ks[4],_Km[4];
179 _Ks[0] = _ks.get(0); _Km[0] = _km.get(0); _Pip[0] = _pip.get(0); _Pip2[0] = _pip2.get(0);
180 _Ks[1] = _ks.get(1); _Km[1] = _km.get(1); _Pip[1] = _pip.get(1); _Pip2[1] = _pip2.get(1);
181 _Ks[2] = _ks.get(2); _Km[2] = _km.get(2); _Pip[2] = _pip.get(2); _Pip2[2] = _pip2.get(2);
182 _Ks[3] = _ks.get(3); _Km[3] = _km.get(3); _Pip[3] = _pip.get(3); _Pip2[3] = _pip2.get(3);
183
184 double _prob;
185 int g0[11]={1,1,1,1,1,1,1,1,1,1,1};
186 int nstates=11;
187
188 calEvaMy(_Ks, _Km, _Pip, _Pip2, mass, width, rho, phi, g0, modetype, nstates, _prob);
189 if(_prob>maxprob) {
190 maxprob=_prob;
191 printf("ir = %d,maxprob = %.10f,prob = %.10f\n\n",ir,maxprob,_prob);
192 }
193 }
194 printf("maxprob = %.10f\n", maxprob);*/
195
197 EvtVector4R ks = p->getDaug(0)->getP4();
198 EvtVector4R km = p->getDaug(1)->getP4();
199 EvtVector4R pip = p->getDaug(2)->getP4();
200 EvtVector4R pip2= p->getDaug(3)->getP4();
201
202 double Pip[4],Pip2[4],Km[4],Ks[4];
203 Ks[0] = ks.get(0); Km[0] = km.get(0); Pip[0] = pip.get(0); Pip2[0] = pip2.get(0);
204 Ks[1] = ks.get(1); Km[1] = km.get(1); Pip[1] = pip.get(1); Pip2[1] = pip2.get(1);
205 Ks[2] = ks.get(2); Km[2] = km.get(2); Pip[2] = pip.get(2); Pip2[2] = pip2.get(2);
206 Ks[3] = ks.get(3); Km[3] = km.get(3); Pip[3] = pip.get(3); Pip2[3] = pip2.get(3);
207
208 //Ks[0] = 0.6096125554; Ks[1] = 0.2036145752; Ks[2] = -0.2256821642; Ks[3] = 0.1884652565 ;
209 //Km[0] = 0.7106732604; Km[1] = -0.2599714599; Km[2] = 0.3643877399; Km[3] = 0.2469330230 ;
210 //Pip[0] = 0.4177584141; Pip[1] = -0.0562004660; Pip[2] = 0.2821155804; Pip[3] =-0.0020089709 ;
211 //Pip2[0] = 0.3945120988; Pip2[1] = -0.1225632042; Pip2[2] = -0.3645928242; Pip2[3] = 0.0521817527 ;
212
213 //Ks[0] = 0.6439438863; Ks[1] = 0.3368405824; Ks[2] = 0.2308441601; Ks[3] = 0.0171298406;
214 //Km[0] = 0.7146391007; Km[1] = 0.4061832461; Km[2] = -0.1841669581; Km[3] = 0.2609401579;
215 //Pip[0] = 0.2817222681; Pip[1] = -0.1005161765;Pip[2] = -0.1172908044;Pip[3] = -0.1898077097;
216 //Pip2[0] = 0.3638598699;Pip2[1] = -0.3106229764;Pip2[2] = 0.0905020525;Pip2[3] = 0.0907574503;
217
218 double value;
219 int g0[11]={1,1,1,1,1,1,1,1,1,1,1};
220 int nstates=11;
221
222 calEvaMy(Ks,Km,Pip,Pip2,mass,width,rho,phi,g0,modetype,nstates,value);
223 setProb(value);
224
225 return;
226}
227
228void EvtDsToKSKmPipPip::Com_Multi(double a1[2], double a2[2], double res[2]){
229 res[0] = a1[0]*a2[0]-a1[1]*a2[1];
230 res[1] = a1[1]*a2[0]+a1[0]*a2[1];
231}
232void EvtDsToKSKmPipPip::Com_Divide(double a1[2], double a2[2], double res[2]){
233 res[0] = (a1[0]*a2[0]+a1[1]*a2[1])/(a2[0]*a2[0]+a2[1]*a2[1]);
234 res[1] = (a1[1]*a2[0]-a1[0]*a2[1])/(a2[0]*a2[0]+a2[1]*a2[1]);
235}
236double EvtDsToKSKmPipPip::SCADot(double a1[4], double a2[4])
237{
238 double _cal = 0.;
239 _cal = a1[0]*a2[0]-a1[1]*a2[1]-a1[2]*a2[2]-a1[3]*a2[3];
240 return _cal;
241}
242double EvtDsToKSKmPipPip::barrier(int l, double sa, double sb, double sc, double r2)
243{
244 double F;
245 double tmp = sa+sb-sc;
246 double q = 0.25*tmp*tmp/sa-sb;
247 if (q < 0) q = 1e-16;
248 double z = q*r2;
249 if (l == 1) {
250 F = sqrt(2.0*z/(1.0+z));
251 }
252 else if (l == 2) {
253 double z2 = z*z;
254 F = sqrt(13.0*z2/(9.0+3.0*z+z2));
255 } else {
256 F = 1.0;
257 }
258 return F;
259}
260void EvtDsToKSKmPipPip::calt1(double daug1[4], double daug0[4], double t1[4]){
261 double p, pq;
262 double pa[4], qa[4];
263 for(int i=0; i<4; i++){
264 pa[i] = daug1[i] + daug0[i];
265 qa[i] = daug1[i] - daug0[i];
266 }
267 p = SCADot(pa,pa);
268 pq = SCADot(pa,qa);
269 for(int i=0; i<4; i++){
270 t1[i] = qa[i] - pq/p*pa[i];
271 }
272}
273void EvtDsToKSKmPipPip::calt2(double daug1[4], double daug0[4], double t2[4][4]){
274 double p, r;
275 double pa[4], t1[4];
276 calt1(daug1,daug0,t1);
277 r = SCADot(t1,t1)/3.0;
278 for(int i=0; i<4; i++) {
279 pa[i] = daug1[i] + daug0[i];
280 }
281 p = SCADot(pa,pa);
282 for(int i=0; i<4; i++) {
283 for(int j=0; j<4; j++) {
284 t2[i][j] = t1[i]*t1[j] - r*(G[i][j]-pa[i]*pa[j]/p);
285 }
286 }
287}
288
289double EvtDsToKSKmPipPip::wid(double mass, double sa, double sb, double sc, double r2, int l){
290 double widm = 0.;
291 double sa0 = mass*mass;
292 double m = sqrt(sa);
293 double tmp = sb-sc;
294 double tmp1 = sa+tmp;
295 double q = 0.25*tmp1*tmp1/sa-sb;
296 if(q<0) q = 1e-16;
297 double tmp2 = sa0+tmp;
298 double q0 = 0.25*tmp2*tmp2/sa0-sb;
299 if(q0<0) q0 = 1e-16;
300 double z = q*r2;
301 double z0 = q0*r2;
302 double t = q/q0;
303 if(l == 0) {widm = sqrt(t)*mass/m;}
304 else if(l == 1) {widm = t*sqrt(t)*mass/m*(1+z0)/(1+z);}
305 else if(l == 2) {widm = t*t*sqrt(t)*mass/m*(9+3*z0+z0*z0)/(9+3*z+z*z);}
306 return widm;
307}
308void EvtDsToKSKmPipPip::Flatte_rhoab(double sa, double sb, double sc, double rho[2])
309{
310 double q = (sa+sb-sc)*(sa+sb-sc)/(4*sa)-sb;
311 if(q>0) {
312 rho[0]=2* sqrt(q/sa);
313 rho[1]=0;
314 }
315 else if(q<0){
316 rho[0]=0;
317 rho[1]=2*sqrt(-q/sa);
318 }
319}
320void EvtDsToKSKmPipPip::propagatora0980(double mass, double sx, double *sb, double *sc, double prop[2])
321{
322 double unit[2]={1.0};
323 double ci[2]={0,1};
324 double rho1[2];
325 Flatte_rhoab(sx,sb[0],sc[0],rho1);
326 double rho2[2];
327 Flatte_rhoab(sx,sb[1],sc[1],rho2);
328
329 double gKK_a980=0.892*0.341, gPiEta_a980=0.341;
330 double tmp1[2]={gKK_a980,0};
331 double tmp11[2];
332 double tmp2[2]={gPiEta_a980,0};
333 double tmp22[2];
334 Com_Multi(tmp1,rho1,tmp11);
335 Com_Multi(tmp2,rho2,tmp22);
336 double tmp3[2]={tmp11[0]+tmp22[0],tmp11[1]+tmp22[1]};
337 double tmp31[2];
338 Com_Multi(tmp3, ci,tmp31);
339 double tmp4[2]={mass*mass-sx-tmp31[0], -1.0*tmp31[1]};
340 Com_Divide( unit,tmp4, prop);
341}
342void EvtDsToKSKmPipPip::propagatorRBW(double mass, double width, double sa, double sb, double sc, double r2, int l, double prop[2]){
343 double a[2], b[2];
344 a[0] = 1;
345 a[1] = 0;
346 b[0] = mass*mass-sa;
347 b[1] = -mass*width*wid(mass,sa,sb,sc,r2,l);
348 Com_Divide(a,b,prop);
349}
350void EvtDsToKSKmPipPip::KPiSLASS(double sa, double sb, double sc, double prop[2]) {
351 const double m1430 = 1.441;
352 const double sa0 = 2.076481; // m1410*m1410;
353 const double w1430 = 0.193;
354 const double Lass1 = 0.25/sa0;
355 double tmp = sb-sc;
356 double tmp1 = sa0+tmp;
357 double q0 = Lass1*tmp1*tmp1-sb;
358 if(q0<0) q0 = 1e-16;
359 double tmp2 = sa+tmp;
360 double qs = 0.25*tmp2*tmp2/sa-sb;
361 double q = sqrt(qs);
362 double width = w1430*q*m1430/sqrt(sa*q0);
363 double temp_R = atan(m1430*width/(sa0-sa));
364 if(temp_R<0) temp_R += math_pi;
365 double deltaR = -1.915 + temp_R; //fiR=-109.7
366 double temp_F = atan(0.226*q/(2.0-3.8194*qs)); // 2.0*0.113=0.226; 0.113*33.8=1.926; a=0.113; r=-33.8
367 if(temp_F<0) temp_F += math_pi;
368 double deltaF = 0.002 + temp_F; //fiF=0.1
369 double deltaS = deltaR + 2.0*deltaF; ///
370 double t1 = 0.96*sin(deltaF); //F= 0.96
371 double t2 = sin(deltaR);
372 double CF[2], CS[2];
373 CF[0] = cos(deltaF);
374 CF[1] = sin(deltaF);
375 CS[0] = cos(deltaS);
376 CS[1] = sin(deltaS);
377 prop[0] = t1*CF[0] + t2*CS[0];
378 prop[1] = t1*CF[1] + t2*CS[1];
379}
380
381void EvtDsToKSKmPipPip::calEvaMy(double* Ks,double* Km, double* Pip, double* Pip2, double *mass1, double *width1, double *amp, double *phase,int* g0, int* modetype, int nstates, double & Result)
382{
383 double pD[4],pKsPip1[4],pKmPip2[4],pKsPip2[4],pKmPip1[4],pKsKmPip1[4],pKsKmPip2[4],pKsKm[4],pKmPipPip[4];
384 for(int i=0;i!=4;i++)
385 {
386 pD[i] = Ks[i]+Pip[i]+Km[i]+Pip2[i];
387 pKsPip1[i] =Ks[i]+Pip[i];
388 pKsPip2[i] =Ks[i]+Pip2[i];
389 pKmPip1[i] =Km[i]+Pip[i];
390 pKmPip2[i] =Km[i]+Pip2[i];
391 pKsKm[i] = Km[i]+Ks[i];
392 pKsKmPip1[i] = pKsPip1[i] + Km[i];
393 pKsKmPip2[i] = pKsPip2[i] + Km[i];
394 pKmPipPip[i] = Km[i]+Pip[i]+Pip2[i];
395 }
396 double sPip2,sPip1,sKs,sKm,sKsPip1,sKmPip2,sD,sKsPip2,sKmPip1,sKsKmPip1,sKsKmPip2,sKsKm,sKmPipPip;
397 sD = SCADot(pD,pD);
398 sKs = SCADot(Ks,Ks);
399 sKm = SCADot(Km,Km);
400 sPip1= SCADot(Pip,Pip);
401 sPip2= SCADot(Pip2,Pip2);
402 sKsPip1 = SCADot(pKsPip1,pKsPip1);
403 sKsPip2 = SCADot(pKsPip2,pKsPip2);
404 sKmPip1 = SCADot(pKmPip1,pKmPip1);
405 sKmPip2 = SCADot(pKmPip2,pKmPip2);
406 sKsKmPip1 = SCADot(pKsKmPip1,pKsKmPip1);
407 sKsKmPip2 = SCADot(pKsKmPip2,pKsKmPip2);
408 sKmPipPip = SCADot(pKmPipPip,pKmPipPip);
409 sKsKm =SCADot(pKsKm,pKsKm);
410
411 double t1KsPip1[4],t1KmPip2[4],t1KsPip2[4],t1KmPip1[4],t2KsKmPip1_kspip[4][4],t2KsKmPip2_kspip[4][4],t2KsKmPip1_kmpip[4][4],t2KsKmPip2_kmpip[4][4],t1KsKmPip1_kspip[4],t1KsKmPip2_kspip[4],t1KsKmPip1_kmpip[4],t1KsKmPip2_kmpip[4],t1KsKmPip1_ksk[4],t1KsKmPip2_ksk[4];
412
413 calt1(Ks,Pip, t1KsPip1);
414 calt1(Ks,Pip2,t1KsPip2);
415 calt1(Km,Pip, t1KmPip1);
416 calt1(Km,Pip2,t1KmPip2);
417 calt1(pKsPip1,Km,t1KsKmPip1_kspip);
418 calt1(pKsPip2,Km,t1KsKmPip2_kspip);
419 calt1(pKmPip1,Ks,t1KsKmPip1_kmpip);
420 calt1(pKmPip2,Ks,t1KsKmPip2_kmpip);
421 calt1(pKsKm,Pip,t1KsKmPip1_ksk);
422 calt1(pKsKm,Pip2,t1KsKmPip2_ksk);
423
424 calt2(pKsPip1,Km,t2KsKmPip1_kspip);
425 calt2(pKsPip2,Km,t2KsKmPip2_kspip);
426 calt2(pKmPip1,Ks,t2KsKmPip1_kmpip);
427 calt2(pKmPip2,Ks,t2KsKmPip2_kmpip);
428
429 double cof[2],amp_tmp[2],amp_PDF[2], PDF[2];
430 amp_PDF[0] = 0;
431 amp_PDF[1] = 0;
432 PDF[0] = 0;
433 PDF[1] = 0;
434 double temp_PDF,tmp1,tmp2,tmp1_a,tmp2_a,temp_PDF_a,tt1,tt2,tt3,tt4;
435 double pro1b[2],pro2b[2],pro1b_a[2],pro2b_a[2],pro_ksk[2],pro3[2],pro4[2],pro1[2],pro2[2];
436 double pro_kspi1[2],pro_kspi2[2],pro_kpi1[2],pro_kpi2[2],pro_kpi2sw[2],pro_kpi1sw[2],pro_kspi2sw[2],pro_kspi1sw[2];
437 double ispro_kspi1=0,ispro_kspi2=0,ispro_kpi1=0,ispro_kpi2=0,ispro_kpi2sw=0,ispro_kpi1sw=0,ispro_kspi2sw=0,ispro_kspi1sw=0,ispro_ksk=0;
438 double barr_kskpi1_ksk=-1.0,barr_kskpi2_ksk=-1.0,barr_kskpi1_kspi=-1.0,barr_kskpi2_kspi=-1.0,barr_kskpi1_kpi=-1.0,barr_kskpi2_kpi=-1.0,barr_kspi1=-1.0,barr_kspi2=-1.0,barr_kpi1=-1.0,barr_kpi2=-1.0,barr_kspi1_kpi2=-1.0,barr_kspi2_kpi1=-1.0,barr_ds_kskpi1=-1.0,barr_ds_kskpi2=-1.0;
439
440 double t1D[4],t2D[4][4],B[3];
441 double sKs2[2]={sKs,mass_Pion*mass_Pion};
442 double sKm2[2]={sKm,mass_Eta*mass_Eta};
443
444 for(int i=0; i<nstates; i++){
445 amp_tmp[0] = 0;
446 amp_tmp[1] = 0;
447 tmp1 = 0;tmp2 = 0;temp_PDF = 0;tmp1_a = 0;tmp2_a = 0;temp_PDF_a = 0;
448
449 cof[0] = amp[i]*cos(phase[i]);
450 cof[1] = amp[i]*sin(phase[i]);
451
452 //printf("phase= %i, %.5f, %.5f\n", i, amp[i],phase[i]);
453
454 //Kst+ Kst0
455 if(modetype[i]==1)
456 {
457 if (ispro_kspi1==0){
458 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
459 ispro_kspi1=1;
460 }
461 if (ispro_kpi2==0) {
462 propagatorRBW(mKst0,GKst0,sKmPip2,sKm,sPip2,rRes2,1,pro_kpi2);
463 ispro_kpi2=1;
464 }
465 Com_Multi(pro_kspi1,pro_kpi2,pro1b);
466
467 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
468 if(barr_kpi2<0.0) barr_kpi2 = barrier(1,sKmPip2,sKm,sPip2,rRes2);
469 for(int i=0; i<4; i++)
470 {
471 temp_PDF += G[i][i]*t1KsPip1[i]*t1KmPip2[i];
472 }
473 tmp1 = barr_kspi1*barr_kpi2*temp_PDF;
474 //printf("%.15f, %.15f , %.15f, %.15f \n", pro_kspi1[0],pro_kpi2[0],pro_kspi1[1],pro_kpi2[1]);
475
476
477 temp_PDF=0;
478 if(ispro_kspi2==0){
479 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
480 ispro_kspi2=1;
481 }
482 if(ispro_kpi1==0){
483 propagatorRBW(mKst0,GKst0,sKmPip1,sKm,sPip1,rRes2,1,pro_kpi1);
484 ispro_kpi1=1;
485 }
486 Com_Multi(pro_kspi2,pro_kpi1,pro2b);
487
488 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
489 if(barr_kpi1<0.0) barr_kpi1 = barrier(1,sKmPip1,sKm,sPip1,rRes2);
490 for(int i=0; i<4; i++)
491 {
492 temp_PDF += G[i][i]*t1KsPip2[i]*t1KmPip1[i];
493 }
494 tmp2 =barr_kspi2*barr_kpi1*temp_PDF;
495 //printf("%.15f, %.15f , %.15f, %.15f \n", pro_kspi2[0],pro_kpi1[0],pro_kspi2[1],pro_kpi1[1]);
496
497 }
498
499 // kstp_kst0_p
500 if(modetype[i]==2)
501 {
502 if (ispro_kspi1==0){
503 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
504 ispro_kspi1=1;
505 }
506 if (ispro_kpi2==0) {
507 propagatorRBW(mKst0,GKst0,sKmPip2,sKm,sPip2,rRes2,1,pro_kpi2);
508 ispro_kpi2=1;
509 }
510 Com_Multi(pro_kspi1,pro_kpi2,pro1b);
511
512 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
513 if(barr_kpi2<0.0) barr_kpi2 = barrier(1,sKmPip2,sKm,sPip2,rRes2);
514 if(barr_kspi1_kpi2<0.0) barr_kspi1_kpi2 = barrier(1,sD,sKsPip1,sKmPip2,rD2);
515 calt1(pKsPip1,pKmPip2,t1D);
516 for(int i=0; i<4; i++)
517 {
518 for(int j=0; j<4; j++)
519 {
520 tt1=G[i][i]*G[j][j]*pD[i]*t1D[j];
521 for(int k=0; k<4; k++)
522 {
523 tt2=t1KsPip1[k]*G[k][k];
524 for(int l=0; l<4; l++)
525 {
526 temp_PDF += E[i][j][k][l]*t1KmPip2[l]*G[l][l]*tt1*tt2;
527 }
528 }
529 }
530 }
531 tmp1 = barr_kspi1*barr_kpi2*barr_kspi1_kpi2*temp_PDF;
532
533 temp_PDF=0;
534 if(ispro_kspi2==0){
535 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
536 ispro_kspi2=1;
537 }
538 if(ispro_kpi1==0){
539 propagatorRBW(mKst0,GKst0,sKmPip1,sKm,sPip1,rRes2,1,pro_kpi1);
540 ispro_kpi1=1;
541 }
542 Com_Multi(pro_kspi2,pro_kpi1,pro2b);
543
544 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
545 if(barr_kpi1<0.0) barr_kpi1 = barrier(1,sKmPip1,sKm,sPip1,rRes2);
546 if(barr_kspi2_kpi1<0.0)barr_kspi2_kpi1 = barrier(1,sD,sKsPip2,sKmPip1,rD2);
547 calt1(pKsPip2,pKmPip1,t1D);
548 for(int i=0; i<4; i++)
549 {
550 for(int j=0; j<4; j++)
551 {
552 tt1 = G[i][i]*G[j][j]*pD[i]*t1D[j];
553 for(int k=0; k<4; k++)
554 {
555 tt2 = t1KsPip2[k]*G[k][k];
556 for(int l=0; l<4; l++)
557 {
558 temp_PDF += E[i][j][k][l]*t1KmPip1[l]*G[l][l]*tt1*tt2;
559 }
560 }
561 }
562 }
563 tmp2 = barr_kspi2*barr_kpi1*barr_kspi2_kpi1*temp_PDF;
564 }
565
566 // kstp_kst0_d
567 if(modetype[i]==3)
568 {
569 if (ispro_kspi1==0){
570 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
571 ispro_kspi1=1;
572 }
573 if (ispro_kpi2==0) {
574 propagatorRBW(mKst0,GKst0,sKmPip2,sKm,sPip2,rRes2,1,pro_kpi2);
575 ispro_kpi2=1;
576 }
577 Com_Multi(pro_kspi1,pro_kpi2,pro1b);
578
579 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
580 if(barr_kpi2<0.0) barr_kpi2 = barrier(1,sKmPip2,sKm,sPip2,rRes2);
581 calt2(pKsPip1,pKmPip2,t2D);
582 for(int i=0; i<4; i++)
583 {
584 for(int j=0; j<4; j++)
585 {
586 temp_PDF += t2D[i][j]*t1KsPip1[i]*t1KmPip2[j]*G[i][i]*G[j][j];
587 }
588 }
589 B[2] = barrier(2,sD,sKsPip1,sKmPip2,rD2);
590 tmp1 = barr_kspi1*barr_kpi2*B[2]*temp_PDF;
591
592 temp_PDF=0;
593 if(ispro_kspi2==0){
594 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
595 ispro_kspi2=1;
596 }
597 if(ispro_kpi1==0){
598 propagatorRBW(mKst0,GKst0,sKmPip1,sKm,sPip1,rRes2,1,pro_kpi1);
599 ispro_kpi1=1;
600 }
601 Com_Multi(pro_kspi2,pro_kpi1,pro2b);
602
603 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
604 if(barr_kpi1<0.0) barr_kpi1 = barrier(1,sKmPip1,sKm,sPip1,rRes2);
605 calt2(pKsPip2,pKmPip1,t2D);
606 for(int i=0; i<4; i++)
607 {
608 for(int j=0; j<4; j++)
609 {
610 temp_PDF += t2D[i][j]*t1KsPip2[i]*t1KmPip1[j]*G[i][i]*G[j][j];
611 }
612 }
613 B[2] = barrier(2,sD,sKsPip2,sKmPip1,rD2);
614 tmp2 = barr_kspi2*barr_kpi1*B[2]*temp_PDF;
615 }
616
617 /// eta1475->a0 pip
618 if(modetype[i]==8)
619 {
620 propagatorRBW(meta1475,Geta1475,sKsKmPip1,sKsKm,sPip1,rRes2,0,pro1);
621 if(ispro_ksk==0){
622 propagatora0980(ma0_980,sKsKm,sKs2,sKm2,pro_ksk);
623 ispro_ksk=1;
624 }
625 //printf("%.15f, %.15f , %.15f, %.15f \n", pro1[0],pro_ksk[0],pro1[1],pro_ksk[1]);
626 Com_Multi(pro1,pro_ksk,pro1b);
627 tmp1=1.0;
628
629 temp_PDF=0;
630 propagatorRBW(meta1475,Geta1475,sKsKmPip2,sKsKm,sPip2,rRes2,0,pro1);
631 if(ispro_ksk==0){
632 propagatora0980(ma0_980,sKsKm,sKs2,sKm2,pro_ksk);
633 ispro_ksk=1;
634 }
635 //printf("%.15f, %.15f , %.15f, %.15f \n", pro1[0],pro_ksk[0],pro1[1],pro_ksk[1]);
636 Com_Multi(pro1,pro_ksk,pro2b);
637 tmp2=1.0;
638 }
639
640 //Ds -> eta(1475)pi+ eta(1475)->K*Ks,K*0>Km pi+
641 if(modetype[i]==4)
642 {
643 propagatorRBW(meta1475,Geta1475,sKsKmPip1,sKmPip1,sKs,rRes2,1,pro1);
644 if(ispro_kpi1==0){
645 propagatorRBW(mKst0,GKst0,sKmPip1,sKm,sPip1,rRes2,1,pro_kpi1);
646 ispro_kpi1=1;
647 }
648 Com_Multi(pro1,pro_kpi1, pro1b);
649
650 if(barr_kpi1<0.0) barr_kpi1 = barrier(1,sKmPip1,sKm,sPip1,rRes2);
651 if(barr_kskpi1_kpi<0.0) barr_kskpi1_kpi = barrier(1,sKsKmPip1,sKmPip1,sKs,rRes2);
652 for(int a=0; a<4; a++)
653 {
654 temp_PDF += Ks[a]*t1KmPip1[a]*G[a][a];
655 }
656 tmp1 = barr_kpi1*barr_kskpi1_kpi*temp_PDF;
657
658 temp_PDF=0;
659 propagatorRBW(meta1475,Geta1475,sKsKmPip2,sKmPip2,sKs,rRes2,1,pro1);
660 if(ispro_kpi2==0){
661 propagatorRBW(mKst0,GKst0,sKmPip2,sKm,sPip2,rRes2,1,pro_kpi2);
662 ispro_kpi2=1;
663 }
664 Com_Multi(pro1,pro_kpi2,pro2b);
665
666 if(barr_kpi2<0.0) barr_kpi2 = barrier(1,sKmPip2,sKm,sPip2,rRes2);
667 if(barr_kskpi2_kpi<0.0) barr_kskpi2_kpi = barrier(1,sKsKmPip2,sKmPip2,sKs,rRes2);
668 for(int a=0; a<4; a++)
669 {
670 temp_PDF += Ks[a]*t1KmPip2[a]*G[a][a];
671 }
672 tmp2 = barr_kpi2*barr_kskpi2_kpi*temp_PDF;
673 }
674
675 //Ds -> eta(1475)pi+ eta(1475)->K*Km,K*+>Ks pi
676 if(modetype[i]==14)
677 {
678 propagatorRBW(meta1475,Geta1475,sKsKmPip1,sKsPip1,sKm,rRes2,1,pro1);
679 if(ispro_kspi1==0){
680 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
681 ispro_kspi1=1;
682 }
683 Com_Multi(pro1,pro_kspi1,pro1b);
684
685 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
686 if(barr_kskpi1_kspi<0.0)barr_kskpi1_kspi= barrier(1,sKsKmPip1,sKsPip1,sKm,rRes2);
687 for(int a=0; a<4; a++)
688 {
689 temp_PDF += Km[a]*t1KsPip1[a]*G[a][a];
690 }
691 tmp1 = barr_kspi1*barr_kskpi1_kspi*temp_PDF;
692
693 temp_PDF=0;
694 propagatorRBW(meta1475,Geta1475,sKsKmPip2,sKsPip2,sKm,rRes2,1,pro1);
695 if(ispro_kspi2==0){
696 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
697 ispro_kspi2=1;
698 }
699 Com_Multi(pro1,pro_kspi2,pro2b);
700
701 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
702 if(barr_kskpi2_kspi<0.0)barr_kskpi2_kspi= barrier(1,sKsKmPip2,sKsPip2,sKm,rRes2);
703 for(int a=0; a<4; a++)
704 {
705 temp_PDF += Km[a]*t1KsPip2[a]*G[a][a];
706 }
707 tmp2 = barr_kspi2*barr_kskpi2_kspi*temp_PDF;
708 }
709
710 //f1285->a0 pip
711 if(modetype[i]==17)
712 {
713 propagatorRBW(mf1285,Gf1285,sKsKmPip1,sKsKm,sPip1,rRes2,1,pro1);
714 if(ispro_ksk==0){
715 propagatora0980(ma0_980,sKsKm,sKs2,sKm2,pro_ksk);
716 ispro_ksk=1;
717 }
718 Com_Multi(pro1,pro_ksk,pro1b);
719 calt1(pKsKmPip1,Pip2,t1D);
720
721 if(barr_kskpi1_ksk<0.0) barr_kskpi1_ksk = barrier(1,sKsKmPip1,sKsKm,sPip1,rRes2);
722 if(barr_ds_kskpi1<0.0) barr_ds_kskpi1 = barrier(1,sD,sKsKmPip1,sPip2,rD2);
723 for(int i=0; i<4; i++)
724 {
725 temp_PDF += G[i][i]*t1D[i]*t1KsKmPip1_ksk[i];
726 }
727 tmp1 = barr_kskpi1_ksk*barr_ds_kskpi1*temp_PDF;
728
729 temp_PDF=0;
730 propagatorRBW(mf1285,Gf1285,sKsKmPip2,sKsKm,sPip2,rRes2,1,pro1);
731 if(ispro_ksk==0){
732 propagatora0980(ma0_980,sKsKm,sKs2,sKm2,pro_ksk);
733 ispro_ksk=1;
734 }
735 Com_Multi(pro1,pro_ksk,pro2b);
736 calt1(pKsKmPip2,Pip,t1D);
737
738 if(barr_kskpi2_ksk<0.0)barr_kskpi2_ksk = barrier(1,sKsKmPip2,sKsKm,sPip2,rRes2);
739 if(barr_ds_kskpi2<0.0) barr_ds_kskpi2 = barrier(1,sD,sKsKmPip2,sPip1,rD2);
740 for(int i=0; i<4; i++)
741 {
742 temp_PDF += G[i][i]*t1D[i]*t1KsKmPip2_ksk[i];
743 }
744 tmp2 = barr_kskpi2_ksk*barr_ds_kskpi2*temp_PDF;
745 }
746
747 //Ds -> 1550 pi+ , 1550 -> K-K*+ ,K*+ ->Ks pi+
748 if(modetype[i]==5){
749 //propagatorRBW(m1510,G1510,sKsKmPip1,sKsPip1,sKm,rRes2,0,pro1);
750 pro1[0]=1;pro1[1]=0;
751 if(ispro_kspi1==0){
752 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
753 ispro_kspi1=1;
754 }
755 //printf("%.15f, %.15f , %.15f, %.15f \n", pro1[0],pro_kspi1[0],pro1[1],pro_kspi1[1]);
756 Com_Multi(pro1,pro_kspi1,pro1b);
757 calt1(pKsKmPip1,Pip2,t1D);
758
759 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
760 if(barr_ds_kskpi1<0.0) barr_ds_kskpi1 = barrier(1,sD,sKsKmPip1,sPip2,rD2);
761 barr_kskpi1_kspi= barrier(2,sKsKmPip1,sKsPip1,sKm,rRes2);
762 for(int i=0; i<4; i++)
763 {
764 for(int j=0; j<4; j++)
765 {
766 temp_PDF += t1D[i]*(pKsKmPip1[i]*pKsKmPip1[j]/sKsKmPip1-G[i][j])*t1KsPip1[j]*G[i][i]*G[j][j];
767 }
768 }
769 tmp1 = barr_kspi1*barr_ds_kskpi1*temp_PDF;
770
771
772 temp_PDF=0;
773 //propagatorRBW(m1510,G1510,sKsKmPip2,sKsPip2,sKm,rRes2,0,pro1);
774 pro1[0]=1;pro1[1]=0;
775 if(ispro_kspi2==0){
776 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
777 ispro_kspi2=1;
778 }
779 Com_Multi(pro1,pro_kspi2,pro2b);
780 //printf("%.15f, %.15f , %.15f, %.15f \n", pro1[0],pro_kspi2[0],pro1[1],pro_kspi2[1]);
781
782 calt1(pKsKmPip2,Pip,t1D);
783 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
784 if(barr_ds_kskpi2<0.0) barr_ds_kskpi2 = barrier(1,sD,sKsKmPip2,sPip1,rD2);
785 barr_kskpi2_kspi = barrier(2,sKsKmPip2,sKsPip2,sKm,rRes2);
786 for(int i=0; i<4; i++)
787 {
788 for(int j=0; j<4; j++)
789 {
790 temp_PDF += t1D[i]*(pKsKmPip2[i]*pKsKmPip2[j]/sKsKmPip2-G[i][j])*t1KsPip2[j]*G[i][i]*G[j][j];
791 }
792 }
793 tmp2 = barr_kspi2*barr_ds_kskpi2*temp_PDF;
794 }
795
796 ///K*+ (Kpi)s_wave
797 if(modetype[i]==21)
798 {
799 if(ispro_kspi1==0){
800 propagatorRBW(mKstp,GKstp,sKsPip1,sKs,sPip1,rRes2,1,pro_kspi1);
801 ispro_kspi1=1;
802 }
803 KPiSLASS(sKmPip2,sKm,sPip2,pro_kpi2sw);
804 Com_Multi(pro_kspi1,pro_kpi2sw,pro1b);
805 calt1(pKsPip1,pKmPip2,t1D);
806
807 if(barr_kspi1<0.0) barr_kspi1 = barrier(1,sKsPip1,sKs,sPip1,rRes2);
808 if(barr_kspi1_kpi2<0.0) barr_kspi1_kpi2 = barrier(1,sD,sKsPip1,sKmPip2,rD2);
809 for(int a=0; a<4; a++)
810 {
811 temp_PDF += t1D[a]*t1KsPip1[a]*G[a][a];
812 }
813 tmp1 = barr_kspi1*barr_kspi1_kpi2*temp_PDF;
814
815 temp_PDF=0;
816 if(ispro_kspi2==0) {
817 propagatorRBW(mKstp,GKstp,sKsPip2,sKs,sPip2,rRes2,1,pro_kspi2);
818 ispro_kspi2=1;
819 }
820 KPiSLASS(sKmPip1,sKm,sPip1,pro_kpi1sw);
821 Com_Multi(pro_kspi2,pro_kpi1sw,pro2b);
822
823 calt1(pKsPip2,pKmPip1,t1D);
824 if(barr_kspi2<0.0) barr_kspi2 = barrier(1,sKsPip2,sKs,sPip2,rRes2);
825 if(barr_kspi2_kpi1<0.0) barr_kspi2_kpi1 = barrier(1,sD,sKsPip2,sKmPip1,rD2);
826 for(int a=0; a<4; a++)
827 {
828 temp_PDF += t1D[a]*t1KsPip2[a]*G[a][a];
829 }
830 tmp2 = barr_kspi2*barr_kspi2_kpi1*temp_PDF;
831 }
832
833 ///Kst0 Kspi_S_wave
834 if(modetype[i]==22)
835 {
836 if(ispro_kpi2==0){
837 propagatorRBW(mKst0,GKst0,sKmPip2,sKm,sPip2,rRes2,1,pro_kpi2);
838 ispro_kpi2=1;
839 }
840 KPiSLASS(sKsPip1,sKs,sPip1,pro_kspi1sw);
841 Com_Multi(pro_kpi2,pro_kspi1sw,pro1b);
842
843 calt1(pKsPip1,pKmPip2,t1D);
844 if(barr_kpi2<0.0) barr_kpi2 = barrier(1,sKmPip2,sKm,sPip2,rRes2);
845 if(barr_kspi1_kpi2<0.0) barr_kspi1_kpi2 = barrier(1,sD,sKsPip1,sKmPip2,rD2);
846 for(int a=0; a<4; a++)
847 {
848 temp_PDF += t1D[a]*t1KmPip2[a]*G[a][a];
849 }
850 tmp1 = barr_kpi2*barr_kspi1_kpi2*temp_PDF;
851
852 temp_PDF=0;
853 if(ispro_kpi1==0){
854 propagatorRBW(mKst0,GKst0,sKmPip1,sKm,sPip1,rRes2,1,pro_kpi1);
855 ispro_kpi1=1;
856 }
857 KPiSLASS(sKsPip2,sKs,sPip2,pro_kspi2sw);
858 Com_Multi(pro_kpi1,pro_kspi2sw,pro2b);
859
860 calt1(pKsPip2,pKmPip1,t1D);
861 if(barr_kpi1<0.0) barr_kpi1 = barrier(1,sKmPip1,sKm,sPip1,rRes2);
862 if(barr_kspi2_kpi1<0.0) barr_kspi2_kpi1 = barrier(1,sD,sKsPip2,sKmPip1,rD2);
863 for(int a=0; a<4; a++)
864 {
865 temp_PDF += t1D[a]*t1KmPip1[a]*G[a][a];
866 }
867 tmp2 = barr_kspi2_kpi1*barr_kpi1*temp_PDF;
868 }
869
870 if(modetype[i]==20)
871 {
872 propagatorRBW(meta1475,Geta1475,sKsKmPip1,sKsPip1,sKm,rRes2,0,pro1);
873 KPiSLASS(sKsPip1,sKs,sPip1,pro_kspi1sw);
874 Com_Multi(pro1,pro_kspi1sw,pro1b);
875 tmp1= 1.0;
876
877 propagatorRBW(meta1475,Geta1475,sKsKmPip2,sKsPip2,sKm,rRes2,0,pro1);
878 KPiSLASS(sKsPip2,sKs,sPip2,pro_kspi2sw);
879 Com_Multi(pro1,pro_kspi2sw,pro2b);
880 tmp2= 1.0;
881 }
882
883
884
885
886 if(modetype[i]==1||modetype[i]==2||modetype[i]==3||modetype[i]==8||modetype[i]==17||modetype[i]==5||modetype[i]==21||modetype[i]==22||modetype[i]==4||modetype[i]==14||modetype[i]==20){
887 amp_tmp[0] = tmp1*pro1b[0] + tmp2*pro2b[0];
888 amp_tmp[1] = tmp1*pro1b[1] + tmp2*pro2b[1];
889 }
890
891 /* if(modetype[i]==20){
892 amp_tmp[0] = tmp1*pro1b[0] + tmp2*pro2b[0];
893 amp_tmp[1] = tmp1*pro1b[1] + tmp2*pro2b[1];
894 } */
895 //printf("pro1b[0],pro2b[0]:: %.15f, %.15f\n", pro1b[0],pro2b[0]);
896 //printf("tmp1,tmp2:: %.15f, %.15f\n", tmp1,tmp2);
897 //printf("amp_tmp[0],amp_tmp[1]:: %.15f, %.15f\n", amp_tmp[0],amp_tmp[1]);
898
899 Com_Multi(amp_tmp,cof,amp_PDF);
900 PDF[0] += amp_PDF[0];
901 PDF[1] += amp_PDF[1];
902
903 //printf("cof= %.15f,%.15f \n",cof[0],cof[1]);
904 //printf("PDF= %.15f,%.15f \n",PDF[0],PDF[1]);
905 //printf("---------------------------------------------------\n");
906}
907 double value = PDF[0]*PDF[0] + PDF[1]*PDF[1];
908 if(value <=0) {value = 1e-20;}
909 Result = value;
910 //printf("value = %.15f\n",value);
911 //printf("---------------------------------------------------\n");
912}
double sin(const BesAngle a)
Definition BesAngle.h:210
double cos(const BesAngle a)
Definition BesAngle.h:213
double mass
*******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 !Latex Output unit
Definition FoamA.h:90
****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
TCrossPart * CS
Definition Mcgpj.cxx:51
TTree * t
Definition binning.cxx:23
void checkSpinDaughter(int d1, EvtSpinType::spintype sp)
void checkSpinParent(EvtSpinType::spintype sp)
void setProbMax(double prbmx)
void checkNDaug(int d1, int d2=-1)
EvtId * getDaugs()
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
void setProb(double prob)
void getName(std::string &name)
void decay(EvtParticle *p)
EvtDecayBase * clone()
const EvtVector4R & getP4() const
EvtParticle * getDaug(int i)
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
double get(int i) const
const double b
Definition slope.cxx:9