9 #include "Pythia8/HiddenValleyFragmentation.h"
21 void HVStringFlav::init() {
24 nFlav = mode(
"HiddenValley:nFlav");
25 probVector = parm(
"HiddenValley:probVector");
29 mT2suppression =
false;
37 FlavContainer HVStringFlav::pick(FlavContainer& flavOld,
double,
double,
41 FlavContainer flavNew;
42 flavNew.rank = flavOld.rank + 1;
45 flavNew.id = 4900100 + min( 1 +
int(nFlav * rndmPtr->flat()), nFlav);
46 if (flavOld.id > 0) flavNew.id = -flavNew.id;
58 int HVStringFlav::combine(FlavContainer& flav1, FlavContainer& flav2) {
63 int idPos = max( flav1.id, flav2.id) - 4900000;
64 int idNeg = -min( flav1.id, flav2.id) - 4900000;
65 if (idPos < 20) idPos = 101;
66 if (idNeg < 20) idNeg = 101;
69 if (idNeg == idPos) idMeson = 4900111;
70 else if (idPos > idNeg) idMeson = 4900211;
71 else idMeson = -4900211;
72 if (rndmPtr->flat() < probVector) idMeson += ((idMeson > 0) ? 2 : -2);
87 void HVStringPT::init() {
90 double sigmamqv = parm(
"HiddenValley:sigmamqv");
91 double sigma = sigmamqv * particleDataPtr->m0( 4900101);
92 sigmaQ = sigma / sqrt(2.);
93 enhancedFraction = 0.;
97 sigma2Had = 2. * pow2( max( SIGMAMIN, sigma) );
100 closePacking =
false;
112 void HVStringZ::init() {
115 aLund = parm(
"HiddenValley:aLund");
116 bmqv2 = parm(
"HiddenValley:bmqv2");
117 rFactqv = parm(
"HiddenValley:rFactqv");
120 mqv2 = pow2( particleDataPtr->m0( 4900101) );
121 bLund = bmqv2 / mqv2;
124 mhvMeson = particleDataPtr->m0( 4900111);
132 double HVStringZ::zFrag(
int ,
int ,
double mT2) {
135 double bShape = bLund * mT2;
136 double cShape = 1. + rFactqv * bmqv2;
137 return zLund( aLund, bShape, cShape);
149 bool HiddenValleyFragmentation::init() {
152 doHVfrag = flag(
"HiddenValley:fragment");
153 if (mode(
"HiddenValley:Ngauge") < 2) doHVfrag =
false;
154 if (!doHVfrag)
return false;
157 nFlav = mode(
"HiddenValley:nFlav");
159 int spinType = particleDataPtr->spinType(4900101);
160 double m0 = particleDataPtr->m0(4900101);
161 for (
int iFlav = 2; iFlav <= nFlav; ++iFlav)
162 particleDataPtr->addParticle( 4900100 + iFlav,
"qv",
"qvbar",
167 mhvMeson = particleDataPtr->m0(4900111);
170 hvEvent.init(
"(Hidden Valley fragmentation)", particleDataPtr);
182 hvColConfig.init(infoPtr, &hvFlavSel);
185 hvStringFrag.init(&hvFlavSel, &hvPTSel, &hvZSel);
186 hvMinistringFrag.init(&hvFlavSel, &hvPTSel, &hvZSel);
197 bool HiddenValleyFragmentation::fragment(
Event& event) {
206 if (!extractHVevent(event))
return true;
209 if (!hvColConfig.insert(ihvParton, hvEvent))
return false;
213 hvColConfig.collect(0, hvEvent,
false);
216 mSys = hvColConfig[0].mass;
219 if (mSys > 3.5 * mhvMeson) {
220 if (!hvStringFrag.fragment( 0, hvColConfig, hvEvent))
return false;
223 }
else if (mSys > 2.1 * mhvMeson) {
224 if (!hvMinistringFrag.fragment( 0, hvColConfig, hvEvent,
true))
228 }
else if (!collapseToMeson())
return false;
231 insertHVevent(event);
242 bool HiddenValleyFragmentation::extractHVevent(
Event& event) {
245 for (
int i = 0; i <
event.size(); ++i) {
246 int idAbs =
event[i].idAbs();
247 bool isHV = (idAbs > 4900000 && idAbs < 4900007)
248 || (idAbs > 4900010 && idAbs < 4900017)
250 || (idAbs > 4900100 && idAbs < 4900109);
252 int iHV = hvEvent.append( event[i]);
254 if (event[i].
id() == 4900021) hvEvent[iHV].id(21);
257 hvEvent[iHV].mothers( 0, i);
258 hvEvent[iHV].daughters( 0, 0);
259 int iMother =
event[i].mother1();
260 for (
int iHVM = 1; iHVM < hvEvent.size(); ++iHVM)
261 if (hvEvent[iHVM].mother2() == iMother) {
262 hvEvent[iHV].mother1( iHVM);
263 if (hvEvent[iHVM].daughter1() == 0) hvEvent[iHVM].daughter1(iHV);
264 else hvEvent[iHVM].daughter2(iHV);
270 hvOldSize = hvEvent.size();
271 if (hvOldSize == 1)
return false;
274 int colBeg = hvEvent.nextColTag();
275 for (
int iHV = 1; iHV < hvOldSize; ++iHV)
276 if (hvEvent[iHV].mother1() == 0) {
277 if (hvEvent[iHV].
id() > 0) hvEvent[iHV].col( colBeg);
278 else hvEvent[iHV].acol( colBeg);
282 for (
int iHV = 1; iHV < hvOldSize; ++iHV) {
283 int dau1 = hvEvent[iHV].daughter1();
284 int dau2 = hvEvent[iHV].daughter2();
285 if (dau1 > 0 && dau2 == 0)
286 hvEvent[dau1].cols( hvEvent[iHV].col(), hvEvent[iHV].acol());
288 int colHV = hvEvent[iHV].col();
289 int acolHV = hvEvent[iHV].acol();
290 int colNew = hvEvent.nextColTag();
292 hvEvent[dau1].cols( colNew, 0);
293 hvEvent[dau2].cols( colHV, colNew);
294 }
else if (colHV == 0) {
295 hvEvent[dau1].cols( 0, colNew);
296 hvEvent[dau2].cols( colNew, acolHV);
298 }
else if (rndmPtr->flat() > 0.5) {
299 hvEvent[dau1].cols( colHV, colNew);
300 hvEvent[dau2].cols( colNew, acolHV);
302 hvEvent[dau1].cols( colNew, acolHV);
303 hvEvent[dau2].cols( colHV, colNew);
310 for (
int iHV = 1; iHV < hvOldSize; ++iHV)
311 if (hvEvent[iHV].isFinal() && hvEvent[iHV].acol() == 0) {
312 ihvParton.push_back( iHV);
313 colNow = hvEvent[iHV].col();
318 for (
int iHV = 1; iHV < hvOldSize; ++iHV)
319 if (hvEvent[iHV].isFinal() && hvEvent[iHV].acol() == colNow) {
320 ihvParton.push_back( iHV);
321 colNow = hvEvent[iHV].col();
335 bool HiddenValleyFragmentation::collapseToMeson() {
338 if (mSys < 1.001 * mhvMeson) {
339 infoPtr->errorMsg(
"Error in HiddenValleyFragmentation::collapseToMeson:"
340 " too low mass to do anything");
345 double mhvGlue = (0.001 + 0.998 * rndmPtr->flat()) * (mSys - mhvMeson);
348 double pAbs = 0.5 * sqrtpos( pow2(mSys*mSys - mhvMeson*mhvMeson
349 - mhvGlue*mhvGlue) - pow2(2. * mhvMeson * mhvGlue) ) / mSys;
350 double pz = (2 * rndmPtr->flat() - 1.) * pAbs;
351 double pT = sqrtpos( pAbs*pAbs - pz*pz);
352 double phi = 2. * M_PI * rndmPtr->flat();
353 double px = pT * cos(phi);
354 double py = pT * sin(phi);
357 Vec4 phvMeson( px, py, pz, sqrt(mhvMeson*mhvMeson + pAbs*pAbs) );
358 Vec4 phvGlue( -px, -py, -pz, sqrt(mhvGlue*mhvGlue + pAbs*pAbs) );
359 phvMeson.bst( hvColConfig[0].pSum );
360 phvGlue.bst( hvColConfig[0].pSum );
363 vector<int> iParton = hvColConfig[0].iParton;
364 int iFirst = hvEvent.append( 4900111, 82, iParton.front(),
365 iParton.back(), 0, 0, 0, 0, phvMeson, mhvMeson);
366 int iLast = hvEvent.append( 4900991, 82, iParton.front(),
367 iParton.back(), 0, 0, 0, 0, phvGlue, mhvGlue);
370 for (
int i = 0; i < int(iParton.size()); ++i) {
371 hvEvent[ iParton[i] ].statusNeg();
372 hvEvent[ iParton[i] ].daughters(iFirst, iLast);
384 bool HiddenValleyFragmentation::insertHVevent(
Event& event) {
387 hvNewSize = hvEvent.size();
388 int nOffset =
event.size() - hvOldSize;
391 int iNew, iMot1, iMot2, iDau1, iDau2;
392 for (
int iHV = hvOldSize; iHV < hvNewSize; ++iHV) {
393 iNew =
event.append( hvEvent[iHV]);
396 if (hvEvent[iHV].
id() == 21)
event[iNew].id(4900021);
397 event[iNew].cols( 0, 0);
400 iMot1 = hvEvent[iHV].mother1();
401 iMot2 = hvEvent[iHV].mother2();
402 iDau1 = hvEvent[iHV].daughter1();
403 iDau2 = hvEvent[iHV].daughter2();
406 if (iMot1 > 0 && iMot1 < hvOldSize) {
407 iMot1 = hvEvent[iMot1].mother2();
408 event[iMot1].statusNeg();
409 event[iMot1].daughter1(iNew);
410 }
else if (iMot1 > 0) iMot1 += nOffset;
411 if (iMot2 > 0 && iMot2 < hvOldSize) {
412 iMot2 = hvEvent[iMot2].mother2();
413 event[iMot2].statusNeg();
414 if (event[iMot2].daughter1() == 0)
event[iMot2].daughter1(iNew);
415 else event[iMot2].daughter2(iNew);
416 }
else if (iMot2 > 0) iMot2 += nOffset;
417 if (iDau1 > 0) iDau1 += nOffset;
418 if (iDau2 > 0) iDau2 += nOffset;
419 event[iNew].mothers( iMot1, iMot2);
420 event[iNew].daughters( iDau1, iDau2);