6 #include "HadronScatter.h"
77 const int SigmaPartialWave::LSHIFT = 1000000;
78 const int SigmaPartialWave::ISHIFT = 1000;
81 const double SigmaPartialWave::CONVERT2MB = 0.389380;
84 const double SigmaPartialWave::WCMBIN = 0.005;
85 const double SigmaPartialWave::CTBIN = 0.2;
87 const int SigmaPartialWave::SUBBIN = 2;
88 const int SigmaPartialWave::ITER = 2;
90 const double SigmaPartialWave::MASSSAFETY = 0.001;
91 const double SigmaPartialWave::GRIDSAFETY = 0.05;
98 bool SigmaPartialWave::init(
int processIn,
string xmlPath,
string filename,
99 Info *infoPtrIn, ParticleData *particleDataPtrIn,
103 particleDataPtr = particleDataPtrIn;
107 if (processIn < 0 || processIn > 2) {
108 infoPtr->errorMsg(
"Error in SigmaPartialWave::init: "
119 if (!readFile(xmlPath, filename))
return false;
123 if (Lmax > 0) PlVec[0] = 1.;
127 if (Lmax > 0) PlpVec[0] = 0.;
128 if (Lmax > 1) PlpVec[1] = 1.;
142 bool SigmaPartialWave::readFile(
string xmlPath,
string filename) {
144 string fullPath = xmlPath + filename;
145 ifstream ifs(fullPath.c_str());
147 infoPtr->errorMsg(
"Error in SigmaPartialWave::init: "
148 "could not read data file");
164 vector < int > Lvec, Ivec, Jvec;
171 for (
unsigned int i = 0; i < line.length(); i++)
172 line[i] = tolower(line[i]);
173 string::size_type startPos = line.find_first_not_of(
" ");
174 if (startPos != string::npos) line = line.substr(startPos);
176 if (line.length() == 0 || line[0] ==
'#')
continue;
180 vector < string > token;
182 startPos = lineT.find_first_of(
" ");
183 token.push_back(lineT.substr(0, startPos));
184 if (startPos == string::npos)
break;
185 startPos = lineT.find_first_not_of(
" ", startPos + 1);
186 if (startPos == string::npos)
break;
187 lineT = lineT.substr(startPos);
191 if (token[0] ==
"set") {
192 bool badSetting =
false;
195 if (token[1] ==
"etype") {
196 if (token[2] ==
"wcm") eType = 0;
197 else if (token[2] ==
"tlab") eType = 1;
198 else badSetting =
true;
201 }
else if (token[1] ==
"eunit") {
202 if (token[2] ==
"gev") eUnit = 0;
203 else if (token[2] ==
"mev") eUnit = 1;
204 else badSetting =
true;
207 }
else if (token[1] ==
"input") {
208 if (token[2] ==
"eta,delta") input = 0;
209 else if (token[2] ==
"sn,delta") input = 1;
210 else if (token[2] ==
"tr,ti") input = 2;
211 else if (token[2] ==
"mod,phi") input = 3;
212 else badSetting =
true;
215 }
else if (token[1] ==
"dunit") {
216 if (token[2] ==
"deg") dUnit = 0;
217 else if (token[2] ==
"rad") dUnit = 1;
218 else badSetting =
true;
221 }
else if (token[1] ==
"norm") {
222 if (token[2] ==
"0") norm = 0;
223 else if (token[2] ==
"1") norm = 1;
224 else badSetting =
true;
229 infoPtr->errorMsg(
"Error in SigmaPartialWave::init: "
237 if (line.substr(0, 1).find_first_of(
"0123456789.") != 0) {
239 Lvec.clear(); Ivec.clear(); Jvec.clear();
242 bool badHeader =
false;
243 for (
unsigned int i = 1; i < token.size(); i++) {
245 startPos = token[i].find_first_of(
",");
246 if (startPos == string::npos) { badHeader =
true;
break; }
247 string Lstr = token[i].substr(0, startPos);
248 token[i] = token[i].substr(startPos + 1);
251 startPos = token[i].find_first_of(
", ");
252 if (startPos == string::npos) {
256 Istr = token[i].substr(0, startPos);
257 token[i] = token[i].substr(startPos + 1);
261 if (token[i].length() != 0) Jstr = token[i];
265 stringstream Lss(Lstr); Lss >> L;
266 stringstream Iss(Istr); Iss >> I;
267 stringstream Jss(Jstr); Jss >> J;
268 if (Lss.fail() || Iss.fail() || Jss.fail())
269 { badHeader =
true;
break; }
277 infoPtr->errorMsg(
"Error in SigmaPartialWave::init: "
278 "malformed header line");
284 bool badData =
false;
287 if (token.size() != 2 * Lvec.size() + 1) badData =
true;
292 stringstream eSS(token[0]);
294 if (eSS.fail()) badData =
true;
296 if (eUnit == 1) eNow *= 1e-3;
298 if (eType == 1) eNow = sqrt(2. * mB * eNow + pow2(mA + mB));
299 binMax = max(binMax, eNow);
304 for (
unsigned int i = 1; i < token.size(); i += 2) {
306 int LIJidx = (i - 1) / 2;
307 int L = Lvec[LIJidx];
308 int I = Ivec[LIJidx];
309 int J = Jvec[LIJidx];
312 stringstream i1SS(token[i]); i1SS >> i1;
313 stringstream i2SS(token[i + 1]); i2SS >> i2;
314 if (i1SS.fail() || i2SS.fail()) { badData =
true;
break; }
317 if (input == 1) i1 = sqrt(1. - i1);
319 if ((input == 0 || input == 1 || input == 3) &&
320 dUnit == 0) i2 *= M_PI / 180.;
324 if (input == 0 || input == 1) {
325 T = (i1 * exp(2. * complex(0., 1.) * i2) - 1.) /
326 2. / complex(0., 1.);
327 }
else if (input == 2) {
329 }
else if (input == 3) {
330 T = i1 * exp(complex(0., 1.) * i2);
334 pwData[L * LSHIFT + I * ISHIFT + J][eNow] = T;
338 infoPtr->errorMsg(
"Error in SigmaPartialWave::init: "
339 "malformed data line");
346 if (!ifs.eof()) { ifs.close();
return false; }
359 void SigmaPartialWave::setupSubprocesses() {
367 sp2in[0] = pair < int, int > ( 211, 211);
368 sp2in[1] = pair < int, int > ( 211, -211);
369 sp2in[2] = pair < int, int > ( 211, 111);
370 sp2in[3] = pair < int, int > ( 111, 111);
371 sp2in[4] = pair < int, int > (-211, 111);
372 sp2in[5] = pair < int, int > (-211, -211);
374 for (
int i = 0; i < subprocessMax; i++) in2sp[sp2in[i]] = i;
376 isoCoeff[0][0] = 0.; isoCoeff[0][2] = 0.; isoCoeff[0][4] = 1.;
377 isoCoeff[1][0] = 1./3.; isoCoeff[1][2] = 1./2.; isoCoeff[1][4] = 1./6.;
378 isoCoeff[2][0] = 0.; isoCoeff[2][2] = 1./2.; isoCoeff[2][4] = 1./2.;
379 isoCoeff[3][0] = 1./3.; isoCoeff[3][2] = 0.; isoCoeff[3][4] = 2./3.;
380 isoCoeff[4][0] = 0.; isoCoeff[4][2] = 1./2.; isoCoeff[4][4] = 1./2.;
381 isoCoeff[5][0] = 0.; isoCoeff[5][2] = 0.; isoCoeff[5][4] = 1.;
388 if (process == 1) { id1 = 321; id2 = 311; }
389 else { id1 = 2212; id2 = 2112; }
393 sp2in[0] = pair < int, int > ( 211, id1);
394 sp2in[1] = pair < int, int > ( 211, id2);
395 sp2in[2] = pair < int, int > ( 111, id1);
396 sp2in[3] = pair < int, int > ( 111, id2);
397 sp2in[4] = pair < int, int > (-211, id1);
398 sp2in[5] = pair < int, int > (-211, id2);
400 isoCoeff[0][1] = 0.; isoCoeff[0][3] = 1.;
401 isoCoeff[1][1] = 2. / 3.; isoCoeff[1][3] = 1. / 3.;
402 isoCoeff[2][1] = 1. / 3.; isoCoeff[2][3] = 2. / 3.;
403 isoCoeff[3][1] = 1. / 3.; isoCoeff[3][3] = 2. / 3.;
404 isoCoeff[4][1] = 2. / 3.; isoCoeff[4][3] = 1. / 3.;
405 isoCoeff[5][1] = 0.; isoCoeff[5][3] = 1.;
407 for (
int i = 0; i < 6; i++) {
408 id1 = ((sp2in[i].first == 111) ? +1 : -1) * sp2in[i].first;
409 sp2in[i + 6] = pair < int, int > (id1, -sp2in[i].second);
410 isoCoeff[i + 6] = isoCoeff[i];
413 for (
int i = 0; i < subprocessMax; i++) in2sp[sp2in[i]] = i;
426 void SigmaPartialWave::setupGrid() {
431 gridMax.resize(subprocessMax);
432 gridNorm.resize(subprocessMax);
433 for (
int sp = 0; sp < subprocessMax; sp++) {
438 int nBin1 = int( (binMax - mA - mB) / WCMBIN );
439 gridMax[subprocess].resize(nBin1);
440 gridNorm[subprocess].resize(nBin1);
441 for (
int n1 = 0; n1 < nBin1; n1++) {
443 double bl1 = mA + mB + double(n1) * WCMBIN;
444 double bu1 = bl1 + WCMBIN;
447 int nBin2 = int( 2. / CTBIN );
448 gridMax[subprocess][n1].resize(nBin2);
449 for (
int n2 = 0; n2 < nBin2; n2++) {
451 double bl2 = -1. + double(n2) * CTBIN;
452 double bu2 = bl2 + CTBIN;
456 double bl3 = bl1, bu3 = bu1, bl4 = bl2, bu4 = bu2;
457 for (
int iter = 0; iter < ITER; iter++) {
458 int i3Save = -1, i4Save = -1;
459 double step3 = (bu3 - bl3) /
double(SUBBIN);
460 double step4 = (bu4 - bl4) /
double(SUBBIN);
461 for (
int i3 = 0; i3 <= SUBBIN; i3++) {
462 double Wcm = bl3 + double(i3) * step3;
463 for (
int i4 = 0; i4 <= SUBBIN; i4++) {
464 double ct = bl4 + double(i4) * step4;
465 double ds = dSigma(Wcm, ct);
474 if (i3Save == -1 && i4Save == -1)
break;
476 bl3 = bl3 + ((i3Save == 0) ? 0. : i3Save - 1.) * step3;
477 bu3 = bl3 + ((i3Save == SUBBIN) ? 1. : 2.) * step3;
480 bl4 = bl4 + ((i4Save == 0) ? 0. : i4Save - 1.) * step4;
481 bu4 = bl4 + ((i4Save == SUBBIN) ? 1. : 2.) * step4;
486 gridMax[subprocess][n1][n2] = maxSig * (1. + GRIDSAFETY);
487 gridNorm[subprocess][n1] += maxSig * (1. + GRIDSAFETY) * CTBIN;
488 sigElMax = max(sigElMax, maxSig);
502 double SigmaPartialWave::pickCosTheta(
double Wcm) {
504 int WcmBin = int((Wcm - mA - mB) / WCMBIN);
505 if (WcmBin < 0) WcmBin = 0;
506 if (WcmBin >=
int(gridMax[subprocess].size()))
507 WcmBin = int(gridMax[subprocess].size() - 1);
514 double y = rndmPtr->flat() * gridNorm[subprocess][WcmBin];
517 for (ctBin = 0; ctBin < int(2. / CTBIN); ctBin++) {
518 if (sum + CTBIN * gridMax[subprocess][WcmBin][ctBin] > y)
break;
519 sum += CTBIN * gridMax[subprocess][WcmBin][ctBin];
523 double x1 = -1. + CTBIN * double(ctBin);
525 double x2 = x1 + CTBIN;
526 double y2 = sum + CTBIN * gridMax[subprocess][WcmBin][ctBin];
527 ct = (x2 - x1) / (y2 - y1) * (y - y1) + x1;
528 wgt = dSigma(Wcm, ct) / gridMax[subprocess][WcmBin][ctBin];
530 infoPtr->errorMsg(
"Warning in SigmaPartialWave::pickCosTheta: "
531 "weight above unity");
534 }
while (wgt <= rndmPtr->flat());
543 bool SigmaPartialWave::setSubprocess(
int spIn) {
544 if (sp2in.find(spIn) == sp2in.end())
return false;
546 pair < int, int > in = sp2in[spIn];
548 mA = particleDataPtr->m0(idA);
550 mB = particleDataPtr->m0(idB);
554 bool SigmaPartialWave::setSubprocess(
int idAin,
int idBin) {
555 pair < int, int > in = pair < int, int > (idAin, idBin);
556 if (in2sp.find(in) == in2sp.end()) {
558 swap(in.first, in.second);
559 if (in2sp.find(in) == in2sp.end())
return false;
561 subprocess = in2sp[in];
563 mA = particleDataPtr->m0(idA);
565 mB = particleDataPtr->m0(idB);
573 double SigmaPartialWave::sigma(
int mode,
double Wcm,
double cTheta) {
575 if (Wcm < (mA + mB + MASSSAFETY))
return 0.;
578 complex amp[2] = { complex(0., 0.) };
582 double s = pow2(Wcm);
583 double k2 = (s - pow2(mB + mA)) * (s - pow2(mB - mA)) / 4. / s;
588 if (process == 2) sTheta = sqrt(1. - pow2(cTheta));
589 legendreP(cTheta, ((process == 2) ?
true :
false));
593 for (
int L = 0; L < Lmax; L++) {
596 complex ampJ[2] = { complex(0., 0.) };
597 int Jstart = (process != 2) ? 0 : 2 * L - 1;
598 int Jend = (process != 2) ? 1 : 2 * L + 2;
599 int Jstep = (process != 2) ? 1 : 2;
601 for (
int J = Jstart; J < Jend; J += Jstep, Jcount++) {
604 for (
int I = 0; I < Imax; I++) {
605 if (isoCoeff[subprocess][I] == 0.)
continue;
608 int LIJ = L * LSHIFT + I * ISHIFT + J;
609 if (pwData.find(LIJ) == pwData.end())
continue;
612 map < double, complex >::iterator it = pwData[LIJ].upper_bound(Wcm);
613 if (it == pwData[LIJ].begin()) ++it;
615 if (it == pwData[LIJ].end()) {
616 ar = (--it)->second.real();
617 ai = it->second.imag();
619 double eA = it->first;
620 complex ampA = (it--)->second;
621 double eB = it->first;
622 complex ampB = it->second;
624 ar = (ampA.real() - ampB.real()) / (eA - eB) *
625 (Wcm - eB) + ampB.real();
626 ai = (ampA.imag() - ampB.imag()) / (eA - eB) *
627 (Wcm - eB) + ampB.imag();
631 ampJ[Jcount] += isoCoeff[subprocess][I] * complex(ar, ai);
637 if (process == 0 || process == 1) {
638 sig += (2. * L + 1.) * (ampJ[0] * conj(ampJ[0])).real();
639 }
else if (process == 2) {
640 sig += ( (L + 0.) * (ampJ[0] * conj(ampJ[0])).real() +
641 (L + 1.) * (ampJ[1] * conj(ampJ[1])).real() );
645 }
else if (mode == 1) {
646 if (process == 0 || process == 1) {
647 sig += (2. * L + 1.) * ampJ[0].imag();
648 }
else if (process == 2) {
649 sig += ( (L + 0.) * ampJ[0].imag() + (L + 1.) * ampJ[1].imag() );
653 }
else if (mode == 2) {
654 if (process == 0 || process == 1) {
655 amp[0] += (2. * L + 1.) * ampJ[0] * PlVec[L];
656 }
else if (process == 2) {
657 amp[0] += ((L + 0.) * ampJ[0] + double(L + 1.) * ampJ[1]) * PlVec[L];
658 amp[1] += complex(0., 1.) * (ampJ[1] - ampJ[0]) * sTheta * PlpVec[L];
665 if (mode == 0 || mode == 1) {
666 if (norm == 0) sig *= 4. * M_PI / k2 * CONVERT2MB;
667 else if (norm == 1) sig *= 64. * M_PI / s * CONVERT2MB;
669 }
else if (mode == 2) {
670 sig = (amp[0] * conj(amp[0])).real() + (amp[1] * conj(amp[1])).real();
671 if (norm == 0) sig *= 2. * M_PI / k2 * CONVERT2MB;
672 else if (norm == 1) sig *= 32. * M_PI / s * CONVERT2MB;
675 return ((idA == idB) ? 0.5 : 1.) * sig;
683 void SigmaPartialWave::legendreP(
double ct,
bool deriv) {
684 if (Lmax > 1) PlVec[1] = ct;
685 for (
int L = 2; L < Lmax; L++) {
686 PlVec[L] = ( (2. * L - 1.) * ct * PlVec[L - 1] -
687 (L - 1.) * PlVec[L - 2] ) /
double(L);
689 PlpVec[L] = ( (2. * L - 1.) * (PlVec[L - 1] + ct * PlpVec[L - 1]) -
690 (L - 1.) * PlpVec[L - 2] ) /
double(L);
704 bool HadronScatter::init(Info* infoPtrIn, Settings& settings,
705 Rndm *rndmPtrIn, ParticleData *particleDataPtr) {
711 doHadronScatter = settings.flag(
"HadronScatter:scatter");
712 afterDecay = settings.flag(
"HadronScatter:afterDecay");
713 allowDecayProd = settings.flag(
"HadronScatter:allowDecayProd");
714 scatterRepeat = settings.flag(
"HadronScatter:scatterRepeat");
716 hadronSelect = settings.mode(
"HadronScatter:hadronSelect");
717 Npar = settings.parm(
"HadronScatter:N");
718 kPar = settings.parm(
"HadronScatter:k");
719 pPar = settings.parm(
"HadronScatter:p");
721 scatterProb = settings.mode(
"HadronScatter:scatterProb");
722 jPar = settings.parm(
"HadronScatter:j");
723 rMax = settings.parm(
"HadronScatter:rMax");
725 doTile = settings.flag(
"HadronScatter:tile");
728 pTsigma = 2.0 * settings.parm(
"StringPT:sigma");
729 pTsigma2 = pTsigma * pTsigma;
730 double pT0ref = settings.parm(
"MultipartonInteractions:pT0ref");
731 double eCMref = settings.parm(
"MultipartonInteractions:eCMref");
732 double eCMpow = settings.parm(
"MultipartonInteractions:eCMpow");
733 double eCMnow = infoPtr->eCM();
734 pT0MPI = pT0ref * pow(eCMnow / eCMref, eCMpow);
737 double mp2 = particleDataPtr->m0(111) * particleDataPtr->m0(111);
738 double eA = infoPtr->eA();
739 double eB = infoPtr->eB();
740 double pzA = sqrt(eA * eA - mp2);
741 double pzB = -sqrt(eB * eB - mp2);
742 yMax = 0.5 * log((eA + pzA) / (eA - pzA));
743 yMin = 0.5 * log((eB + pzB) / (eB - pzB));
746 ytMax = int((yMax - yMin) / rMax);
747 ytSize = (yMax - yMin) /
double(ytMax);
748 ptMax = int(2. * M_PI / rMax);
749 ptSize = 2. * M_PI / double(ptMax);
752 ytSize = yMax - yMin;
758 for (
int yt = 0; yt < ytMax; yt++) tile[yt].resize(ptMax);
765 const char* PYTHIA8DATA =
"PYTHIA8DATA";
766 char* envPath = getenv(PYTHIA8DATA);
767 if (envPath != 0 && *envPath !=
'\0') {
769 while (*(envPath+i) !=
'\0') xmlPath += *(envPath+(i++));
770 }
else xmlPath =
"../xmldoc";
771 if (xmlPath[ xmlPath.length() - 1 ] !=
'/') xmlPath +=
"/";
774 if ( !sigmaPW[0].init(0, xmlPath,
"pipi-Froggatt.dat",
775 infoPtr, particleDataPtr, rndmPtr) )
return false;
776 if ( !sigmaPW[1].init(1, xmlPath,
"piK-Estabrooks.dat",
777 infoPtr, particleDataPtr, rndmPtr) )
return false;
778 if ( !sigmaPW[2].init(2, xmlPath,
"piN-SAID-WI08.dat",
779 infoPtr, particleDataPtr, rndmPtr) )
return false;
781 sigElMax = max(sigElMax, sigmaPW[0].getSigmaElMax());
782 sigElMax = max(sigElMax, sigmaPW[1].getSigmaElMax());
783 sigElMax = max(sigElMax, sigmaPW[2].getSigmaElMax());
796 void HadronScatter::debugOutput() {
798 cout <<
"Hadron scattering:" << endl
799 <<
" scatter = " << ((doHadronScatter) ?
"on" :
"off") << endl
800 <<
" afterDecay = " << ((afterDecay) ?
"on" :
"off") << endl
801 <<
" allowDecayProd = " << ((allowDecayProd) ?
"on" :
"off") << endl
802 <<
" scatterRepeat = " << ((scatterRepeat) ?
"on" :
"off") << endl
803 <<
" tile = " << ((doTile) ?
"on" :
"off") << endl
804 <<
" yMin = " << yMin << endl
805 <<
" yMax = " << yMax << endl
806 <<
" ytMax = " << ytMax << endl
807 <<
" ytSize = " << ytSize << endl
808 <<
" ptMax = " << ptMax << endl
809 <<
" ptSize = " << ptSize << endl
811 <<
" hadronSelect = " << hadronSelect << endl
812 <<
" N = " << Npar << endl
813 <<
" k = " << kPar << endl
814 <<
" p = " << pPar << endl
816 <<
" scatterProb = " << scatterProb << endl
817 <<
" j = " << jPar << endl
818 <<
" rMax = " << rMax << endl
820 <<
" pTsigma = " << pTsigma2 << endl
821 <<
" pT0MPI = " << pT0MPI << endl
823 <<
" sigElMax = " << sigElMax << endl << endl;
833 void HadronScatter::scatter(
Event& event) {
835 for (
int yt = 0; yt < ytMax; yt++)
836 for (
int pt = 0; pt < ptMax; pt++)
837 tile[yt][pt].clear();
840 for (
int i = 0; i <
event.size(); i++)
841 if (event[i].isFinal() &&
event[i].isHadron() && canScatter(event, i))
842 tile[yTile(event, i)][pTile(event, i)].insert(HSIndex(i, i));
845 vector < HadronScatterPair > scatterList;
847 for (
int pt1 = 0; pt1 < ptMax; pt1++)
848 for (
int yt1 = 0; yt1 < ytMax; yt1++)
849 for (set < HSIndex >::iterator si1 = tile[yt1][pt1].begin();
850 si1 != tile[yt1][pt1].end(); si1++)
851 tileIntProb(scatterList, event, *si1, yt1, pt1,
false);
853 sort(scatterList.rbegin(), scatterList.rend());
856 if (scatterRepeat) scattered.clear();
859 while (scatterList.size() > 0) {
861 HadronScatterPair &hsp = scatterList[0];
862 if (!event[hsp.i1.second].isFinal() || !
event[hsp.i2.second].isFinal()) {
863 scatterList.erase(scatterList.begin());
867 tile[hsp.yt1][hsp.pt1].erase(hsp.i1);
868 tile[hsp.yt2][hsp.pt2].erase(hsp.i2);
869 hadronScatter(event, hsp);
871 HSIndex iNew1 = hsp.i1, iNew2 = hsp.i2;
877 HSIndex iNew = iNew1;
878 for (
int i = 0; i < 2; i++) {
879 if (canScatter(event, iNew.second)) {
882 if (i == 1) scattered.insert(HSIndex(min(iNew1.first, iNew2.first),
883 max(iNew1.first, iNew2.first)));
884 int yt = yTile(event, iNew.second);
885 int pt = pTile(event, iNew.second);
886 resort = tileIntProb(scatterList, event, iNew, yt, pt,
true);
887 tile[yt][pt].insert(iNew);
895 scatterList.erase(scatterList.begin());
897 if (resort) sort(scatterList.rbegin(), scatterList.rend());
910 bool HadronScatter::canScatter(
Event& event,
int i) {
915 if (scatterProb == 1 || scatterProb == 2)
916 if (event[i].idAbs() != 111 && event[i].idAbs() != 211 &&
917 event[i].idAbs() != 321 && event[i].idAbs() != 2212)
921 switch (hadronSelect) {
923 double t1 = exp( - event[i].pT2() / 2. / pTsigma2);
924 double t2 = pow(pT0MPI, pPar) /
925 pow(pT0MPI * pT0MPI + event[i].pT2(), pPar / 2.);
926 p = Npar * t1 / ( (1 - kPar) * t1 + kPar * t2 );
931 if (p > rndmPtr->flat())
return true;
939 bool HadronScatter::doesScatter(
Event& event,
const HSIndex &i1,
941 Particle &p1 =
event[i1.second];
942 Particle &p2 =
event[i2.second];
945 if (!allowDecayProd && event[i1.first].mother1() ==
946 event[i2.first].mother1() &&
event[
event[i1.first].mother1()].isHadron())
951 scattered.find(HSIndex(min(i1.first, i2.first), max(i1.first, i2.first)))
952 != scattered.end())
return false;
955 int id1 = min(p1.idAbs(), p2.idAbs());
956 int id2 = max(p1.idAbs(), p2.idAbs());
957 if (scatterProb == 1 || scatterProb == 2) {
958 if ((id1 == 321 || id1 == 2212) && id1 == id2)
return false;
959 if (id1 == 321 && id2 == 2212)
return false;
963 double dy = p1.y() - p2.y();
964 double dp = abs(p1.phi() - p2.phi());
965 if (dp > M_PI) dp = 2 * M_PI - dp;
966 double dr2 = dy * dy + dp * dp;
967 double p = max(0., 1. - dr2 / rMax2);
970 if (scatterProb == 0 || scatterProb == 1) {
974 }
else if (scatterProb == 2) {
976 double Wcm = (p1.p() + p2.p()).mCalc();
978 if ((id1 == 111 || id1 == 211) && (id2 == 111 || id2 == 211)) pw = 0;
979 else if ((id1 == 111 || id1 == 211) && id2 == 321) pw = 1;
980 else if ((id1 == 111 || id1 == 211) && id2 == 2212) pw = 2;
982 infoPtr->errorMsg(
"Error in HadronScatter::doesScatter:"
983 "unknown subprocess");
985 if (!sigmaPW[pw].setSubprocess(p1.id(), p2.id())) {
986 infoPtr->errorMsg(
"Error in HadronScatter::doesScatter:"
987 "setSubprocess failed");
989 p *= 1 - exp(-jPar * sigmaPW[pw].sigmaEl(Wcm));
994 return (p > rndmPtr->flat());
1002 double HadronScatter::measure(
Event& event,
int idx1,
int idx2) {
1003 Particle &p1 =
event[idx1];
1004 Particle &p2 =
event[idx2];
1005 return abs(p1.pT() / p1.mT() - p2.pT() / p2.mT());
1013 bool HadronScatter::hadronScatter(
Event& event, HadronScatterPair &hsp) {
1014 bool doSwap = (0.5 < rndmPtr->flat()) ?
true :
false;
1015 if (doSwap) swap(hsp.i1, hsp.i2);
1017 Particle &p1 =
event[hsp.i1.second];
1018 Particle &p2 =
event[hsp.i2.second];
1021 double ct = 0., phi = 2 * M_PI * rndmPtr->flat();
1024 if (scatterProb == 0 || scatterProb == 1) {
1025 ct = 2. * rndmPtr->flat() - 1.;
1028 }
else if (scatterProb == 2) {
1030 int id1 = min(p1.idAbs(), p2.idAbs());
1031 int id2 = max(p1.idAbs(), p2.idAbs());
1032 double Wcm = (p1.p() + p2.p()).mCalc();
1034 if ((id1 == 111 || id1 == 211) && (id2 == 111 || id2 == 211)) pw = 0;
1035 else if ((id1 == 111 || id1 == 211) && id2 == 321) pw = 1;
1036 else if ((id1 == 111 || id1 == 211) && id2 == 2212) pw = 2;
1038 infoPtr->errorMsg(
"Error in HadronScatter::hadronScatter:"
1039 "unknown subprocess");
1041 sigmaPW[pw].setSubprocess(p1.id(), p2.id());
1042 ct = sigmaPW[pw].pickCosTheta(Wcm);
1047 sMat.toCMframe(p1.p(), p2.p());
1048 sMat.rot(acos(ct), phi);
1049 sMat.fromCMframe(p1.p(), p2.p());
1050 Vec4 v1 = p1.p(), v2 = p2.p();
1055 int iNew1 =
event.copy(hsp.i1.second, 98);
1057 event[iNew1].e(event[iNew1].eCalc());
1058 event[hsp.i1.second].statusNeg();
1059 int iNew2 =
event.copy(hsp.i2.second, 98);
1061 event[iNew2].e(event[iNew2].eCalc());
1062 event[hsp.i2.second].statusNeg();
1065 hsp.i1.second = iNew1;
1066 hsp.i2.second = iNew2;
1068 if (doSwap) swap(hsp.i1, hsp.i2);
1083 bool HadronScatter::tileIntProb(vector < HadronScatterPair > &hsp,
1084 Event &event,
const HSIndex &i1,
int yt1,
int pt1,
bool doAll) {
1086 bool pairAdded =
false;
1090 set < HSIndex >::iterator si2 = tile[yt1][pt1].find(i1);
1091 while (++si2 != tile[yt1][pt1].end())
1093 if (doesScatter(event, i1, *si2)) {
1094 double m = measure(event, i1.second, si2->second);
1095 hsp.push_back(HadronScatterPair(i1, yt1, pt1, *si2, yt1, pt1, m));
1100 int tileMax = (doAll) ? 9 : 4;
1101 for (
int tileNow = 0; tileNow < tileMax; tileNow++) {
1103 int yt2 = yt1, pt2 = pt1;
1105 case 0: ++pt2;
break;
1106 case 1: ++yt2; ++pt2;
break;
1107 case 2: ++yt2;
break;
1108 case 3: ++yt2; --pt2;
break;
1109 case 4: --pt2;
break;
1110 case 5: --yt2; --pt2;
break;
1111 case 6: --yt2;
break;
1112 case 7: --yt2; ++pt2;
break;
1116 if (yt2 < 0 || yt2 >= ytMax)
continue;
1120 if (pt2 == pt1 || pt2 == pt1 + 1)
continue;
1122 }
else if (pt2 >= ptMax) {
1124 if (pt2 == pt1 || pt2 == pt1 - 1)
continue;
1128 for (set < HSIndex >::iterator si2 = tile[yt2][pt2].begin();
1129 si2 != tile[yt2][pt2].end(); si2++) {
1131 if (doesScatter(event, i1, *si2)) {
1132 double m = measure(event, i1.second, si2->second);
1133 hsp.push_back(HadronScatterPair(i1, yt1, pt1, *si2, yt2, pt2, m));