StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HiddenValleyFragmentation.cc
1 // HiddenValleyFragmentation.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the
7 // HiddenValleyFragmentation class and its helper classes.
8 
9 #include "Pythia8/HiddenValleyFragmentation.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The HVStringFlav class is used to select HV-quark and HV-hadron flavours.
16 
17 //--------------------------------------------------------------------------
18 
19 // Initialize data members of the flavour generation.
20 
21 void HVStringFlav::init(Settings& settings, ParticleData* particleDataPtrIn,
22  Rndm* rndmPtrIn, Info* infoPtrIn) {
23 
24  // Save pointers.
25  particleDataPtr = particleDataPtrIn;
26  rndmPtr = rndmPtrIn;
27  infoPtr = infoPtrIn;
28 
29  // Read in data from Settings.
30  nFlav = settings.mode("HiddenValley:nFlav");
31  probVector = settings.parm("HiddenValley:probVector");
32  thermalModel = false;
33  useWidthPre = false;
34  closePacking = false;
35  mT2suppression = false;
36 
37 }
38 
39 //--------------------------------------------------------------------------
40 
41 // Pick a new HV-flavour given an incoming one.
42 
43 FlavContainer HVStringFlav::pick(FlavContainer& flavOld, double, double) {
44 
45  // Initial values for new flavour.
46  FlavContainer flavNew;
47  flavNew.rank = flavOld.rank + 1;
48 
49  // Pick new HV-flavour at random; keep track of sign.
50  flavNew.id = 4900100 + min( 1 + int(nFlav * rndmPtr->flat()), nFlav);
51  if (flavOld.id > 0) flavNew.id = -flavNew.id;
52 
53  // Done.
54  return flavNew;
55 
56 }
57 
58 //--------------------------------------------------------------------------
59 
60 // Combine two HV-flavours to produce an HV-hadron.
61 // This is simplified procedure, assuming only two HV mesons defined.
62 
63 int HVStringFlav::combine(FlavContainer& flav1, FlavContainer& flav2) {
64 
65  // Positive and negative flavour. Note that with kinetic mixing
66  // the Fv are really intended to represent qv, so remap.
67  int idMeson = 0;
68  int idPos = max( flav1.id, flav2.id) - 4900000;
69  int idNeg = -min( flav1.id, flav2.id) - 4900000;
70  if (idPos < 20) idPos = 101;
71  if (idNeg < 20) idNeg = 101;
72 
73  // Pick HV-meson code, spin either 0 or 1.
74  if (idNeg == idPos) idMeson = 4900111;
75  else if (idPos > idNeg) idMeson = 4900211;
76  else idMeson = -4900211;
77  if (rndmPtr->flat() < probVector) idMeson += ((idMeson > 0) ? 2 : -2);
78 
79  // Done.
80  return idMeson;
81 
82 }
83 
84 //==========================================================================
85 
86 // The HVStringPT class is used to select pT in HV fragmentation.
87 
88 //--------------------------------------------------------------------------
89 
90 // Initialize data members of the string pT selection.
91 
92 void HVStringPT::init(Settings& settings, ParticleData* particleDataPtrIn,
93  Rndm* rndmPtrIn, Info* infoPtrIn) {
94 
95  // Save pointer.
96  particleDataPtr = particleDataPtrIn;
97  rndmPtr = rndmPtrIn;
98  infoPtr = infoPtrIn;
99 
100  // Parameter of the pT width. No enhancement, since this is finetuning.
101  double sigmamqv = settings.parm("HiddenValley:sigmamqv");
102  double sigma = sigmamqv * particleDataPtrIn->m0( 4900101);
103  sigmaQ = sigma / sqrt(2.);
104  enhancedFraction = 0.;
105  enhancedWidth = 0.;
106 
107  // Parameter for pT suppression in MiniStringFragmentation.
108  sigma2Had = 2. * pow2( max( SIGMAMIN, sigma) );
109  thermalModel = false;
110  useWidthPre = false;
111  closePacking = false;
112 
113 }
114 
115 //==========================================================================
116 
117 // The HVStringZ class is used to select z in HV fragmentation.
118 
119 //--------------------------------------------------------------------------
120 
121 // Initialize data members of the string z selection.
122 
123 void HVStringZ::init(Settings& settings, ParticleData& particleData,
124  Rndm* rndmPtrIn, Info* infoPtrIn) {
125 
126  // Save pointer.
127  rndmPtr = rndmPtrIn;
128  infoPtr = infoPtrIn;
129 
130  // Paramaters of Lund/Bowler symmetric fragmentation function.
131  aLund = settings.parm("HiddenValley:aLund");
132  bmqv2 = settings.parm("HiddenValley:bmqv2");
133  rFactqv = settings.parm("HiddenValley:rFactqv");
134 
135  // Use qv mass to set scale of bEff = b * m^2;
136  mqv2 = pow2( particleData.m0( 4900101) );
137  bLund = bmqv2 / mqv2;
138 
139  // Mass of qv meson used to set stop scale for fragmentation iteration.
140  mhvMeson = particleData.m0( 4900111);
141 
142 }
143 
144 //--------------------------------------------------------------------------
145 
146 // Generate the fraction z that the next hadron will take using Lund/Bowler.
147 
148 double HVStringZ::zFrag( int , int , double mT2) {
149 
150  // Shape parameters of Lund symmetric fragmentation function.
151  double bShape = bLund * mT2;
152  double cShape = 1. + rFactqv * bmqv2;
153  return zLund( aLund, bShape, cShape);
154 
155 }
156 
157 //==========================================================================
158 
159 // The HiddenValleyFragmentation class.
160 
161 //--------------------------------------------------------------------------
162 
163 // Initialize and save pointers.
164 
165 bool HiddenValleyFragmentation::init(Info* infoPtrIn, Settings& settings,
166  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn) {
167 
168  // Save pointers.
169  infoPtr = infoPtrIn;
170  particleDataPtr = particleDataPtrIn;
171  rndmPtr = rndmPtrIn;
172 
173  // Check whether Hidden Valley fragmentation switched on, and SU(N).
174  doHVfrag = settings.flag("HiddenValley:fragment");
175  if (settings.mode("HiddenValley:Ngauge") < 2) doHVfrag = false;
176  if (!doHVfrag) return false;
177 
178  // Several copies of qv may be needed. Taken to have same mass.
179  nFlav = settings.mode("HiddenValley:nFlav");
180  if (nFlav > 1) {
181  int spinType = particleDataPtr->spinType(4900101);
182  double m0 = particleDataPtr->m0(4900101);
183  for (int iFlav = 2; iFlav <= nFlav; ++iFlav)
184  particleDataPtr->addParticle( 4900100 + iFlav, "qv", "qvbar",
185  spinType, 0, 0, m0);
186  }
187 
188  // Hidden Valley meson mass used to choose hadronization mode.
189  mhvMeson = particleDataPtr->m0(4900111);
190 
191  // Initialize the hvEvent instance of an event record.
192  hvEvent.init( "(Hidden Valley fragmentation)", particleDataPtr);
193 
194  // Create HVStringFlav instance for HV-flavour selection.
195  hvFlavSelPtr = new HVStringFlav();
196  hvFlavSelPtr->init( settings, particleDataPtr, rndmPtr, infoPtr);
197 
198  // Create HVStringPT instance for pT selection in HV fragmentation.
199  hvPTSelPtr = new HVStringPT();
200  hvPTSelPtr->init( settings, particleDataPtr, rndmPtr, infoPtr);
201 
202  // Create HVStringZ instance for z selection in HV fragmentation.
203  hvZSelPtr = new HVStringZ();
204  hvZSelPtr->init( settings, *particleDataPtr, rndmPtr, infoPtr);
205 
206  // Initialize auxiliary administrative class.
207  hvColConfig.init(infoPtr, settings, hvFlavSelPtr);
208 
209  // Initialize HV-string and HV-ministring fragmentation.
210  hvStringFrag.init(infoPtr, settings, particleDataPtr, rndmPtr,
211  hvFlavSelPtr, hvPTSelPtr, hvZSelPtr);
212  hvMinistringFrag.init(infoPtr, settings, particleDataPtr, rndmPtr,
213  hvFlavSelPtr, hvPTSelPtr, hvZSelPtr);
214 
215  // Done.
216  return true;
217 
218 }
219 
220 //--------------------------------------------------------------------------
221 
222 // Perform the fragmentation.
223 
224 bool HiddenValleyFragmentation::fragment(Event& event) {
225 
226  // Reset containers for next event.
227  hvEvent.reset();
228  hvColConfig.clear();
229  ihvParton.resize(0);
230 
231  // Extract HV-particles from event to hvEvent. Assign HV-colours.
232  // Done if no HV-particles found.
233  if (!extractHVevent(event)) return true;
234 
235  // Store found string system. Analyze its properties.
236  if (!hvColConfig.insert(ihvParton, hvEvent)) return false;
237 
238  // Collect sequentially all partons in the HV subsystem.
239  // Copy also if already in order, or else history tracing may fail.
240  hvColConfig.collect(0, hvEvent, false);
241 
242  // Mass used to decide how to fragment system.
243  mSys = hvColConfig[0].mass;
244 
245  // HV-string fragmentation when enough mass to produce >= 3 HV-mesons.
246  if (mSys > 3.5 * mhvMeson) {
247  if (!hvStringFrag.fragment( 0, hvColConfig, hvEvent)) return false;
248 
249  // HV-ministring fragmentation when enough mass to produce 2 HV-mesons.
250  } else if (mSys > 2.1 * mhvMeson) {
251  if (!hvMinistringFrag.fragment( 0, hvColConfig, hvEvent, true))
252  return false;
253 
254  // If only enough mass for one HV-meson assume HV-glueballs emitted.
255  } else if (!collapseToMeson()) return false;
256 
257  // Insert HV particles from hvEvent to event.
258  insertHVevent(event);
259 
260  // Done.
261  return true;
262 
263 }
264 
265 //--------------------------------------------------------------------------
266 
267 // Extract HV-particles from event to hvEvent. Assign HV-colours.
268 
269 bool HiddenValleyFragmentation::extractHVevent(Event& event) {
270 
271  // Copy Hidden-Valley particles to special event record.
272  for (int i = 0; i < event.size(); ++i) {
273  int idAbs = event[i].idAbs();
274  bool isHV = (idAbs > 4900000 && idAbs < 4900007)
275  || (idAbs > 4900010 && idAbs < 4900017)
276  || idAbs == 4900021
277  || (idAbs > 4900100 && idAbs < 4900109);
278  if (isHV) {
279  int iHV = hvEvent.append( event[i]);
280  // Convert HV-gluons into normal ones so as to use normal machinery.
281  if (event[i].id() == 4900021) hvEvent[iHV].id(21);
282  // Second mother points back to position in complete event;
283  // otherwise construct the HV history inside hvEvent.
284  hvEvent[iHV].mothers( 0, i);
285  hvEvent[iHV].daughters( 0, 0);
286  int iMother = event[i].mother1();
287  for (int iHVM = 1; iHVM < hvEvent.size(); ++iHVM)
288  if (hvEvent[iHVM].mother2() == iMother) {
289  hvEvent[iHV].mother1( iHVM);
290  if (hvEvent[iHVM].daughter1() == 0) hvEvent[iHVM].daughter1(iHV);
291  else hvEvent[iHVM].daughter2(iHV);
292  }
293  }
294  }
295 
296  // Done if no HV particles found.
297  hvOldSize = hvEvent.size();
298  if (hvOldSize == 1) return false;
299 
300  // Initial colour - anticolour parton pair.
301  int colBeg = hvEvent.nextColTag();
302  for (int iHV = 1; iHV < hvOldSize; ++iHV)
303  if (hvEvent[iHV].mother1() == 0) {
304  if (hvEvent[iHV].id() > 0) hvEvent[iHV].col( colBeg);
305  else hvEvent[iHV].acol( colBeg);
306  }
307 
308  // Then trace colours down to daughters; new colour if two daughters.
309  for (int iHV = 1; iHV < hvOldSize; ++iHV) {
310  int dau1 = hvEvent[iHV].daughter1();
311  int dau2 = hvEvent[iHV].daughter2();
312  if (dau1 > 0 && dau2 == 0)
313  hvEvent[dau1].cols( hvEvent[iHV].col(), hvEvent[iHV].acol());
314  else if (dau2 > 0) {
315  int colHV = hvEvent[iHV].col();
316  int acolHV = hvEvent[iHV].acol();
317  int colNew = hvEvent.nextColTag();
318  if (acolHV == 0) {
319  hvEvent[dau1].cols( colNew, 0);
320  hvEvent[dau2].cols( colHV, colNew);
321  } else if (colHV == 0) {
322  hvEvent[dau1].cols( 0, colNew);
323  hvEvent[dau2].cols( colNew, acolHV);
324  // Temporary: should seek recoiling dipole end!??
325  } else if (rndmPtr->flat() > 0.5) {
326  hvEvent[dau1].cols( colHV, colNew);
327  hvEvent[dau2].cols( colNew, acolHV);
328  } else {
329  hvEvent[dau1].cols( colNew, acolHV);
330  hvEvent[dau2].cols( colHV, colNew);
331  }
332  }
333  }
334 
335  // Pick up the colour end.
336  int colNow = 0;
337  for (int iHV = 1; iHV < hvOldSize; ++iHV)
338  if (hvEvent[iHV].isFinal() && hvEvent[iHV].acol() == 0) {
339  ihvParton.push_back( iHV);
340  colNow = hvEvent[iHV].col();
341  }
342 
343  // Trace colour by colour until reached anticolour end.
344  while (colNow > 0) {
345  for (int iHV = 1; iHV < hvOldSize; ++iHV)
346  if (hvEvent[iHV].isFinal() && hvEvent[iHV].acol() == colNow) {
347  ihvParton.push_back( iHV);
348  colNow = hvEvent[iHV].col();
349  break;
350  }
351  }
352 
353  // Done.
354  return true;
355 
356 }
357 
358 //--------------------------------------------------------------------------
359 
360 // Collapse of light system to one HV-meson, by the emission of HV-glueballs.
361 
362 bool HiddenValleyFragmentation::collapseToMeson() {
363 
364  // If too low mass then cannot do anything. Should not happen.
365  if (mSys < 1.001 * mhvMeson) {
366  infoPtr->errorMsg("Error in HiddenValleyFragmentation::collapseToMeson:"
367  " too low mass to do anything");
368  return false;
369  }
370 
371  // Choose mass of collective HV-glueball states flat between limits.
372  double mhvGlue = (0.001 + 0.998 * rndmPtr->flat()) * (mSys - mhvMeson);
373 
374  // Find momentum in rest frame, with isotropic "decay" angles.
375  double pAbs = 0.5 * sqrtpos( pow2(mSys*mSys - mhvMeson*mhvMeson
376  - mhvGlue*mhvGlue) - pow2(2. * mhvMeson * mhvGlue) ) / mSys;
377  double pz = (2 * rndmPtr->flat() - 1.) * pAbs;
378  double pT = sqrtpos( pAbs*pAbs - pz*pz);
379  double phi = 2. * M_PI * rndmPtr->flat();
380  double px = pT * cos(phi);
381  double py = pT * sin(phi);
382 
383  // Construct four-vectors and boost them to event frame.
384  Vec4 phvMeson( px, py, pz, sqrt(mhvMeson*mhvMeson + pAbs*pAbs) );
385  Vec4 phvGlue( -px, -py, -pz, sqrt(mhvGlue*mhvGlue + pAbs*pAbs) );
386  phvMeson.bst( hvColConfig[0].pSum );
387  phvGlue.bst( hvColConfig[0].pSum );
388 
389  // Add produced particles to the event record.
390  vector<int> iParton = hvColConfig[0].iParton;
391  int iFirst = hvEvent.append( 4900111, 82, iParton.front(),
392  iParton.back(), 0, 0, 0, 0, phvMeson, mhvMeson);
393  int iLast = hvEvent.append( 4900991, 82, iParton.front(),
394  iParton.back(), 0, 0, 0, 0, phvGlue, mhvGlue);
395 
396  // Mark original partons as hadronized and set their daughter range.
397  for (int i = 0; i < int(iParton.size()); ++i) {
398  hvEvent[ iParton[i] ].statusNeg();
399  hvEvent[ iParton[i] ].daughters(iFirst, iLast);
400  }
401 
402  // Done.
403  return true;
404 
405 }
406 
407 //--------------------------------------------------------------------------
408 
409 // Insert HV-particles from hvEvent to event.
410 
411 bool HiddenValleyFragmentation::insertHVevent(Event& event) {
412 
413  // Offset for mother/daughter indices.
414  hvNewSize = hvEvent.size();
415  int nOffset = event.size() - hvOldSize;
416 
417  // Copy back HV-particles.
418  int iNew, iMot1, iMot2, iDau1, iDau2;
419  for (int iHV = hvOldSize; iHV < hvNewSize; ++iHV) {
420  iNew = event.append( hvEvent[iHV]);
421 
422  // Restore HV-gluon codes. Do not keep HV-colours, to avoid confusion.
423  if (hvEvent[iHV].id() == 21) event[iNew].id(4900021);
424  event[iNew].cols( 0, 0);
425 
426  // Begin history construction.
427  iMot1 = hvEvent[iHV].mother1();
428  iMot2 = hvEvent[iHV].mother2();
429  iDau1 = hvEvent[iHV].daughter1();
430  iDau2 = hvEvent[iHV].daughter2();
431  // Special mother for partons copied from event, else simple offset.
432  // Also set daughters of mothers in original record.
433  if (iMot1 > 0 && iMot1 < hvOldSize) {
434  iMot1 = hvEvent[iMot1].mother2();
435  event[iMot1].statusNeg();
436  event[iMot1].daughter1(iNew);
437  } else if (iMot1 > 0) iMot1 += nOffset;
438  if (iMot2 > 0 && iMot2 < hvOldSize) {
439  iMot2 = hvEvent[iMot2].mother2();
440  event[iMot2].statusNeg();
441  if (event[iMot2].daughter1() == 0) event[iMot2].daughter1(iNew);
442  else event[iMot2].daughter2(iNew);
443  } else if (iMot2 > 0) iMot2 += nOffset;
444  if (iDau1 > 0) iDau1 += nOffset;
445  if (iDau2 > 0) iDau2 += nOffset;
446  event[iNew].mothers( iMot1, iMot2);
447  event[iNew].daughters( iDau1, iDau2);
448  }
449 
450  // Done.
451  return true;
452 
453 }
454 
455 //==========================================================================
456 
457 } // end namespace Pythia8
Definition: AgUStep.h:26