StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MiniStringFragmentation.cc
1 // MiniStringFragmentation.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, 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 // MiniStringFragmentation class
8 
9 #include "Pythia8/MiniStringFragmentation.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The MiniStringFragmentation class.
16 
17 //--------------------------------------------------------------------------
18 
19 // Constants: could be changed here if desired, but normally should not.
20 // These are of technical nature, as described for each.
21 
22 // Since diffractive by definition is > 1 particle, try hard.
23 const int MiniStringFragmentation::NTRYDIFFRACTIVE = 200;
24 
25 // After one-body fragmentation failed, try two-body once more.
26 const int MiniStringFragmentation::NTRYLASTRESORT = 100;
27 
28 // Loop try to combine available endquarks to valid hadron.
29 const int MiniStringFragmentation::NTRYFLAV = 10;
30 
31 //--------------------------------------------------------------------------
32 
33 // Initialize and save pointers.
34 
35 void MiniStringFragmentation::init(Info* infoPtrIn, Settings& settings,
36  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
37  StringFlav* flavSelPtrIn, StringPT* pTSelPtrIn, StringZ* zSelPtrIn) {
38 
39  // Save pointers.
40  infoPtr = infoPtrIn;
41  particleDataPtr = particleDataPtrIn;
42  rndmPtr = rndmPtrIn;
43  flavSelPtr = flavSelPtrIn;
44  pTSelPtr = pTSelPtrIn;
45  zSelPtr = zSelPtrIn;
46 
47  // Initialize the MiniStringFragmentation class proper.
48  nTryMass = settings.mode("MiniStringFragmentation:nTry");
49 
50  // Initialize the b parameter of the z spectrum, used when joining jets.
51  bLund = zSelPtr->bAreaLund();
52 
53 }
54 
55 //--------------------------------------------------------------------------
56 
57 // Do the fragmentation: driver routine.
58 
59 bool MiniStringFragmentation::fragment(int iSub, ColConfig& colConfig,
60  Event& event, bool isDiff) {
61 
62  // Junction topologies not yet considered - is very rare.
63  iParton = colConfig[iSub].iParton;
64  if (iParton.front() < 0) {
65  infoPtr->errorMsg("Error in MiniStringFragmentation::fragment: "
66  "very low-mass junction topologies not yet handled");
67  return false;
68  }
69 
70  // Read in info on system to be treated.
71  flav1 = FlavContainer( event[ iParton.front() ].id() );
72  flav2 = FlavContainer( event[ iParton.back() ].id() );
73  pSum = colConfig[iSub].pSum;
74  mSum = colConfig[iSub].mass;
75  m2Sum = mSum*mSum;
76  isClosed = colConfig[iSub].isClosed;
77 
78  // Do not want diffractive systems to easily collapse to one particle.
79  int nTryFirst = (isDiff) ? NTRYDIFFRACTIVE : nTryMass;
80 
81  // First try to produce two particles from the system.
82  if (ministring2two( nTryFirst, event)) return true;
83 
84  // If this fails, then form one hadron and shuffle momentum.
85  if (ministring2one( iSub, colConfig, event)) return true;
86 
87  // If also this fails, then try harder to produce two particles.
88  if (ministring2two( NTRYLASTRESORT, event)) return true;
89 
90  // Else complete failure.
91  infoPtr->errorMsg("Error in MiniStringFragmentation::fragment: "
92  "no 1- or 2-body state found above mass threshold");
93  return false;
94 
95 }
96 
97 //--------------------------------------------------------------------------
98 
99  // Attempt to produce two particles from the ministring.
100 
101 bool MiniStringFragmentation::ministring2two( int nTry, Event& event) {
102 
103  // Properties of the produced hadrons.
104  int idHad1 = 0;
105  int idHad2 = 0;
106  double mHad1 = 0.;
107  double mHad2 = 0.;
108  double mHadSum = 0.;
109 
110  // Allow a few attempts to find a particle pair with low enough masses.
111  for (int iTry = 0; iTry < nTry; ++iTry) {
112 
113  // For closed gluon loop need to pick an initial flavour.
114  if (isClosed) do {
115  int idStart = flavSelPtr->pickLightQ();
116  FlavContainer flavStart(idStart, 1);
117  flavStart = flavSelPtr->pick( flavStart);
118  flav1 = flavSelPtr->pick( flavStart);
119  flav2.anti(flav1);
120  } while (flav1.id == 0 || flav1.nPop > 0);
121 
122  // Create a new q qbar flavour to form two hadrons.
123  // Start from a diquark, if any.
124  do {
125  FlavContainer flav3 =
126  (flav1.isDiquark() || (!flav2.isDiquark() && rndmPtr->flat() < 0.5) )
127  ? flavSelPtr->pick( flav1) : flavSelPtr->pick( flav2).anti();
128  idHad1 = flavSelPtr->combine( flav1, flav3);
129  idHad2 = flavSelPtr->combine( flav2, flav3.anti());
130  } while (idHad1 == 0 || idHad2 == 0);
131 
132  // Check whether the mass sum fits inside the available phase space.
133  mHad1 = particleDataPtr->mSel(idHad1);
134  mHad2 = particleDataPtr->mSel(idHad2);
135  mHadSum = mHad1 + mHad2;
136  if (mHadSum < mSum) break;
137  }
138  if (mHadSum >= mSum) return false;
139 
140  // Define an effective two-parton string, by splitting intermediate
141  // gluon momenta in proportion to their closeness to either endpoint.
142  Vec4 pSum1 = event[ iParton.front() ].p();
143  Vec4 pSum2 = event[ iParton.back() ].p();
144  if (iParton.size() > 2) {
145  Vec4 pEnd1 = pSum1;
146  Vec4 pEnd2 = pSum2;
147  Vec4 pEndSum = pEnd1 + pEnd2;
148  for (int i = 1; i < int(iParton.size()) - 1 ; ++i) {
149  Vec4 pNow = event[ iParton[i] ].p();
150  double ratio = (pEnd2 * pNow) / (pEndSum * pNow);
151  pSum1 += ratio * pNow;
152  pSum2 += (1. - ratio) * pNow;
153  }
154  }
155 
156  // Set up a string region based on the two effective endpoints.
157  StringRegion region;
158  region.setUp( pSum1, pSum2);
159 
160  // Generate an isotropic decay in the ministring rest frame,
161  // suppressed at large pT by a fragmentation pT Gaussian.
162  double pAbs2 = 0.25 * ( pow2(m2Sum - mHad1*mHad1 - mHad2*mHad2)
163  - pow2(2. * mHad1 * mHad2) ) / m2Sum;
164  double pT2 = 0.;
165  do {
166  double cosTheta = rndmPtr->flat();
167  pT2 = (1. - pow2(cosTheta)) * pAbs2;
168  } while (pTSelPtr->suppressPT2(pT2) < rndmPtr->flat() );
169 
170  // Construct the forward-backward asymmetry of the two particles.
171  double mT21 = mHad1*mHad1 + pT2;
172  double mT22 = mHad2*mHad2 + pT2;
173  double lambda = sqrtpos( pow2(m2Sum - mT21 - mT22) - 4. * mT21 * mT22 );
174  double probReverse = 1. / (1. + exp( min( 50., bLund * lambda) ) );
175 
176  // Construct kinematics, as viewed in the transverse rest frame.
177  double xpz1 = 0.5 * lambda/ m2Sum;
178  if (probReverse > rndmPtr->flat()) xpz1 = -xpz1;
179  double xmDiff = (mT21 - mT22) / m2Sum;
180  double xe1 = 0.5 * (1. + xmDiff);
181  double xe2 = 0.5 * (1. - xmDiff );
182 
183  // Distribute pT isotropically in angle.
184  double phi = 2. * M_PI * rndmPtr->flat();
185  double pT = sqrt(pT2);
186  double px = pT * cos(phi);
187  double py = pT * sin(phi);
188 
189  // Translate this into kinematics in the string frame.
190  Vec4 pHad1 = region.pHad( xe1 + xpz1, xe1 - xpz1, px, py);
191  Vec4 pHad2 = region.pHad( xe2 - xpz1, xe2 + xpz1, -px, -py);
192 
193  // Add produced particles to the event record.
194  int iFirst = event.append( idHad1, 82, iParton.front(), iParton.back(),
195  0, 0, 0, 0, pHad1, mHad1);
196  int iLast = event.append( idHad2, 82, iParton.front(), iParton.back(),
197  0, 0, 0, 0, pHad2, mHad2);
198 
199  // Set decay vertex when this is displaced.
200  if (event[iParton.front()].hasVertex()) {
201  Vec4 vDec = event[iParton.front()].vDec();
202  event[iFirst].vProd( vDec );
203  event[iLast].vProd( vDec );
204  }
205 
206  // Set lifetime of hadrons.
207  event[iFirst].tau( event[iFirst].tau0() * rndmPtr->exp() );
208  event[iLast].tau( event[iLast].tau0() * rndmPtr->exp() );
209 
210  // Mark original partons as hadronized and set their daughter range.
211  for (int i = 0; i < int(iParton.size()); ++i) {
212  event[ iParton[i] ].statusNeg();
213  event[ iParton[i] ].daughters(iFirst, iLast);
214  }
215 
216  // Successfully done.
217  return true;
218 
219 }
220 
221 //--------------------------------------------------------------------------
222 
223 // Attempt to produce one particle from a ministring.
224 // Current algorithm: find the system with largest invariant mass
225 // relative to the existing one, and boost that system appropriately.
226 // Try more sophisticated alternatives later?? (Z0 mass shifted??)
227 // Also, if problems, attempt several times to obtain closer mass match??
228 
229 bool MiniStringFragmentation::ministring2one( int iSub,
230  ColConfig& colConfig, Event& event) {
231 
232  // Cannot handle qq + qbarqbar system.
233  if (abs(flav1.id) > 100 && abs(flav2.id) > 100) return false;
234 
235  // For closed gluon loop need to pick an initial flavour.
236  if (isClosed) do {
237  int idStart = flavSelPtr->pickLightQ();
238  FlavContainer flavStart(idStart, 1);
239  flav1 = flavSelPtr->pick( flavStart);
240  flav2 = flav1.anti();
241  } while (abs(flav1.id) > 100);
242 
243  // Select hadron flavour from available quark flavours.
244  int idHad = 0;
245  for (int iTryFlav = 0; iTryFlav < NTRYFLAV; ++iTryFlav) {
246  idHad = flavSelPtr->combine( flav1, flav2);
247  if (idHad != 0) break;
248  }
249  if (idHad == 0) return false;
250 
251  // Find mass.
252  double mHad = particleDataPtr->mSel(idHad);
253 
254  // Find the untreated parton system which combines to the largest
255  // squared mass above mimimum required.
256  int iMax = -1;
257  double deltaM2 = mHad*mHad - mSum*mSum;
258  double delta2Max = 0.;
259  for (int iRec = iSub + 1; iRec < colConfig.size(); ++iRec) {
260  double delta2Rec = 2. * (pSum * colConfig[iRec].pSum) - deltaM2
261  - 2. * mHad * colConfig[iRec].mass;
262  if (delta2Rec > delta2Max) { iMax = iRec; delta2Max = delta2Rec;}
263  }
264  if (iMax == -1) return false;
265 
266  // Construct kinematics of the hadron and recoiling system.
267  Vec4& pRec = colConfig[iMax].pSum;
268  double mRec = colConfig[iMax].mass;
269  double vecProd = pSum * pRec;
270  double coefOld = mSum*mSum + vecProd;
271  double coefNew = mHad*mHad + vecProd;
272  double coefRec = mRec*mRec + vecProd;
273  double coefSum = coefOld + coefNew;
274  double sHat = coefOld + coefRec;
275  double root = sqrtpos( (pow2(coefSum) - 4. * sHat * mHad*mHad)
276  / (pow2(vecProd) - pow2(mSum * mRec)) );
277  double k2 = 0.5 * (coefOld * root - coefSum) / sHat;
278  double k1 = (coefRec * k2 + 0.5 * deltaM2) / coefOld;
279  Vec4 pHad = (1. + k1) * pSum - k2 * pRec;
280  Vec4 pRecNew = (1. + k2) * pRec - k1 * pSum;
281 
282  // Add the produced particle to the event record.
283  int iHad = event.append( idHad, 81, iParton.front(), iParton.back(),
284  0, 0, 0, 0, pHad, mHad);
285 
286  // Set decay vertex when this is displaced.
287  if (event[iParton.front()].hasVertex()) {
288  Vec4 vDec = event[iParton.front()].vDec();
289  event[iHad].vProd( vDec );
290  }
291 
292  // Set lifetime of hadron.
293  event[iHad].tau( event[iHad].tau0() * rndmPtr->exp() );
294 
295  // Mark original partons as hadronized and set their daughter range.
296  for (int i = 0; i < int(iParton.size()); ++i) {
297  event[ iParton[i] ].statusNeg();
298  event[ iParton[i] ].daughters(iHad, iHad);
299  }
300 
301  // Copy down recoiling system, with boosted momentum. Update current partons.
302  RotBstMatrix M;
303  M.bst(pRec, pRecNew);
304  for (int i = 0; i < colConfig[iMax].size(); ++i) {
305  int iOld = colConfig[iMax].iParton[i];
306  // Do not touch negative iOld = beginning of new junction leg.
307  if (iOld >= 0) {
308  int iNew = event.copy(iOld, 72);
309  event[iNew].rotbst(M);
310  colConfig[iMax].iParton[i] = iNew;
311  }
312  }
313  colConfig[iMax].pSum = pRecNew;
314  colConfig[iMax].isCollected = true;
315 
316  // Successfully done.
317  return true;
318 
319 }
320 
321 //==========================================================================
322 
323 } // end namespace Pythia8
Definition: AgUStep.h:26