10 #include "Pythia8/Ropewalk.h"
23 OverlappingRopeDipole::OverlappingRopeDipole(RopeDipole* d,
double m0,
24 RotBstMatrix& r) : dipole(d), dir(1) {
27 b1 = d->d1Ptr()->getParticlePtr()->vProd() * MM2FM;
29 b2 = d->d2Ptr()->getParticlePtr()->vProd() * MM2FM;
31 y1 = d->d1Ptr()->rap(m0,r);
32 y2 = d->d2Ptr()->rap(m0,r);
33 if (y1 < y2) dir = -1;
41 bool OverlappingRopeDipole::overlap(
double y, Vec4 ba,
double r0) {
43 if (y < min(y1, y2) || y > max(y1, y2))
return false;
44 Vec4 bb = b1 + (b2 - b1) * (y - y1) / (y2 - y1);
46 return (tmp.pT() <= 2 * r0);
54 bool OverlappingRopeDipole::hadronized() {
56 return dipole->hadronized();
70 RopeDipole::RopeDipole(RopeDipoleEnd d1In, RopeDipoleEnd d2In,
int iSubIn,
72 : d1(d1In), d2(d2In), iSub(iSubIn), hasRotFrom(false), hasRotTo(false),
73 isHadronized(false), infoPtr(infoPtrIn) {
76 if (d1In.getParticlePtr()->col() == d2In.getParticlePtr()->acol()
77 && d1In.getParticlePtr()->col() != 0) { }
78 else { d2 = d1In, d1 = d2In; }
86 void RopeDipole::addExcitation(
double ylab, Particle* ex) {
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;
92 if (ex == itr->second)
return;
93 excitations.insert( make_pair(ylab,ex) );
101 void RopeDipole::propagateInit(
double deltat) {
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");
113 double mTc = sqrt(mT2c);
114 double mTa = sqrt(mT2a);
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);
120 d1.getParticlePtr()->vProdAdd(newv1 * FM2MM);
121 d2.getParticlePtr()->vProdAdd(newv2 * FM2MM);
129 void RopeDipole::propagate(
double deltat,
double m0) {
132 propagateInit(deltat);
133 for (map<double, Particle*>::iterator eItr = excitations.begin();
134 eItr != excitations.end(); ++eItr) {
136 Vec4 em = eItr->second->p();
137 em.rotbst(getDipoleLabFrame());
141 Vec4 newVert = Vec4( deltat * em.px() / em.pT(),
142 deltat * em.py() / em.pT(), 0, 0);
143 eItr->second->vProdAdd(newVert * FM2MM);
145 else eItr->second->vProd( bInterpolateLab(eItr->first,m0) * FM2MM);
154 void RopeDipole::excitationsToString(
double m0,
Event& event) {
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;
162 excitations.erase(eraseMe);
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.");
176 vector<int> daughters;
180 if (d1.rap(m0) == minRapidity(m0)) {
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))));
191 d2.getParticlePtr()->acol(acol);
192 event[d2.getNe()].acol(acol);
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))));
205 d2.getParticlePtr()->acol(acol);
206 event[d2.getNe()].acol(acol);
208 bool stringEnd =
false;
209 if (d2.getParticlePtr()->col() == 0) stringEnd =
true;
212 event[d1.getNe()].statusNeg();
213 Particle cc1 = *d1.getParticlePtr();
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 ] );
219 event[d2.getNe()].statusNeg();
220 Particle cc2 = *d2.getParticlePtr();
222 cc2.mothers(d2.getNe(),d2.getNe());
223 int did =
event.append(cc2);
224 event[d2.getNe()].daughters(did,did);
233 void RopeDipole::splitMomentum(Vec4 mom, Particle* p1, Particle* p2,
236 Vec4 p1new = p1->p() + frac * mom;
237 Vec4 p2new = p2->p() + (1. - frac) * mom;
250 bool RopeDipole::recoil(Vec4& pg,
bool dummy) {
254 if (d1.rap(1.0) > d2.rap(1.0)) sign = -1;
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();
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;
275 double sqarg = pow2(pplus * pminus - mta2 - mtc2) - 4. * mta2 * mtc2;
276 if (sqarg <= 0.0 )
return false;
278 ppa = 0.5 * (pplus * pminus + mta2 - mtc2 + sqrt(sqarg)) / pminus;
283 if ( ppa * mtc < ppc * mta )
return false;
286 pma = 0.5 * (pplus * pminus + mta2 - mtc2 + sqrt(sqarg)) / pplus;
291 if ( ppa*mtc > ppc*mta )
return false;
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));
311 RotBstMatrix RopeDipole::getDipoleRestFrame() {
313 if (hasRotTo)
return rotTo;
316 r.toCMframe(d1.getParticlePtr()->p(),d2.getParticlePtr()->p());
327 RotBstMatrix RopeDipole::getDipoleLabFrame() {
329 if(hasRotFrom)
return rotFrom;
332 r.fromCMframe(d1.getParticlePtr()->p(),d2.getParticlePtr()->p());
342 Vec4 RopeDipole::dipoleMomentum() {
344 Vec4 ret = d1.getParticlePtr()->p() + d2.getParticlePtr()->p();
355 Vec4 RopeDipole::bInterpolateDip(
double y,
double m0) {
356 if(!hasRotTo) getDipoleRestFrame();
357 Vec4 bb1 = d1.getParticlePtr()->vProd() * MM2FM;
359 Vec4 bb2 = d2.getParticlePtr()->vProd() * MM2FM;
361 double y1 = d1.rap(m0,rotTo);
362 double y2 = d2.rap(m0,rotTo);
363 return bb1 + y * (bb2 - bb1) / (y2 - y1);
371 Vec4 RopeDipole::bInterpolateLab(
double y,
double m0) {
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);
386 Vec4 RopeDipole::bInterpolate(
double y, RotBstMatrix
rb,
double m0) {
388 Vec4 bb1 = d1.getParticlePtr()->vProd() * MM2FM;
389 Vec4 bb2 = d2.getParticlePtr()->vProd() * MM2FM;
392 double y1 = d1.rap(m0);
393 double y2 = d2.rap(m0);
394 return bb1 + y * (bb2 - bb1) / (y2 - y1);
403 pair<int, int> RopeDipole::getOverlaps(
double yfrac,
double m0,
double r0) {
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;
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;
417 return make_pair(m,n);
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) { }
438 void setParticlePtrs(Particle* p1, Particle* p2) {
442 dip1->addExcitation(y, pp1);
443 dip2->addExcitation(y, pp2);
447 void shove(
double dpx,
double dpy) {
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);
461 Vec4 deltap1 = p1new - p1;
462 Vec4 deltap2 = p2new - p2;
464 if ( dip2->recoil(deltap2) ) {
465 if ( dip1->recoil(deltap1) ) {
477 return dip1->bInterpolateDip(y,m0) -
478 dip2->bInterpolateDip(y,m0);
503 bool Ropewalk::init() {
505 StringInteractions::init();
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");
527 if ( flag(
"Ropewalk:doShoving") ) {
529 if (deltat > tShove) {
530 infoPtr->errorMsg(
"Error in Ropewalk::init: "
531 "deltat cannot be larger than tShove");
534 if ( !flag(
"PartonVertex:setVertex") ) {
535 infoPtr->errorMsg(
"Error in Ropewalk::init: "
536 "Shoving enabled, but no vertex information.");
539 stringrepPtr = make_shared<RopewalkShover>(*this);
540 registerSubObject(*stringrepPtr);
541 if ( !stringrepPtr->init() )
return false;
544 if ( flag(
"Ropewalk:doFlavour") ) {
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");
555 fragmodPtr = make_shared<FlavourRope>(*this);
556 registerSubObject(*fragmodPtr);
557 if ( !fragmodPtr->init() )
return false;
569 double Ropewalk::averageKappa() {
573 for (DMap::iterator itr = dipoles.begin(); itr != dipoles.end(); ++itr) {
576 pair<int,int> overlap = itr->second.getOverlaps( rndmPtr->flat(), m0, r0);
580 pair<double, double> pq = select( overlap.first + 1, overlap.second,
582 double enh = 0.25 * (2. + 2. * pq.first + pq.second);
583 kap += (enh > 1.0 ? enh : 1.0);
595 double Ropewalk::getKappaHere(
int e1,
int e2,
double yfrac) {
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);
605 pair<int, int> overlap = d->getOverlaps(yfrac, m0, r0);
606 pair<double, double> pq;
610 pq = make_pair(overlap.first + 1, overlap.second);
614 pq = select(overlap.first + 1, overlap.second, rndmPtr);
617 double enh = 0.25 * (2. + 2. * pq.first + pq.second);
618 return (enh > 1.0 ? enh : 1.0);
626 bool Ropewalk::calculateOverlaps() {
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;
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;
641 for (multimap< pair<int,int>, RopeDipole>::iterator itr2 = dipoles.begin();
642 itr2 != dipoles.end(); ++itr2) {
643 RopeDipole* d2 = &(itr2->second);
646 if (d1 == d2)
continue;
647 if (d2->dipoleMomentum().m2Calc() < pow2(m0))
continue;
650 OverlappingRopeDipole od(d2, m0, dipoleRestFrame);
651 if (min(od.y1, od.y2) > yc1 || max(od.y1, od.y2) < ya1 || od.y1 == od.y2)
654 d1->addOverlappingDipole(od);
665 pair<int, int> Ropewalk::select(
int m,
int n, Rndm* rndm) {
672 while (mm + nn > 0) {
675 if (rndm->flat() < 0.5 && mm > 0) {
679 double p1 = multiplicity(p + 1, q);
680 double p2 = multiplicity(p, q - 1);
681 double p3 = multiplicity(p - 1, q + 1);
684 double sum = p1 + p2 + p3;
685 p1 /= sum, p2 /= sum, p3 /= sum;
688 double pick = rndm->flat();
690 else if (pick < p1 + p2) --q;
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();
704 else if (pick < p1 + p2) --p;
710 return make_pair( (p < 0 ? 0 : p), (q < 0 ? 0 : q) );
718 void Ropewalk::shoveTheDipoles(
Event& event) {
722 for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr)
723 dItr->second.propagateInit(tInit);
727 multimap<double, RopeDipole *> rapDipoles;
732 for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr) {
733 RopeDipole* dip = &(dItr->second);
735 rapDipoles.insert( make_pair(dip->maxRapidity(m0), dip) );
737 if (dip->minRapidity(m0) < ymin) ymin = dip->minRapidity(m0);
738 if (dip->maxRapidity(m0) > ymax) ymax = dip->maxRapidity(m0);
742 vector<double> rapidities;
743 for (
double y = ymin; y < ymax; y += deltay) rapidities.push_back(y);
746 map<double, vector<Exc> > exPairs;
747 for (
int i = 0, N = eParticles.size(); i < N; ++i) eParticles[i].clear();
749 for (
int i = 0, N = rapidities.size(); i < N; ++i) {
751 eParticles.push_back( vector<Particle>() );
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);
763 vector<int> eraseDipoles;
765 for (
int j = 0, M = tmp.size(); j < M; ++j) {
768 if (!tmp[j]->recoil(ex,
true) ) {
769 eraseDipoles.push_back(j);
773 for (
int j = 0, M = eraseDipoles.size(); j < M; ++j) {
774 tmp.erase( tmp.begin() + (eraseDipoles[j]-j) );
777 if(
int(tmp.size()) > 1)
778 for (
int j = 0, M = tmp.size(); j < M; ++j) {
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);
787 exPairs[ySample] = vector<Exc>();
788 for (
int j = 0, M = tmp.size(); j < M; ++j)
789 for (
int k = 0; k < M; ++k) {
791 if(j != k && tmp[j]->index() != tmp[k]->index() )
792 exPairs[ySample].push_back( Exc(ySample, m0, i, j, k, tmp[j],
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] );
807 for (
double t = tInit; t < tShove + tInit; t += deltat) {
809 for (map<
double, vector<Exc> >::iterator slItr = exPairs.begin();
810 slItr != exPairs.end(); ++slItr)
812 for (
int i = 0, N = slItr->second.size(); i < N; ++i) {
813 Exc& ep = slItr->second[i];
815 Vec4 direction = ep.direction();
820 double rt = max(t, 1. / showerCut / 5.068);
821 rt = min(rt, r0 * gExponent);
822 double dist = direction.pT();
824 if (dist < rCutOff * rt) {
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;
835 for (DMap::iterator dItr = dipoles.begin(); dItr != dipoles.end(); ++dItr)
836 dItr->second.propagate(deltat, m0);
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);
850 bool Ropewalk::extractDipoles(
Event& event, ColConfig& colConfig) {
854 for (
int iSub = 0; iSub < colConfig.size(); ++iSub) {
857 if (colConfig[iSub].hasJunction && !shoveJunctionStrings)
continue;
858 if (colConfig[iSub].isClosed && !shoveGluonLoops)
continue;
859 if (colConfig[iSub].massExcess <= mStringMin && !shoveMiniStrings)
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) {
870 RopeDipoleEnd next( &event, stringPartons[iPar]);
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)) );
881 dipoles.insert( pair< pair<int, int>, RopeDipole>(dipoleER,
882 RopeDipole( previous, next, iSub, infoPtr)));
905 const double RopeFragPars::DELTAA = 0.1;
908 const double RopeFragPars::ACONV = 0.001;
911 const double RopeFragPars::ZCUT = 1.0e-4;
917 bool RopeFragPars::init() {
920 beta = parm(
"Ropewalk:beta");
924 string params [len] = {
"StringPT:sigma",
"StringZ:aLund",
925 "StringZ:aExtraDiquark",
"StringZ:bLund",
"StringFlav:probStoUD",
926 "StringFlav:probSQtoQQ",
"StringFlav:probQQ1toQQ0",
"StringFlav:probQQtoQ",
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]);
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.");
949 map<string,double> RopeFragPars::getEffectiveParameters(
double h) {
951 map<double, map<string, double> >::iterator parItr = parameters.find(h);
954 if ( parItr != parameters.end())
return parItr->second;
957 if (!calculateEffectiveParameters(h))
958 infoPtr->errorMsg(
"Error in RopeFragPars::getEffectiveParameters:"
959 " calculating effective parameters.");
962 if (!insertEffectiveParameters(h))
963 infoPtr->errorMsg(
"Error in RopeFragPars::getEffectiveParameters:"
964 " inserting effective parameters.");
967 return getEffectiveParameters(h);
975 double RopeFragPars::getEffectiveA(
double thisb,
double mT2,
bool isDiquark) {
978 if (thisb == bIn)
return (isDiquark ? aIn + adiqIn : aIn);
981 map<double, double>* aMPtr = (isDiquark ? &aDiqMap : &aMap);
982 double bmT2 = mT2 * thisb;
985 map<double,double>::iterator aItr = aMPtr->find(bmT2);
986 if (aItr != aMPtr->end())
return aItr->second;
989 double ae = ( isDiquark ? aEffective(aIn + adiqIn, thisb, mT2)
990 : aEffective(aIn, thisb, mT2) );
992 double suba = getEffectiveA(thisb, mT2,
false);
993 aMPtr->insert( make_pair(bmT2,
ae - suba) );
995 else aMPtr->insert(make_pair(bmT2,
ae));
1004 bool RopeFragPars::calculateEffectiveParameters(
double h) {
1006 if (h <= 0)
return false;
1007 double hinv = 1.0 / h;
1011 kappaEff = kappaIn * h;
1013 rhoEff = pow(rhoIn, hinv);
1015 xEff = pow(xIn, hinv);
1017 yEff = pow(yIn, hinv);
1019 sigmaEff = sigmaIn * sqrt(h);
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)
1029 xiEff = alphaEff * beta * pow( xiIn / alpha / beta, hinv);
1030 if (xiEff > 1.0) xiEff = 1.0;
1031 if (xiEff < xiIn) xiEff = xiIn;
1034 bEff = (2. + rhoEff) / (2. + rhoIn) * bIn;
1035 if (bEff < bIn) bEff = bIn;
1036 if (bEff > 2.0) bEff = 2.0;
1039 aEff = getEffectiveA( bEff, 1.0,
false);
1040 adiqEff = getEffectiveA( bEff, 1.0,
true) - aEff;
1050 bool RopeFragPars::insertEffectiveParameters(
double h) {
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;
1063 return (parameters.insert( make_pair(h,p)).second );
1071 double RopeFragPars::aEffective(
double aOrig,
double thisb,
double mT2) {
1074 double N = integrateFragFun(aOrig, bIn, mT2);
1075 double NEff = integrateFragFun(aOrig, thisb, mT2);
1076 int s = (N < NEff) ? -1 : 1;
1078 double aNew = aOrig - s * da;
1083 NEff = integrateFragFun(aNew, thisb, mT2);
1084 if ( ((N < NEff) ? -1 : 1) != s ) {
1085 s = (N < NEff) ? -1 : 1;
1091 if (aNew < 0.0) {aNew = 0.1;
break;}
1092 if (aNew > 2.0) {aNew = 2.0;
break;}
1093 }
while (da > ACONV);
1102 double RopeFragPars::fragf(
double z,
double a,
double b,
double mT2) {
1104 if (z < ZCUT)
return 0.0;
1105 return pow(1 - z, a) * exp(-b * mT2 / z) / z;
1113 double RopeFragPars::integrateFragFun(
double a,
double b,
double mT2) {
1116 double nextIter, nextComb;
1117 double thisComb = 0.0, thisIter = 0.0;
1119 double error = 1.0e-2;
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))
1127 thisIter = nextIter;
1128 thisComb = nextComb;
1130 infoPtr->errorMsg(
"RopeFragPars::integrateFragFun:"
1131 "No convergence of frag fun integral.");
1140 double RopeFragPars::trapIntegrate(
double a,
double b,
double mT2,
1141 double sOld,
int n) {
1145 if (n == 1)
return 0.5 * (fragf(0.0, a, b, mT2) + fragf(1.0, a, b, mT2));
1149 double deltaz = 1.0 / double(intp);
1150 double z = 0.5 * deltaz;
1153 for (
int i = 0; i < intp; ++i, z += deltaz) sum += fragf( z, a, b, mT2);
1154 return 0.5 * (sOld + sum / double(intp));
1167 bool FlavourRope::doChangeFragPar(StringFlav* flavPtr, StringZ* zPtr,
1168 StringPT * pTPtr,
double m2Had, vector<int> iParton,
int endId) {
1171 map<string, double> newPar;
1173 newPar = fetchParametersBuffon(m2Had, iParton, endId);
1175 newPar = fetchParameters(m2Had, iParton, endId);
1177 for (map<string, double>::iterator itr = newPar.begin(); itr!=newPar.end();
1178 ++itr) settingsPtr->parm( itr->first, itr->second);
1191 map<string, double> FlavourRope::fetchParametersBuffon(
double m2Had,
1192 vector<int> iParton,
int endId) {
1194 if (fixedKappa)
return fp.getEffectiveParameters(h);
1196 infoPtr->errorMsg(
"Error in FlavourRope::fetchParametersBuffon:"
1197 " Event pointer not set in FlavourRope");
1198 return fp.getEffectiveParameters(1.0);
1200 if(find(hadronized.begin(),hadronized.end(),*iParton.begin()) ==
1202 hadronized.reserve(hadronized.size() + iParton.size());
1203 hadronized.insert(hadronized.end(),iParton.begin(),iParton.end());
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);
1216 if(ePtr->at(*(iParton.begin())).id() != endId)
1217 reverse(iParton.begin(),iParton.end());
1220 Vec4 hadronic4Momentum(0,0,0,0);
1222 double dipFrac = -1.0;
1223 vector<int>::iterator dipItr;
1225 for(dipItr = iParton.begin(); dipItr != iParton.end(); ++dipItr){
1226 double m2Big = hadronic4Momentum.m2Calc();
1227 if( m2Had <= m2Big){
1235 else if(dipItr - 1 == iParton.begin()){
1236 dipFrac = sqrt(m2Had/m2Big);
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);
1245 hadronic4Momentum -= 0.5*ePtr->at(*(dipItr -1)).p();
1246 double m2Small = hadronic4Momentum.m2Calc();
1248 dipFrac = (sqrt(m2Had) - sqrt(m2Small)) /
1249 (sqrt(m2Big) - sqrt(m2Small));
1253 hadronic4Momentum += ePtr->at(*dipItr).id() == 21 ?
1254 0.5*ePtr->at(*dipItr).p() : ePtr->at(*dipItr).p();
1258 if(dipItr == iParton.end())
1259 return fp.getEffectiveParameters(1.0);
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);
1271 yBreak = ePtr->at(*dipItr).y();
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);
1279 double dy = ePtr->at(*dipItr).y() - ePtr->at(*(dipItr - 1)).y();
1280 yBreak = ePtr->at(*(dipItr - 1)).y() + dipFrac*dy;
1287 for(
int i = 0; i < ePtr->size(); ++i){
1290 if(find(iParton.begin(),iParton.end(), i) != iParton.end())
1293 if(find(hadronized.begin(),hadronized.end(),i) != hadronized.end())
1295 double pRap = ePtr->at(i).y();
1296 if(pRap > yBreak - rapiditySpan && pRap < yBreak + rapiditySpan ){
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();
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;
1313 enh = 0.25*(2.0*p+q+2.0);
1315 return fp.getEffectiveParameters(enh);
1320 return fp.getEffectiveParameters(1.0);
1322 return fp.getEffectiveParameters(1.0);
1328 map<string, double> FlavourRope::fetchParameters(
double m2Had,
1329 vector<int> iParton,
int endId) {
1331 if (fixedKappa)
return fp.getEffectiveParameters(h);
1333 infoPtr->errorMsg(
"Error in FlavourRope::fetchParameters:"
1334 " Event pointer not set in FlavourRope");
1335 return fp.getEffectiveParameters(1.0);
1338 int eventIndex = -1;
1341 if( ePtr->at(iParton[0]).id() == endId) dirPos =
true;
1342 else if( ePtr->at(iParton[iParton.size() - 1]).
id() == endId) dirPos =
false;
1344 infoPtr->errorMsg(
"Error in FlavourRope::fetchParameters:"
1345 " Could not get string direction");
1346 return fp.getEffectiveParameters(1.0);
1349 for (
int i = 0, N = iParton.size(); i < N; ++i) {
1351 int j = (dirPos ? i : N - 1 - i);
1353 if ( iParton[j] < 0)
continue;
1354 mom += ePtr->at(iParton[j]).p();
1355 if ( mom.m2Calc() > m2Had) {
1363 double m2Here = mom.m2Calc();
1365 eventIndex = max(eventIndex, 1);
1367 double dipFrac = 0.;
1369 if (eventIndex > 1) {
1370 mom-=ePtr->at(iParton[eventIndex]).p();
1372 double m2Small = max(mom.m2Calc(), 0.);
1373 dipFrac = (sqrt(m2Had) - sqrt(m2Small)) / (sqrt(m2Here) - sqrt(m2Small));
1375 else dipFrac = sqrt(m2Had / m2Here);
1376 double enh = rwPtr->getKappaHere( iParton[eventIndex - 1],
1377 iParton[eventIndex], dipFrac);
1378 return fp.getEffectiveParameters(enh);
1385 bool FlavourRope::initEvent(
Event& event, ColConfig& colConfig) {
1388 if (flag(
"PartonVertex:setVertex") && !flag(
"Ropewalk:doBuffon")) {
1389 rwPtr->extractDipoles(event, colConfig);
1390 rwPtr->calculateOverlaps();
1401 bool RopewalkShover::stringRepulsion(
Event & event, ColConfig & colConfig) {
1402 rwPtr->extractDipoles(event, colConfig);
1403 rwPtr->shoveTheDipoles(event);
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...