StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ropewalk.cc
1 // Ropewalk.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 // RopeRandState, RopeDipoleEnd, OverlappingRopeDipole, RopeDipole, Ropewalk,
8 // RopeFragPars and FlavourRope classes.
9 
10 #include "Pythia8/Ropewalk.h"
11 
12 namespace Pythia8 {
13 
14 //==========================================================================
15 
16 // OverlappingRopeDipole class.
17 // This class describes dipoles overlapping with a given dipole.
18 
19 //--------------------------------------------------------------------------
20 
21 // Constructor sets up coordinates in other dipole's rest frame.
22 
23 OverlappingRopeDipole::OverlappingRopeDipole(RopeDipole* d, double m0,
24  RotBstMatrix& r) : dipole(d), dir(1) {
25 
26  // Coordinates in other dipole's rest frame
27  b1 = d->d1Ptr()->getParticlePtr()->vProd() * MM2FM;
28  b1.rotbst(r);
29  b2 = d->d2Ptr()->getParticlePtr()->vProd() * MM2FM;
30  b2.rotbst(r);
31  y1 = d->d1Ptr()->rap(m0,r);
32  y2 = d->d2Ptr()->rap(m0,r);
33  if (y1 < y2) dir = -1;
34 
35 }
36 
37 //--------------------------------------------------------------------------
38 
39 // Calculate the overlap at given y and b.
40 
41 bool OverlappingRopeDipole::overlap(double y, Vec4 ba, double r0) {
42 
43  if (y < min(y1, y2) || y > max(y1, y2)) return false;
44  Vec4 bb = b1 + (b2 - b1) * (y - y1) / (y2 - y1);
45  Vec4 tmp = ba - bb;
46  return (tmp.pT() <= 2 * r0);
47 
48 }
49 
50 //--------------------------------------------------------------------------
51 
52 // Has the dipole been hadronized?
53 
54 bool OverlappingRopeDipole::hadronized() {
55 
56  return dipole->hadronized();
57 
58 }
59 
60 //==========================================================================
61 
62 // RopeDipole class.
63 // This class describes a dipoles in impact parameter space.
64 
65 //--------------------------------------------------------------------------
66 
67 // The RopeDipole constructor makes sure that d1 is always the colored
68 // end and d2 the anti-colored.
69 
70 RopeDipole::RopeDipole(RopeDipoleEnd d1In, RopeDipoleEnd d2In, int iSubIn,
71  Info* infoPtrIn)
72  : d1(d1In), d2(d2In), iSub(iSubIn), hasRotFrom(false), hasRotTo(false),
73  isHadronized(false), infoPtr(infoPtrIn) {
74 
75  // Test if d1 is colored end and d2 anti-colored.
76  if (d1In.getParticlePtr()->col() == d2In.getParticlePtr()->acol()
77  && d1In.getParticlePtr()->col() != 0) { }
78  else { d2 = d1In, d1 = d2In; }
79 
80 }
81 
82 //--------------------------------------------------------------------------
83 
84 // Insert an excitation on dipole, if not already there.
85 
86 void RopeDipole::addExcitation(double ylab, Particle* ex) {
87 
88  pair<map<double, Particle*>::iterator, map<double, Particle*>::iterator>
89  ret = excitations.equal_range(ylab);
90  for (map<double, Particle*>::iterator itr = ret.first; itr != ret.second;
91  ++itr)
92  if (ex == itr->second) return;
93  excitations.insert( make_pair(ylab,ex) );
94 
95 }
96 
97 //--------------------------------------------------------------------------
98 
99 // Propagate the dipole itself.
100 
101 void RopeDipole::propagateInit(double deltat) {
102 
103  // Dipole end momenta.
104  Vec4 pcm = d1.getParticlePtr()->p();
105  Vec4 pam = d2.getParticlePtr()->p();
106  double mT2c = pcm.pT2() + pcm.m2Calc();
107  double mT2a = pam.pT2() + pam.m2Calc();
108  if (mT2c <= 0 || mT2a <= 0) {
109  infoPtr->errorMsg("Error in RopeDipole::propagateInit: Tried to"
110  "propagate a RopeDipoleEnd with mT <= 0");
111  return;
112  }
113  double mTc = sqrt(mT2c);
114  double mTa = sqrt(mT2a);
115  // New vertices in the lab frame.
116  Vec4 newv1 = Vec4( deltat * pcm.px() / mTc, deltat * pcm.py() / mTc, 0, 0);
117  Vec4 newv2 = Vec4( deltat * pam.px() / mTa, deltat * pam.py() / mTa, 0, 0);
118 
119  // Set the new vertices deep.
120  d1.getParticlePtr()->vProdAdd(newv1 * FM2MM);
121  d2.getParticlePtr()->vProdAdd(newv2 * FM2MM);
122 
123 }
124 
125 //--------------------------------------------------------------------------
126 
127 // Propagate both dipole ends as well as all excitations.
128 
129 void RopeDipole::propagate(double deltat, double m0) {
130 
131  // First propagate the dipole ends.
132  propagateInit(deltat);
133  for (map<double, Particle*>::iterator eItr = excitations.begin();
134  eItr != excitations.end(); ++eItr) {
135 
136  Vec4 em = eItr->second->p();
137  em.rotbst(getDipoleLabFrame());
138  // Propagate excitations.
139 
140  if (em.pT() > 0.0){
141  Vec4 newVert = Vec4( deltat * em.px() / em.pT(),
142  deltat * em.py() / em.pT(), 0, 0);
143  eItr->second->vProdAdd(newVert * FM2MM);
144  }
145  else eItr->second->vProd( bInterpolateLab(eItr->first,m0) * FM2MM);
146  }
147 
148 }
149 
150 //--------------------------------------------------------------------------
151 
152 // Put gluon excitations on the dipole.
153 
154 void RopeDipole::excitationsToString(double m0, Event& event) {
155 
156  // Erase excitations below cut-off.
157  map<double, Particle*>::iterator pItr = excitations.begin();
158  while (pItr != excitations.end() ) {
159  if (pItr->second->pAbs() < 1e-6) {
160  map<double, Particle*>::iterator eraseMe = pItr;
161  ++pItr;
162  excitations.erase(eraseMe);
163  }
164  else ++pItr;
165  }
166 
167  // We now colour connect the excitations to the dipole.
168  // The dipole is connected like this sketch:
169  // acol (d1) col ------ acol (d2) col.
170  int oldcol = d1.getParticlePtr()->col();
171  if (oldcol != d2.getParticlePtr()->acol()) {
172  infoPtr->errorMsg("Error in Ropewalk::RopeDipole::excitationsToString: "
173  "color indices do not match.");
174  return;
175  }
176  vector<int> daughters;
177 
178  // Two cases depending on which end we should start at.
179  // We always start at min rapidity and connect from there.
180  if (d1.rap(m0) == minRapidity(m0)) {
181  int acol = oldcol;
182  for (map<double, Particle*>::iterator itr = excitations.begin();
183  itr != excitations.end(); ++itr) {
184  int col = event.nextColTag();
185  itr->second->status(51);
186  itr->second->mothers(d1.getNe(),d1.getNe());
187  itr->second->cols(col,acol);
188  daughters.push_back(event.append(Particle(*(itr->second))));
189  acol = col;
190  }
191  d2.getParticlePtr()->acol(acol);
192  event[d2.getNe()].acol(acol);
193  }
194  else {
195  int acol = oldcol;
196  for (map<double, Particle*>::reverse_iterator itr = excitations.rbegin();
197  itr != excitations.rend(); ++itr) {
198  int col = event.nextColTag();
199  itr->second->status(51);
200  itr->second->mothers(d1.getNe(),d1.getNe());
201  itr->second->cols(col,acol);
202  daughters.push_back(event.append(Particle(*(itr->second))));
203  acol = col;
204  }
205  d2.getParticlePtr()->acol(acol);
206  event[d2.getNe()].acol(acol);
207  }
208  bool stringEnd = false;
209  if (d2.getParticlePtr()->col() == 0) stringEnd = true;
210 
211  // Update status codes and mother/daughter indexing.
212  event[d1.getNe()].statusNeg();
213  Particle cc1 = *d1.getParticlePtr();
214  cc1.statusPos();
215  cc1.mothers(d1.getNe(),d1.getNe());
216  daughters.push_back(event.append(cc1));
217  event[d1.getNe()].daughters( daughters[0], daughters[daughters.size() -1 ] );
218  if (stringEnd) {
219  event[d2.getNe()].statusNeg();
220  Particle cc2 = *d2.getParticlePtr();
221  cc2.statusPos();
222  cc2.mothers(d2.getNe(),d2.getNe());
223  int did = event.append(cc2);
224  event[d2.getNe()].daughters(did,did);
225  }
226 
227 }
228 
229 //--------------------------------------------------------------------------
230 
231 // Redistribute momentum to two particles.
232 
233 void RopeDipole::splitMomentum(Vec4 mom, Particle* p1, Particle* p2,
234  double frac) {
235 
236  Vec4 p1new = p1->p() + frac * mom;
237  Vec4 p2new = p2->p() + (1. - frac) * mom;
238  p1->p(p1new);
239  p2->p(p2new);
240  return;
241 
242 }
243 
244 //--------------------------------------------------------------------------
245 
246 // Recoil the dipole from adding a gluon. If the "dummy" option is set,
247 // the recoil will not be added, but only checked.
248 // Note: the gluon will not actually be added, only the recoil (if possible).
249 
250 bool RopeDipole::recoil(Vec4& pg, bool dummy) {
251 
252  // Keep track of direction.
253  int sign = 1;
254  if (d1.rap(1.0) > d2.rap(1.0)) sign = -1;
255 
256  // Lightcone momenta after inserting the gluon.
257  Particle* epaPtr = d1.getParticlePtr();
258  Particle* epcPtr = d2.getParticlePtr();
259  double pplus = epcPtr->pPos() + epaPtr->pPos() - pg.pPos();
260  double pminus = epcPtr->pNeg() + epaPtr->pNeg() - pg.pNeg();
261 
262  // The new lightcone momenta of the dipole ends.
263  double ppa = 0.0;
264  double ppc = 0.0;
265  double pma = 0.0;
266  double pmc = 0.0;
267  double mta2 = epaPtr->mT2();
268  double mtc2 = epcPtr->mT2();
269  double mta = sqrt(mta2);
270  double mtc = sqrt(mtc2);
271  if ( pplus * pminus <= pow2(mta + mtc)
272  || pplus <= 0.0 || pminus <= 0.0 ) return false;
273 
274  // Calculate the new momenta.
275  double sqarg = pow2(pplus * pminus - mta2 - mtc2) - 4. * mta2 * mtc2;
276  if (sqarg <= 0.0 ) return false;
277  if (sign > 0) {
278  ppa = 0.5 * (pplus * pminus + mta2 - mtc2 + sqrt(sqarg)) / pminus;
279  pma = mta2 / ppa;
280  pmc = pminus - pma;
281  ppc = mtc2 / pmc;
282  // Check rapidity after recoil.
283  if ( ppa * mtc < ppc * mta ) return false;
284  }
285  else {
286  pma = 0.5 * (pplus * pminus + mta2 - mtc2 + sqrt(sqarg)) / pplus;
287  ppa = mta2 / pma;
288  ppc = pplus - ppa;
289  pmc = mtc2 / ppc;
290  // Check rapidity after recoil.
291  if ( ppa*mtc > ppc*mta ) return false;
292  }
293 
294  // Set up and store the new momenta.
295  Vec4 shifta = Vec4( epaPtr->px(), epaPtr->py(),
296  0.5 * (ppa - pma), 0.5 * (ppa + pma));
297  Vec4 shiftc = Vec4( epcPtr->px(), epcPtr->py(),
298  0.5 * (ppc - pmc), 0.5 * (ppc + pmc));
299  if (!dummy) {
300  epaPtr->p(shifta);
301  epcPtr->p(shiftc);
302  }
303  return true;
304 
305 }
306 
307 //--------------------------------------------------------------------------
308 
309 // Get the Lorentz matrix to go to the dipole rest frame.
310 
311 RotBstMatrix RopeDipole::getDipoleRestFrame() {
312 
313  if (hasRotTo) return rotTo;
314 
315  RotBstMatrix r;
316  r.toCMframe(d1.getParticlePtr()->p(),d2.getParticlePtr()->p());
317  rotTo = r;
318  hasRotTo = true;
319  return rotTo;
320 
321 }
322 
323 //--------------------------------------------------------------------------
324 
325 // Get the Lorentz matrix to go from the dipole rest frame to lab frame.
326 
327 RotBstMatrix RopeDipole::getDipoleLabFrame() {
328 
329  if(hasRotFrom) return rotFrom;
330 
331  RotBstMatrix r;
332  r.fromCMframe(d1.getParticlePtr()->p(),d2.getParticlePtr()->p());
333  rotFrom = r;
334  hasRotFrom = true;
335  return rotFrom;
336 
337 }
338 //--------------------------------------------------------------------------
339 
340 // The dipole four-momentum.
341 
342 Vec4 RopeDipole::dipoleMomentum() {
343 
344  Vec4 ret = d1.getParticlePtr()->p() + d2.getParticlePtr()->p();
345  return ret;
346 
347 }
348 
349 //--------------------------------------------------------------------------
350 
351 // Interpolate (linear) between dipole ends to get b position at given y.
352 // Here y must be given in dipole rest frame, and the resulting b-position
353 // will also be in the dipole rest frame.
354 
355 Vec4 RopeDipole::bInterpolateDip(double y, double m0) {
356  if(!hasRotTo) getDipoleRestFrame();
357  Vec4 bb1 = d1.getParticlePtr()->vProd() * MM2FM;
358  bb1.rotbst(rotTo);
359  Vec4 bb2 = d2.getParticlePtr()->vProd() * MM2FM;
360  bb2.rotbst(rotTo);
361  double y1 = d1.rap(m0,rotTo);
362  double y2 = d2.rap(m0,rotTo);
363  return bb1 + y * (bb2 - bb1) / (y2 - y1);
364 
365 }
366 
367 //--------------------------------------------------------------------------
368 
369 // Interpolate (linear) between dipole ends to get b position at given y.
370 
371 Vec4 RopeDipole::bInterpolateLab(double y, double m0) {
372 
373  Vec4 bb1 = d1.getParticlePtr()->vProd() * MM2FM;
374  Vec4 bb2 = d2.getParticlePtr()->vProd() * MM2FM;
375  double y1 = d1.rap(m0);
376  double y2 = d2.rap(m0);
377  return bb1 + y * (bb2 - bb1) / (y2 - y1);
378 
379 }
380 
381 //--------------------------------------------------------------------------
382 
383 // Interpolate (linear) between dipole ends to get b position in
384 // a given frame at given y.
385 
386 Vec4 RopeDipole::bInterpolate(double y, RotBstMatrix rb, double m0) {
387 
388  Vec4 bb1 = d1.getParticlePtr()->vProd() * MM2FM;
389  Vec4 bb2 = d2.getParticlePtr()->vProd() * MM2FM;
390  bb1.rotbst(rb);
391  bb2.rotbst(rb);
392  double y1 = d1.rap(m0);
393  double y2 = d2.rap(m0);
394  return bb1 + y * (bb2 - bb1) / (y2 - y1);
395 
396 }
397 
398 //--------------------------------------------------------------------------
399 
400 // Calculate the amount of overlapping dipoles at a given rapidity.
401 // Return the number of overlapping dipoles to the "left" and "right".
402 
403 pair<int, int> RopeDipole::getOverlaps(double yfrac, double m0, double r0) {
404  // Transform yfrac to y in dipole rest frame
405  if (!hasRotTo) getDipoleRestFrame();
406  double yL = d1.rap(m0,rotTo);
407  double yS = d2.rap(m0,rotTo);
408  double yH = yS + (yL - yS) * yfrac;
409  int m = 0, n = 0;
410  for (size_t i = 0; i < overlaps.size(); ++i) {
411  if (overlaps[i].overlap( yfrac, bInterpolateDip(yH,m0), r0)
412  && !overlaps[i].hadronized()) {
413  if (overlaps[i].dir > 0) ++m;
414  else ++n;
415  }
416  }
417  return make_pair(m,n);
418 
419 }
420 
421 //==========================================================================
422 
423 // Exc class.
424 // It is a helper class to Ropewalk, used to describe a pair of excitations
425 // needed for shoving. It is kept away from users, as there are a lot of
426 // raw pointers floating around.
427 
428 //--------------------------------------------------------------------------
429 
430 struct Exc {
431 
432 // The constructor.
433 Exc(double yIn, double m0In, int iIn, int jIn, int kIn, RopeDipole* dip1In,
434  RopeDipole* dip2In) : y(yIn), m0(m0In), i(iIn), j(jIn), k(kIn), pp1(NULL),
435  pp2(NULL), dip1(dip1In), dip2(dip2In) { }
436 
437 // Set pointers to the two excitation particles.
438 void setParticlePtrs(Particle* p1, Particle* p2) {
439  pp1 = p1;
440  pp2 = p2;
441  // Set a pointer to the excitation in the dipole.
442  dip1->addExcitation(y, pp1);
443  dip2->addExcitation(y, pp2);
444 }
445 
446 // Give the excitation a kick in the x and y direction,
447 void shove(double dpx, double dpy) {
448  // The old momenta.
449  Vec4 p2 = pp2->p();
450  Vec4 p1 = pp1->p();
451  // The new momenta, after the shove.
452  double mt2new = sqrt(pow2(p2.px() - dpx) + pow2(p2.py() - dpy));
453  double e2new = mt2new * cosh(y);
454  double p2znew = mt2new * sinh(y);
455  Vec4 p2new(p2.px() - dpx, p2.py() - dpy, p2znew, e2new);
456  double mt1new = sqrt(pow2(p1.px() + dpx) + pow2(p1.py() + dpy));
457  double e1new = mt1new * cosh(y);
458  double p1znew = mt1new * sinh(y);
459  Vec4 p1new(p1.px() + dpx, p1.py() + dpy, p1znew, e1new);
460  // The differences between the two.
461  Vec4 deltap1 = p1new - p1;
462  Vec4 deltap2 = p2new - p2;
463  // Now check if we can add these two excitations to the dipoles.
464  if ( dip2->recoil(deltap2) ) {
465  if ( dip1->recoil(deltap1) ) {
466  pp1->p(p1new);
467  pp2->p(p2new);
468  } else {
469  Vec4 dp2 = -deltap2;
470  dip2->recoil(dp2);
471  }
472  }
473 }
474 
475 // The push direction as a four vector.
476 Vec4 direction() {
477  return dip1->bInterpolateDip(y,m0) -
478  dip2->bInterpolateDip(y,m0);
479 }
480 
481 // Member variables, slice rapidity and small cut-off mass.
482 double y;
483 double m0;
484 
485 // Local particle indices.
486 int i, j, k;
487 Particle* pp1;
488 Particle* pp2;
489 RopeDipole* dip1;
490 RopeDipole* dip2;
491 };
492 
493 //==========================================================================
494 
495 // Ropewalk class.
496 // This class keeps track of all the strings making up ropes for shoving
497 // as well as flavour enhancement.
498 
499 //--------------------------------------------------------------------------
500 
501 // The Ropewalk init function sets parameters and pointers.
502 
503 bool Ropewalk::init() {
504 
505  StringInteractions::init();
506 
507  // Parameters of the ropewalk.
508  shoveMiniStrings = flag("Ropewalk:shoveMiniStrings");
509  shoveJunctionStrings = flag("Ropewalk:shoveJunctionStrings");
510  shoveGluonLoops = flag("Ropewalk:shoveGluonLoops");
511  limitMom = flag("Ropewalk:limitMom");
512  mStringMin = parm("HadronLevel:mStringMin");
513  r0 = parm("Ropewalk:r0");
514  m0 = parm("Ropewalk:m0");
515  pTcut = parm("Ropewalk:pTcut");
516  rCutOff = parm("Ropewalk:rCutOff");
517  gAmplitude = parm("Ropewalk:gAmplitude");
518  gExponent = parm("Ropewalk:gExponent");
519  deltay = parm("Ropewalk:deltay");
520  deltat = parm("Ropewalk:deltat");
521  tShove = parm("Ropewalk:tShove");
522  tInit = parm("Ropewalk:tInit");
523  showerCut = parm("TimeShower:pTmin");
524  alwaysHighest = flag("Ropewalk:alwaysHighest");
525 
526  // Creat the interface objects.
527  if ( flag("Ropewalk:doShoving") ) {
528  // Check consistency.
529  if (deltat > tShove) {
530  infoPtr->errorMsg("Error in Ropewalk::init: "
531  "deltat cannot be larger than tShove");
532  return false;
533  }
534  if ( !flag("PartonVertex:setVertex") ) {
535  infoPtr->errorMsg("Error in Ropewalk::init: "
536  "Shoving enabled, but no vertex information.");
537  return false;
538  }
539  stringrepPtr = make_shared<RopewalkShover>(*this);
540  registerSubObject(*stringrepPtr);
541  if ( !stringrepPtr->init() ) return false;
542  }
543 
544  if ( flag("Ropewalk:doFlavour") ) {
545  // Sanity check of flavour rope and vertex information.
546  // Flavour ropes requires vertex information, unless an effective
547  // string tension is supplied by hand or Buffon treatment.
548  if (!flag("PartonVertex:setVertex") &&
549  (!flag("Ropewalk:setFixedKappa") &&
550  !flag("Ropewalk:doBuffon")) ) {
551  infoPtr->errorMsg("Error in Ropewalk::init: "
552  "failed initialization of flavour ropes");
553  return false;
554  }
555  fragmodPtr = make_shared<FlavourRope>(*this);
556  registerSubObject(*fragmodPtr);
557  if ( !fragmodPtr->init() ) return false;
558  }
559 
560  return true;
561 
562 }
563 
564 //--------------------------------------------------------------------------
565 
566 // Calculate the average string tension of the event, in units of the default
567 // string tension (ie. 1 GeV/fm), using random walk in colour space.
568 
569 double Ropewalk::averageKappa() {
570 
571  double kap = 0.;
572  double nd = 0;
573  for (DMap::iterator itr = dipoles.begin(); itr != dipoles.end(); ++itr) {
574 
575  // Getting the overlaps: m, n.
576  pair<int,int> overlap = itr->second.getOverlaps( rndmPtr->flat(), m0, r0);
577 
578  // Overlaps define the number of steps taken in the random walk.
579  // We define the present dipole to always point in the p-direction.
580  pair<double, double> pq = select( overlap.first + 1, overlap.second,
581  rndmPtr);
582  double enh = 0.25 * (2. + 2. * pq.first + pq.second);
583  kap += (enh > 1.0 ? enh : 1.0);
584  nd += 1.0;
585  }
586  return kap / nd;
587 
588 }
589 
590 //--------------------------------------------------------------------------
591 
592 // Calculate the effective string tension a fraction yfrac in on the dipole,
593 // given by indices e1 and e2.
594 
595 double Ropewalk::getKappaHere(int e1, int e2, double yfrac) {
596 
597  multimap< pair<int,int>, RopeDipole >::iterator
598  itr = dipoles.find( make_pair(e1,e2) );
599  if (itr == dipoles.end()) itr = dipoles.find( make_pair(e2,e1) );
600  if (itr == dipoles.end()) return 1.0;
601  RopeDipole* d = &(itr->second);
602  d->hadronized(true);
603 
604  // Get quantum numbers m and n.
605  pair<int, int> overlap = d->getOverlaps(yfrac, m0, r0);
606  pair<double, double> pq;
607  // If we are always in the highest multiplet, we need not do
608  // a random walk
609  if (alwaysHighest) {
610  pq = make_pair(overlap.first + 1, overlap.second);
611  }
612  // Invoke default random walk procedure.
613  else {
614  pq = select(overlap.first + 1, overlap.second, rndmPtr);
615  }
616  // Calculate enhancement factor.
617  double enh = 0.25 * (2. + 2. * pq.first + pq.second);
618  return (enh > 1.0 ? enh : 1.0);
619 
620 }
621 
622 //--------------------------------------------------------------------------
623 
624 // Calculate all overlaps of all dipoles and store as OverlappingRopeDipoles.
625 
626 bool Ropewalk::calculateOverlaps() {
627 
628  // Go through all dipoles.
629  for (multimap< pair<int,int>, RopeDipole>::iterator itr = dipoles.begin();
630  itr != dipoles.end(); ++itr ) {
631  RopeDipole* d1 = &(itr->second);
632  if (d1->dipoleMomentum().m2Calc() < pow2(m0)) continue;
633 
634  // RopeDipoles rapidities in dipole rest frame.
635  RotBstMatrix dipoleRestFrame = d1->getDipoleRestFrame();
636  double yc1 = d1->d1Ptr()->rap(m0, dipoleRestFrame);
637  double ya1 = d1->d2Ptr()->rap(m0, dipoleRestFrame);
638  if (yc1 <= ya1) continue;
639 
640  // Go through all possible overlapping dipoles.
641  for (multimap< pair<int,int>, RopeDipole>::iterator itr2 = dipoles.begin();
642  itr2 != dipoles.end(); ++itr2) {
643  RopeDipole* d2 = &(itr2->second);
644 
645  // Skip self and overlaps with miniscule dipoles.
646  if (d1 == d2) continue;
647  if (d2->dipoleMomentum().m2Calc() < pow2(m0)) continue;
648 
649  // Ignore if not overlapping in rapidity.
650  OverlappingRopeDipole od(d2, m0, dipoleRestFrame);
651  if (min(od.y1, od.y2) > yc1 || max(od.y1, od.y2) < ya1 || od.y1 == od.y2)
652  continue;
653 
654  d1->addOverlappingDipole(od);
655 
656  }
657  }
658  return true;
659 
660 }
661 //--------------------------------------------------------------------------
662 
663 // Invoke the random walk of colour states.
664 
665 pair<int, int> Ropewalk::select(int m, int n, Rndm* rndm) {
666 
667  // Initial valuesM mm and nn are step counters.
668  int p = 0, q = 0;
669  int mm = m, nn = n;
670 
671  // We must take all steps before terminating.
672  while (mm + nn > 0) {
673 
674  // Take randomly a step in one direction.
675  if (rndm->flat() < 0.5 && mm > 0) {
676  --mm;
677 
678  // Calculate the step probabilities.
679  double p1 = multiplicity(p + 1, q);
680  double p2 = multiplicity(p, q - 1);
681  double p3 = multiplicity(p - 1, q + 1);
682 
683  // Normalization.
684  double sum = p1 + p2 + p3;
685  p1 /= sum, p2 /= sum, p3 /= sum;
686 
687  // Select a state.
688  double pick = rndm->flat();
689  if (pick < p1) ++p;
690  else if (pick < p1 + p2) --q;
691  else --p, ++q;
692  }
693 
694  // As above, but for nn.
695  else if (nn > 0) {
696  --nn;
697  double p1 = multiplicity(p, q + 1);
698  double p2 = multiplicity(p -1, q);
699  double p3 = multiplicity(p + 1, q - 1);
700  double sum = p1 + p2 + p3;
701  p1 /= sum, p2 /= sum, p3 /= sum;
702  double pick = rndm->flat();
703  if (pick < p1) ++q;
704  else if (pick < p1 + p2) --p;
705  else --q, ++p;
706  }
707  }
708 
709  // Done.
710  return make_pair( (p < 0 ? 0 : p), (q < 0 ? 0 : q) );
711 
712 }
713 
714 //--------------------------------------------------------------------------
715 
716 // Shove all dipoles in the event.
717 
718 void Ropewalk::shoveTheDipoles(Event& event) {
719 
720  // Possibility of some initial propagation.
721  if ( tInit > 0.0) {
722  for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr)
723  dItr->second.propagateInit(tInit);
724  }
725 
726  // The rapidity slices.
727  multimap<double, RopeDipole *> rapDipoles;
728 
729  // Order the dipoles in max rapidity.
730  double ymin = 0;
731  double ymax = 0;
732  for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr) {
733  RopeDipole* dip = &(dItr->second);
734  // Order the dipoles in max rapidity.
735  rapDipoles.insert( make_pair(dip->maxRapidity(m0), dip) );
736  // Find maximal and minimal rapidity to sample.
737  if (dip->minRapidity(m0) < ymin) ymin = dip->minRapidity(m0);
738  if (dip->maxRapidity(m0) > ymax) ymax = dip->maxRapidity(m0);
739  }
740 
741  // Do the sampling from flat distribution.
742  vector<double> rapidities;
743  for (double y = ymin; y < ymax; y += deltay) rapidities.push_back(y);
744 
745  // For each value of ySample, we have a vector of excitation pairs.
746  map<double, vector<Exc> > exPairs;
747  for (int i = 0, N = eParticles.size(); i < N; ++i) eParticles[i].clear();
748  eParticles.clear();
749  for (int i = 0, N = rapidities.size(); i < N; ++i) {
750  // Construct an empty vector of excitation particles.
751  eParticles.push_back( vector<Particle>() );
752 
753  // Find dipoles sampled in this slice, and store them temporarily.
754  double ySample = rapidities[i];
755  vector<RopeDipole*> tmp;
756  for (multimap<double, RopeDipole*>::iterator
757  rItr = rapDipoles.lower_bound(ySample); rItr != rapDipoles.end(); ++rItr) {
758  if (rItr->second->minRapidity(m0) < ySample)
759  tmp.push_back(rItr->second);
760  }
761 
762  // Construct excitation particles, one for each sampled dipole in this slice.
763  vector<int> eraseDipoles;
764 
765  for (int j = 0, M = tmp.size(); j < M; ++j) {
766  Vec4 ex;
767  // Test if the dipole can bear the excitation.
768  if (!tmp[j]->recoil(ex,true) ) {
769  eraseDipoles.push_back(j);
770  }
771  }
772  // Erase dipoles which could not bear an excitation.
773  for (int j = 0, M = eraseDipoles.size(); j < M; ++j) {
774  tmp.erase( tmp.begin() + (eraseDipoles[j]-j) );
775  }
776  // Add the actual excitations, but only if we can create pairs.
777  if( int(tmp.size()) > 1)
778  for (int j = 0, M = tmp.size(); j < M; ++j) {
779  Vec4 ex;
780  // We boost the excitation back from dipole rest frame.
781  tmp[j]->recoil(ex,false);
782  Particle pp = Particle(21, 22, 0, 0, 0, 0, 0, 0, ex);
783  pp.vProd( tmp[j]->bInterpolateLab(ySample,m0) * FM2MM);
784  eParticles[i].push_back(pp);
785  }
786  // Construct all pairs of possible excitations in this slice.
787  exPairs[ySample] = vector<Exc>();
788  for (int j = 0, M = tmp.size(); j < M; ++j)
789  for (int k = 0; k < M; ++k) {
790  // Don't allow a string to shove itself.
791  if(j != k && tmp[j]->index() != tmp[k]->index() )
792  exPairs[ySample].push_back( Exc(ySample, m0, i, j, k, tmp[j],
793  tmp[k]) );
794  }
795  }
796 
797  // Give the excitations pointers to the excitation particles.
798  for (map<double, vector<Exc> >::iterator slItr = exPairs.begin();
799  slItr != exPairs.end(); ++slItr) {
800  for (int i = 0, N = slItr->second.size(); i < N; ++i) {
801  Exc& ep = slItr->second[i];
802  ep.setParticlePtrs( &eParticles[ep.i][ep.j], &eParticles[ep.i][ep.k] );
803  }
804  }
805 
806  // Shoving loop.
807  for (double t = tInit; t < tShove + tInit; t += deltat) {
808  // For all slices.
809  for (map<double, vector<Exc> >::iterator slItr = exPairs.begin();
810  slItr != exPairs.end(); ++slItr)
811  // For all excitation pairs.
812  for (int i = 0, N = slItr->second.size(); i < N; ++i) {
813  Exc& ep = slItr->second[i];
814  // The direction vector is a space-time four-vector.
815  Vec4 direction = ep.direction();
816  // The string radius is time dependent,
817  // growing with the speed of light.
818  // Minimal string size is 1 / shower cut-off
819  // converted to fm.
820  double rt = max(t, 1. / showerCut / 5.068);
821  rt = min(rt, r0 * gExponent);
822  double dist = direction.pT();
823  // Calculate the push, its direction and do the shoving.
824  if (dist < rCutOff * rt) {
825  // Gain function.
826  double gain = 0.5 * deltay * deltat * gAmplitude * dist / rt / rt
827  * exp( -0.25 * dist * dist / rt / rt);
828  double dpx = dist > 0.0 ? gain * direction.px() / dist: 0.0;
829  double dpy = dist > 0.0 ? gain * direction.py() / dist: 0.0;
830  ep.shove(dpx, dpy);
831  }
832  }
833 
834  // Propagate the dipoles.
835  for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr)
836  dItr->second.propagate(deltat, m0);
837  }
838 
839  // Add the excitations to the dipoles.
840  for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr) {
841  RopeDipole* dip = &(dItr->second);
842  if (dip->nExcitations() > 0) dip->excitationsToString( m0, event);
843  }
844 }
845 
846 //--------------------------------------------------------------------------
847 
848 // Extract all dipoles from an event.
849 
850 bool Ropewalk::extractDipoles(Event& event, ColConfig& colConfig) {
851 
852  // Go through all strings in the event.
853  dipoles.clear();
854  for (int iSub = 0; iSub < colConfig.size(); ++iSub) {
855 
856  // We can exclude loops, junctions and ministrings from the Ropewalk.
857  if (colConfig[iSub].hasJunction && !shoveJunctionStrings) continue;
858  if (colConfig[iSub].isClosed && !shoveGluonLoops) continue;
859  if (colConfig[iSub].massExcess <= mStringMin && !shoveMiniStrings)
860  continue;
861 
862  colConfig.collect(iSub,event);
863  vector<int> stringPartons = colConfig[iSub].iParton;
864  vector<RopeDipole> stringDipole;
865  bool stringStart = true;
866  RopeDipoleEnd previous;
867  for (int iPar = int(stringPartons.size() - 1); iPar > -1; --iPar) {
868  if (stringPartons[iPar] > 0) {
869  // Ordinary particle.
870  RopeDipoleEnd next( &event, stringPartons[iPar]);
871  // If we are at the first parton, no dipole.
872  if ( !stringStart) {
873  // Get the parton placement in Event Record.
874  pair<int,int> dipoleER = make_pair( stringPartons[iPar + 1],
875  stringPartons[iPar] );
876  RopeDipole test(previous, next, iSub, infoPtr);
877  if ( limitMom && test.dipoleMomentum().pT() < pTcut)
878  dipoles.insert( pair< pair<int, int>, RopeDipole>(dipoleER,
879  RopeDipole( previous, next, iSub, infoPtr)) );
880  else if (!limitMom)
881  dipoles.insert( pair< pair<int, int>, RopeDipole>(dipoleER,
882  RopeDipole( previous, next, iSub, infoPtr)));
883  }
884  previous = next;
885  stringStart = false;
886  }
887  else continue;
888  }
889  // We have created all dipoles.
890  }
891  return true;
892 
893 }
894 
895 //==========================================================================
896 
897 // RopeFragPars recalculates fragmentation parameters according to a
898 // changed string tension. Helper class to FlavourRope.
899 
900 //--------------------------------------------------------------------------
901 
902 // Constants: could be changed here if desired, but normally should not.
903 
904 // Initial step size for a calculation.
905 const double RopeFragPars::DELTAA = 0.1;
906 
907 // Convergence criterion for a calculation.
908 const double RopeFragPars::ACONV = 0.001;
909 
910 // Low z cut-off in fragmentation function.
911 const double RopeFragPars::ZCUT = 1.0e-4;
912 
913 //--------------------------------------------------------------------------
914 
915 // The init function sets up initial parameters from settings.
916 
917 bool RopeFragPars::init() {
918 
919  // The junction parameter.
920  beta = parm("Ropewalk:beta");
921 
922  // Initialize default values from input settings.
923  const int len = 9;
924  string params [len] = {"StringPT:sigma", "StringZ:aLund",
925  "StringZ:aExtraDiquark","StringZ:bLund", "StringFlav:probStoUD",
926  "StringFlav:probSQtoQQ", "StringFlav:probQQ1toQQ0", "StringFlav:probQQtoQ",
927  "StringFlav:kappa"};
928  double* variables[len] = {&sigmaIn, &aIn, &adiqIn, &bIn, &rhoIn, &xIn,
929  &yIn, &xiIn, &kappaIn};
930  for (int i = 0; i < len; ++i) *variables[i] = parm(params[i]);
931 
932  // Insert the h = 1 case immediately.
933  sigmaEff = sigmaIn, aEff = aIn, adiqEff = adiqIn, bEff = bIn,
934  rhoEff = rhoIn, xEff = xIn, yEff = yIn, xiEff = xiIn, kappaEff = kappaIn;
935  if (!insertEffectiveParameters(1.0)) { infoPtr->errorMsg(
936  "Error in RopeFragPars::init: failed to insert defaults.");
937  return false;
938  }
939 
940  return true;
941 
942 }
943 
944 //--------------------------------------------------------------------------
945 
946 // Return parameters at given string tension, ordered by their
947 // name for easy insertion in settings.
948 
949 map<string,double> RopeFragPars::getEffectiveParameters(double h) {
950 
951  map<double, map<string, double> >::iterator parItr = parameters.find(h);
952 
953  // If the parameters are already calculated, return them.
954  if ( parItr != parameters.end()) return parItr->second;
955 
956  // Otherwise calculate them.
957  if (!calculateEffectiveParameters(h))
958  infoPtr->errorMsg("Error in RopeFragPars::getEffectiveParameters:"
959  " calculating effective parameters.");
960 
961  // And insert them.
962  if (!insertEffectiveParameters(h))
963  infoPtr->errorMsg("Error in RopeFragPars::getEffectiveParameters:"
964  " inserting effective parameters.");
965 
966  // And recurse.
967  return getEffectiveParameters(h);
968 
969 }
970 
971 //--------------------------------------------------------------------------
972 
973 // Get the Fragmentation function a parameter from cache or calculate it.
974 
975 double RopeFragPars::getEffectiveA(double thisb, double mT2, bool isDiquark) {
976 
977  // Check for the trivial case.
978  if (thisb == bIn) return (isDiquark ? aIn + adiqIn : aIn);
979 
980  // We order by b*mT2
981  map<double, double>* aMPtr = (isDiquark ? &aDiqMap : &aMap);
982  double bmT2 = mT2 * thisb;
983 
984  // Check if we have already calculated this a value before.
985  map<double,double>::iterator aItr = aMPtr->find(bmT2);
986  if (aItr != aMPtr->end()) return aItr->second;
987 
988  // Otherwise calculate it.
989  double ae = ( isDiquark ? aEffective(aIn + adiqIn, thisb, mT2)
990  : aEffective(aIn, thisb, mT2) );
991  if (isDiquark) {
992  double suba = getEffectiveA(thisb, mT2, false);
993  aMPtr->insert( make_pair(bmT2, ae - suba) );
994  }
995  else aMPtr->insert(make_pair(bmT2, ae));
996  return ae;
997 
998 }
999 
1000 //--------------------------------------------------------------------------
1001 
1002 // Calculate the effective parameters.
1003 
1004 bool RopeFragPars::calculateEffectiveParameters(double h) {
1005 
1006  if (h <= 0) return false;
1007  double hinv = 1.0 / h;
1008 
1009  // Start with the easiest transformations.
1010  // The string tension kappa.
1011  kappaEff = kappaIn * h;
1012  // Strangeness.
1013  rhoEff = pow(rhoIn, hinv);
1014  // Strange diquarks.
1015  xEff = pow(xIn, hinv);
1016  // Spin.
1017  yEff = pow(yIn, hinv);
1018  // pT distribution width.
1019  sigmaEff = sigmaIn * sqrt(h);
1020 
1021  // Derived quantity alpha.
1022  double alpha = (1 + 2. * xIn * rhoIn + 9. * yIn + 6. * xIn * rhoIn * yIn
1023  + 3. * yIn * xIn * xIn * rhoIn * rhoIn) / (2. + rhoIn);
1024  double alphaEff = (1. + 2. * xEff * rhoEff + 9. * yEff
1025  + 6. * xEff * rhoEff * yEff + 3. * yEff * xEff * xEff * rhoEff * rhoEff)
1026  / (2. + rhoEff);
1027 
1028  // Baryons.
1029  xiEff = alphaEff * beta * pow( xiIn / alpha / beta, hinv);
1030  if (xiEff > 1.0) xiEff = 1.0;
1031  if (xiEff < xiIn) xiEff = xiIn;
1032 
1033  // Fragmentation function b.
1034  bEff = (2. + rhoEff) / (2. + rhoIn) * bIn;
1035  if (bEff < bIn) bEff = bIn;
1036  if (bEff > 2.0) bEff = 2.0;
1037 
1038  // We calculate a for a typical particle with mT2 = 1 GeV^2.
1039  aEff = getEffectiveA( bEff, 1.0, false);
1040  adiqEff = getEffectiveA( bEff, 1.0, true) - aEff;
1041 
1042  return true;
1043 
1044 }
1045 
1046 //--------------------------------------------------------------------------
1047 
1048 // Insert calculated parameters in cache for later (re-)use.
1049 
1050 bool RopeFragPars::insertEffectiveParameters(double h) {
1051 
1052  map<string,double> p;
1053  p["StringPT:sigma"] = sigmaEff;
1054  p["StringZ:bLund"] = bEff;
1055  p["StringFlav:probStoUD"] = rhoEff;
1056  p["StringFlav:probSQtoQQ"] = xEff;
1057  p["StringFlav:probQQ1toQQ0"] = yEff;
1058  p["StringFlav:probQQtoQ"] = xiEff;
1059  p["StringZ:aLund"] = aEff;
1060  p["StringZ:aExtraDiquark"] = adiqEff;
1061  p["StringFlav:kappa"] = kappaEff;
1062 
1063  return (parameters.insert( make_pair(h,p)).second );
1064 
1065 }
1066 
1067 //--------------------------------------------------------------------------
1068 
1069 // Calculate the a parameter.
1070 
1071 double RopeFragPars::aEffective(double aOrig, double thisb, double mT2) {
1072 
1073  // Calculate initial normalization constants.
1074  double N = integrateFragFun(aOrig, bIn, mT2);
1075  double NEff = integrateFragFun(aOrig, thisb, mT2);
1076  int s = (N < NEff) ? -1 : 1;
1077  double da = DELTAA;
1078  double aNew = aOrig - s * da;
1079 
1080  // Iterate until we meet preset convergence criterion.
1081  do {
1082  // Calculate normalization with current a.
1083  NEff = integrateFragFun(aNew, thisb, mT2);
1084  if ( ((N < NEff) ? -1 : 1) != s ) {
1085  s = (N < NEff) ? -1 : 1;
1086  // If we have crossed over the solution, decrease
1087  // the step size and turn around.
1088  da /= 10.0;
1089  }
1090  aNew -= s * da;
1091  if (aNew < 0.0) {aNew = 0.1; break;}
1092  if (aNew > 2.0) {aNew = 2.0; break;}
1093  } while (da > ACONV);
1094  return aNew;
1095 
1096 }
1097 
1098 //--------------------------------------------------------------------------
1099 
1100 // The Lund fragmentation function.
1101 
1102 double RopeFragPars::fragf(double z, double a, double b, double mT2) {
1103 
1104  if (z < ZCUT) return 0.0;
1105  return pow(1 - z, a) * exp(-b * mT2 / z) / z;
1106 
1107 }
1108 
1109 //--------------------------------------------------------------------------
1110 
1111 // Integral of the Lund fragmentation function with given parameter values.
1112 
1113 double RopeFragPars::integrateFragFun(double a, double b, double mT2) {
1114 
1115  // Using Simpson's rule to integrate the Lund fragmentation function.
1116  double nextIter, nextComb;
1117  double thisComb = 0.0, thisIter = 0.0;
1118  // The target error on the integral should never be changed.
1119  double error = 1.0e-2;
1120 
1121  // 20 is the max number of iterations, 3 is min. Should not be changed.
1122  for (int i = 1; i < 20; ++i) {
1123  nextIter = trapIntegrate( a, b, mT2, thisIter, i);
1124  nextComb = (4.0 * nextIter - thisIter) / 3.0;
1125  if (i > 3 && abs(nextComb - thisComb) < error * abs(nextComb))
1126  return nextComb;
1127  thisIter = nextIter;
1128  thisComb = nextComb;
1129  }
1130  infoPtr->errorMsg("RopeFragPars::integrateFragFun:"
1131  "No convergence of frag fun integral.");
1132  return 0.0;
1133 
1134 }
1135 
1136 //--------------------------------------------------------------------------
1137 
1138 // Helper function for integration.
1139 
1140 double RopeFragPars::trapIntegrate( double a, double b, double mT2,
1141  double sOld, int n) {
1142 
1143  // Compute the nth correction to the integral of fragfunc between 0 and 1
1144  // using extended trapezoidal rule.
1145  if (n == 1) return 0.5 * (fragf(0.0, a, b, mT2) + fragf(1.0, a, b, mT2));
1146  // We want 2^(n-2) interior points (intp). Use bitwise shift to speed up.
1147  int intp = 1;
1148  intp <<= n - 2;
1149  double deltaz = 1.0 / double(intp);
1150  double z = 0.5 * deltaz;
1151  double sum = 0.0;
1152  // Do the integral.
1153  for (int i = 0; i < intp; ++i, z += deltaz) sum += fragf( z, a, b, mT2);
1154  return 0.5 * (sOld + sum / double(intp));
1155 
1156 }
1157 
1158 //==========================================================================
1159 
1160 // The FlavourRope class takes care of placing a string breakup in
1161 // the event, and assigning the string breakup effective parameters.
1162 
1163 //--------------------------------------------------------------------------
1164 
1165 // Change the fragmentation parameters.
1166 
1167 bool FlavourRope::doChangeFragPar(StringFlav* flavPtr, StringZ* zPtr,
1168  StringPT * pTPtr, double m2Had, vector<int> iParton, int endId) {
1169 
1170  // The new parameters.
1171  map<string, double> newPar;
1172  if (doBuffon)
1173  newPar = fetchParametersBuffon(m2Had, iParton, endId);
1174  else
1175  newPar = fetchParameters(m2Had, iParton, endId);
1176  // Change settings to new settings.
1177  for (map<string, double>::iterator itr = newPar.begin(); itr!=newPar.end();
1178  ++itr) settingsPtr->parm( itr->first, itr->second);
1179  // Re-initialize flavour, z, and pT selection with new settings.
1180  flavPtr->init();
1181  zPtr->init();
1182  pTPtr->init();
1183  return true;
1184 
1185 }
1186 
1187 //--------------------------------------------------------------------------
1188 
1189 // Find breakup placement and fetch effective parameters using Buffon.
1190 
1191 map<string, double> FlavourRope::fetchParametersBuffon(double m2Had,
1192  vector<int> iParton, int endId) {
1193  // If effective string tension is set manually, use that.
1194  if (fixedKappa) return fp.getEffectiveParameters(h);
1195  if (!ePtr) {
1196  infoPtr->errorMsg("Error in FlavourRope::fetchParametersBuffon:"
1197  " Event pointer not set in FlavourRope");
1198  return fp.getEffectiveParameters(1.0);
1199  }
1200  if(find(hadronized.begin(),hadronized.end(),*iParton.begin()) ==
1201  hadronized.end()){
1202  hadronized.reserve(hadronized.size() + iParton.size());
1203  hadronized.insert(hadronized.end(),iParton.begin(),iParton.end());
1204  }
1205  // Quark string ends, default mode
1206  if (endId != 21){
1207  // Test consistency
1208  if(ePtr->at(*(iParton.begin())).id() != endId &&
1209  ePtr->at(*(iParton.end() - 1)).id() != endId) {
1210  infoPtr->errorMsg("Error in FlavourRope::fetchParametersBuffon:"
1211  " Quark end inconsistency.");
1212  return fp.getEffectiveParameters(1.0);
1213  }
1214 
1215  // First we must let the string vector point in the right direction
1216  if(ePtr->at(*(iParton.begin())).id() != endId)
1217  reverse(iParton.begin(),iParton.end());
1218 
1219  // Initialize a bit
1220  Vec4 hadronic4Momentum(0,0,0,0);
1221  double enh = 1.0;
1222  double dipFrac = -1.0;
1223  vector<int>::iterator dipItr;
1224  // Find out when invariant mass exceeds m2Had
1225  for(dipItr = iParton.begin(); dipItr != iParton.end(); ++dipItr){
1226  double m2Big = hadronic4Momentum.m2Calc();
1227  if( m2Had <= m2Big){
1228  // Approximate the fraction we are in on the dipole, this goes
1229  // in three cases.
1230  // We are at the beginning.
1231  if(m2Had == 0){
1232  dipFrac = 0;
1233  }
1234  // We are somewhere in the first dipole
1235  else if(dipItr - 1 == iParton.begin()){
1236  dipFrac = sqrt(m2Had/m2Big);
1237  }
1238  else{
1239  if(ePtr->at(*(dipItr - 1)).id() != 21) {
1240  infoPtr->errorMsg("Error in FlavourRope::fetchParameters"
1241  "Buffon: Connecting partons should always be gluons.");
1242  return fp.getEffectiveParameters(1.0);
1243  }
1244 
1245  hadronic4Momentum -= 0.5*ePtr->at(*(dipItr -1)).p();
1246  double m2Small = hadronic4Momentum.m2Calc();
1247 
1248  dipFrac = (sqrt(m2Had) - sqrt(m2Small)) /
1249  (sqrt(m2Big) - sqrt(m2Small));
1250  }
1251  break;
1252  }
1253  hadronic4Momentum += ePtr->at(*dipItr).id() == 21 ?
1254  0.5*ePtr->at(*dipItr).p() : ePtr->at(*dipItr).p();
1255  }
1256  // If we reached the end
1257  // we are in a small string that should just be collapsed
1258  if(dipItr == iParton.end())
1259  return fp.getEffectiveParameters(1.0);
1260  // Sanity check
1261  if(dipFrac < 0 || dipFrac > 1) {
1262  infoPtr->errorMsg("Error in FlavourRope::fetchParametersBuffon:"
1263  " Dipole exceed with fraction less than 0 or greater than 1.");
1264  return fp.getEffectiveParameters(1.0);
1265  }
1266  // We now figure out at what rapidity value,
1267  // in the lab system, the string is breaking
1268  double yBreak;
1269  // Trivial case, just inherit
1270  if(dipFrac == 0)
1271  yBreak = ePtr->at(*dipItr).y();
1272  else{
1273  // Sanity check
1274  if(dipItr == iParton.begin()) {
1275  infoPtr->errorMsg("Error in FlavourRope::fetchParametersBuffon:"
1276  " We are somehow before the first dipole on a string.");
1277  return fp.getEffectiveParameters(1.0);
1278  }
1279  double dy = ePtr->at(*dipItr).y() - ePtr->at(*(dipItr - 1)).y();
1280  yBreak = ePtr->at(*(dipItr - 1)).y() + dipFrac*dy;
1281  }
1282  // Count the number of partons in the whole
1283  // event within deltay of breaking point
1284  double p = 1;
1285  double q = 0;
1286 
1287  for(int i = 0; i < ePtr->size(); ++i){
1288  // Don't double count partons from this
1289  // string (no self-overlap)
1290  if(find(iParton.begin(),iParton.end(), i) != iParton.end())
1291  continue;
1292  // Don't count strings that are already hadronized
1293  if(find(hadronized.begin(),hadronized.end(),i) != hadronized.end())
1294  continue;
1295  double pRap = ePtr->at(i).y();
1296  if(pRap > yBreak - rapiditySpan && pRap < yBreak + rapiditySpan ){
1297  // Do a "Buffon" selection to decide whether
1298  // two strings overlap in impact parameter space
1299  // given ratio of string diameter to collision diameter
1300  double r1 = rndmPtr->flat();
1301  double r2 = rndmPtr->flat();
1302  double theta1 = 2*M_PI*rndmPtr->flat();
1303  double theta2 = 2*M_PI*rndmPtr->flat();
1304  // Overlap?
1305  if(4*pow2(stringProtonRatio) > pow2(sqrt(r1)*cos(theta1) -
1306  sqrt(r2)*cos(theta2)) + pow2(sqrt(r1)*sin(theta1) -
1307  sqrt(r2)*sin(theta2))) {
1308  if(rndmPtr->flat() < 0.5) p += 0.5;
1309  else q += 0.5;
1310  }
1311  }
1312  }
1313  enh = 0.25*(2.0*p+q+2.0);
1314 
1315  return fp.getEffectiveParameters(enh);
1316  }
1317  // For closed gluon loops we cannot distinguish the ends.
1318  // Do nothing instead
1319  else{
1320  return fp.getEffectiveParameters(1.0);
1321  }
1322  return fp.getEffectiveParameters(1.0);
1323 }
1324 
1325 //--------------------------------------------------------------------------
1326 // Find breakup placement and fetch effective parameters using Ropewalk.
1327 
1328 map<string, double> FlavourRope::fetchParameters(double m2Had,
1329  vector<int> iParton, int endId) {
1330  // If effective string tension is set manually, use that.
1331  if (fixedKappa) return fp.getEffectiveParameters(h);
1332  if (!ePtr) {
1333  infoPtr->errorMsg("Error in FlavourRope::fetchParameters:"
1334  " Event pointer not set in FlavourRope");
1335  return fp.getEffectiveParameters(1.0);
1336  }
1337  Vec4 mom;
1338  int eventIndex = -1;
1339  // Set direction
1340  bool dirPos;
1341  if( ePtr->at(iParton[0]).id() == endId) dirPos = true;
1342  else if( ePtr->at(iParton[iParton.size() - 1]).id() == endId) dirPos = false;
1343  else {
1344  infoPtr->errorMsg("Error in FlavourRope::fetchParameters:"
1345  " Could not get string direction");
1346  return fp.getEffectiveParameters(1.0);
1347  }
1348 
1349  for (int i = 0, N = iParton.size(); i < N; ++i) {
1350  // Change to right direction
1351  int j = (dirPos ? i : N - 1 - i);
1352  // Skip the junction entry.
1353  if ( iParton[j] < 0) continue;
1354  mom += ePtr->at(iParton[j]).p();
1355  if ( mom.m2Calc() > m2Had) {
1356  eventIndex = j;
1357  break;
1358  }
1359  }
1360 
1361  // We figure out where we are on the dipole.
1362  // Set some values.
1363  double m2Here = mom.m2Calc();
1364  // Catch event indices less than first dipole.
1365  eventIndex = max(eventIndex, 1);
1366  // Dipfrac is the fraction moved in on a dipole.
1367  double dipFrac = 0.;
1368  // Trace one step back if possible.
1369  if (eventIndex > 1) {
1370  mom-=ePtr->at(iParton[eventIndex]).p();
1371  // Starting point can never be less than 0.
1372  double m2Small = max(mom.m2Calc(), 0.);
1373  dipFrac = (sqrt(m2Had) - sqrt(m2Small)) / (sqrt(m2Here) - sqrt(m2Small));
1374  }
1375  else dipFrac = sqrt(m2Had / m2Here);
1376  double enh = rwPtr->getKappaHere( iParton[eventIndex - 1],
1377  iParton[eventIndex], dipFrac);
1378  return fp.getEffectiveParameters(enh);
1379 
1380 }
1381 
1382 //--------------------------------------------------------------------------
1383 // Inteface to he Ropewalk class.
1384 
1385 bool FlavourRope::initEvent(Event& event, ColConfig& colConfig) {
1386 
1387  setEventPtr(event);
1388  if (flag("PartonVertex:setVertex") && !flag("Ropewalk:doBuffon")) {
1389  rwPtr->extractDipoles(event, colConfig);
1390  rwPtr->calculateOverlaps();
1391  }
1392  return true;
1393 
1394 }
1395 
1396 //==========================================================================
1397 
1398 // The RopewalkShover class is an interface to shoving functions in
1399 // the Ropewalk class.
1400 
1401 bool RopewalkShover::stringRepulsion(Event & event, ColConfig & colConfig) {
1402  rwPtr->extractDipoles(event, colConfig);
1403  rwPtr->shoveTheDipoles(event);
1404  return true;
1405 }
1406 
1407 } // End namespace Pythia8
void ae(int tracks=-1, int hits=-1)
This function is to search for the next non-empty event and draw it by looping over StBFChain (readin...
Definition: Ed.C:102
Definition: rb.hh:21
Definition: AgUStep.h:26