15 #include "Pythia8/Pythia.h"
17 using namespace Pythia8;
33 double normTheta2qqbar, normTheta2llbar, normTheta2ggg;
36 virtual void initConstants();
43 virtual void calcWidth(
bool =
false);
51 void ResonanceTheta::initConstants() {
54 normTheta2qqbar = 0.0001;
55 normTheta2llbar = 0.0001;
56 normTheta2ggg = 0.001;
63 void ResonanceTheta::calcWidth(
bool) {
66 if (id1Abs < 6) widNow = 3. * normTheta2qqbar * mHat;
69 else if (id1Abs == 11 || id1Abs == 13 || id1Abs == 15)
70 widNow = normTheta2llbar * mHat;
73 else if (id1Abs == 21) widNow = 8. * normTheta2ggg * mHat;
89 virtual void initProc();
92 virtual void sigmaKin();
95 virtual double sigmaHat() {
return sigma;}
98 virtual void setIdColAcol();
101 virtual double weightDecay(
Event& process,
int iResBeg,
int iResEnd);
104 virtual string name()
const {
return "q qbar -> Theta";}
105 virtual int code()
const {
return 621;}
106 virtual string inFlux()
const {
return "qqbarSame";}
107 virtual int resonanceA()
const {
return 663;}
113 double mRes, GammaRes, m2Res, GamMRat, normTheta2qqbar, sigma;
124 void Sigma1qqbar2Theta::initProc() {
128 mRes = particleDataPtr->m0(idTheta);
129 GammaRes = particleDataPtr->mWidth(idTheta);
131 GamMRat = GammaRes / mRes;
134 normTheta2qqbar = 0.0001;
137 particlePtr = particleDataPtr->particleDataEntryPtr(idTheta);
145 void Sigma1qqbar2Theta::sigmaKin() {
148 double widthIn = normTheta2qqbar * mH / 3.;
151 double sigBW = 12. * M_PI / ( pow2(sH - m2Res) + pow2(sH * GamMRat) );
154 double widthOut = particlePtr->resWidthOpen(663, mH);
157 sigma = widthIn * sigBW * widthOut;
165 void Sigma1qqbar2Theta::setIdColAcol() {
168 setId( id1, id2, idTheta);
171 setColAcol( 1, 0, 0, 1, 0, 0);
172 if (id1 < 0) swapColAcol();
180 double Sigma1qqbar2Theta::weightDecay(
Event& process,
int iResBeg,
184 if (iResEnd != iResBeg || process[iResBeg].idAbs() != idTheta)
188 int i1 = process[iResBeg].daughter1();
191 if (i3 != process[iResBeg].daughter2() || process[i1].
id() != 21)
195 double x1 = 2. * process[i1].p() * process[iResBeg].p()
196 / process[iResBeg].m2();
197 double x2 = 2. * process[i2].p() * process[iResBeg].p()
198 / process[iResBeg].m2();
199 double x3 = 2. * process[i3].p() * process[iResBeg].p()
200 / process[iResBeg].m2();
203 double wtME = pow2( (1. - x1) / (x2 * x3) )
204 + pow2( (1. - x2) / (x1 * x3) ) + pow2( (1. - x3) / (x1 * x2) );
206 return wtME / wtMEmax;
228 pythia.readString(
"663:new = Theta void 3 0 0 342.0 0.2 300. 400. 0.");
229 pythia.readString(
"663:addChannel = 1 0. 0 1 -1");
230 pythia.readString(
"663:addChannel = 1 0. 0 2 -2");
231 pythia.readString(
"663:addChannel = 1 0. 0 3 -3");
232 pythia.readString(
"663:addChannel = 1 0. 0 4 -4");
233 pythia.readString(
"663:addChannel = 1 0. 0 5 -5");
234 pythia.readString(
"663:addChannel = 1 0. 0 11 -11");
235 pythia.readString(
"663:addChannel = 1 0. 0 13 -13");
236 pythia.readString(
"663:addChannel = 1 0. 0 15 -15");
237 pythia.readString(
"663:addChannel = 1 0. 0 21 21 21");
244 pythia.setResonancePtr(resonanceTheta);
249 pythia.setSigmaPtr(sigma1Theta);
253 pythia.readString(
"Check:nErrList = 2");
259 Hist mTheta(
"Theta mass", 100, 300., 400.);
263 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
266 if (!pythia.next()) {
267 if (++iAbort < nAbort)
continue;
268 cout <<
" Event generation aborted prematurely, owing to error!\n";
273 mTheta.fill( pythia.process[5].m() );