11 #ifndef Pythia8_ColourReconnection_H 12 #define Pythia8_ColourReconnection_H 14 #include "Pythia8/Basics.h" 15 #include "Pythia8/BeamParticle.h" 16 #include "Pythia8/Event.h" 17 #include "Pythia8/FragmentationFlavZpT.h" 18 #include "Pythia8/Info.h" 19 #include "Pythia8/ParticleData.h" 20 #include "Pythia8/StringFragmentation.h" 21 #include "Pythia8/PartonDistributions.h" 22 #include "Pythia8/PartonSystems.h" 23 #include "Pythia8/PythiaStdlib.h" 24 #include "Pythia8/Settings.h" 25 #include "Pythia8/StringLength.h" 26 #include "Pythia8/StringInteractions.h" 42 int colReconnectionIn = 0,
bool isJunIn =
false,
bool isAntiJunIn =
false,
43 bool isActiveIn =
true,
bool isRealIn =
false) :
col(colIn), iCol(iColIn),
44 iAcol(iAcolIn), colReconnection(colReconnectionIn), isJun(isJunIn),
45 isAntiJun(isAntiJunIn),isActive(isActiveIn), isReal(isRealIn)
46 {iColLeg = 0; iAcolLeg = 0; printed =
false; p1p2 = 0.;}
48 double mDip(
Event & event) {
49 if (isJun || isAntiJun)
return 1E9;
50 else return m(event[iCol].p(),event[iAcol].p());
54 int col, iCol, iAcol, iColLeg, iAcolLeg, colReconnection;
55 bool isJun, isAntiJun, isActive, isReal, printed;
56 weak_ptr<ColourDipole> leftDip{}, rightDip{};
57 vector<weak_ptr<ColourDipole> > colDips, acolDips;
62 int ciCol{-1}, ciAcol{-1};
63 bool pCalculated{
false};
72 inline bool operator<(
const ColourDipolePtr& d1,
const ColourDipolePtr& d2) {
73 return ( d1 && d2 ? d1->index < d2->index : !d1 );}
87 dipsOrig() {
for(
int i = 0;i < 3;++i) {
88 dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];}
91 this->Junction::operator=(ju);
92 for(
int i = 0;i < 3;++i) {
93 dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];
98 ColourDipolePtr getColDip(
int i) {
return dips[i];}
99 void setColDip(
int i, ColourDipolePtr dip) {dips[i] = dip;}
100 ColourDipolePtr dips[3];
101 ColourDipolePtr dipsOrig[3];
117 ColourDipolePtr dip3In = 0, ColourDipolePtr dip4In = 0,
int modeIn = 0,
118 double lambdaDiffIn = 0) {
119 dips.push_back(dip1In); dips.push_back(dip2In);
120 dips.push_back(dip3In); dips.push_back(dip4In);
121 mode = modeIn; lambdaDiff = lambdaDiffIn;
125 cout <<
"mode: " << mode <<
" " <<
"lambdaDiff: " << lambdaDiff << endl;
126 for (
int i = 0;i < int(dips.size()) && dips[i] != 0;++i) {
127 cout <<
" "; dips[i]->list(); }
130 vector<ColourDipolePtr> dips;
148 vector<vector<ColourDipolePtr> > dips;
149 vector<bool> colEndIncluded, acolEndIncluded;
150 vector<ColourDipolePtr> activeDips;
156 void listActiveDips();
173 singleReconOnly(), lowerLambdaOnly(), allowDiqJunCR(), nSys(),
174 nReconCols(), swap1(), swap2(), reconnectMode(), flipMode(),
175 timeDilationMode(), eCM(), sCM(), pT0(), pT20Rec(), pT0Ref(), ecmRef(),
176 ecmPow(), reconnectRange(), m0(), mPseudo(), m2Lambda(), fracGluon(),
177 dLambdaCut(), timeDilationPar(), timeDilationParGeV(), tfrag(), blowR(),
178 blowT(), rHadron(), kI(), dipMaxDist(), nColMove() {}
185 BeamParticlePtr beamBPtrIn) {
186 beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;}
189 bool next(
Event & event,
int oldSize);
194 static const double MINIMUMGAIN, MINIMUMGAINJUN, TINYP1P2;
195 static const int MAXRECONNECTIONS;
198 bool allowJunctions, sameNeighbourCol, singleReconOnly, lowerLambdaOnly;
200 int nSys, nReconCols, swap1, swap2, reconnectMode, flipMode,
202 double eCM, sCM, pT0, pT20Rec, pT0Ref, ecmRef, ecmPow, reconnectRange,
203 m0, mPseudo, m2Lambda, fracGluon, dLambdaCut, timeDilationPar,
204 timeDilationParGeV, tfrag, blowR, blowT, rHadron, kI;
208 vector<ColourDipolePtr> dipoles, usedDipoles;
214 void addDipole(
int colIn = 0,
int iColIn = 0,
int iAcolIn = 0,
215 int colReconnectionIn = 0,
bool isJunIn =
false,
bool isAntiJunIn =
false,
216 bool isActiveIn =
true,
bool isRealIn =
false) {
217 dipoles.push_back(make_shared<ColourDipole>(colIn, iColIn, iAcolIn,
218 colReconnectionIn, isJunIn, isAntiJunIn, isActiveIn, isRealIn));
219 dipoles.back()->index = ++dipoleIndex;
224 dipoles.push_back(make_shared<ColourDipole>(dipole));
225 dipoles.back()->index = ++dipoleIndex;
229 vector<ColourJunction> junctions;
230 vector<ColourParticle> particles;
231 vector<TrialReconnection> junTrials, dipTrials;
232 vector<vector<int> > iColJun;
233 vector<double> formationTimes;
242 bool nextNew(
Event & event,
int oldSize);
245 void swapDipoles(ColourDipolePtr& dip1, ColourDipolePtr& dip2,
249 void setupDipoles(
Event& event,
int iFirst = 0);
252 void makePseudoParticle(ColourDipolePtr& dip,
int status,
253 bool setupDone =
false);
257 bool getJunctionIndices(
const ColourDipolePtr& dip,
int &iJun,
int &i0,
258 int &i1,
int &i2,
int &junLeg0,
int &junLeg1,
int &junLeg2)
const;
261 void makeAllPseudoParticles(
Event & event,
int iFirst = 0);
264 void updateEvent(
Event& event,
int iFirst = 0);
266 double calculateStringLength(
const ColourDipolePtr& dip,
267 vector<ColourDipolePtr> & dips);
270 double calculateStringLength(
int i,
int j)
const;
274 double calculateJunctionLength(
const int i,
const int j,
const int k)
const;
279 double calculateDoubleJunctionLength(
const int i,
const int j,
const int k,
285 bool findJunctionParticles(
int iJun, vector<int>& iParticles,
286 vector<bool> &usedJuns,
int &nJuns, vector<ColourDipolePtr> &dips);
289 void singleReconnection( ColourDipolePtr& dip1, ColourDipolePtr& dip2);
292 void singleJunction(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
295 void singleJunction(
const ColourDipolePtr& dip1,
const ColourDipolePtr& dip2,
296 const ColourDipolePtr& dip3);
299 void listChain(ColourDipolePtr& dip);
302 void listAllChains();
305 void listDipoles(
bool onlyActive =
false,
bool onlyReal =
false)
const;
308 void listParticles()
const;
311 void listJunctions()
const;
317 void checkRealDipoles(
Event& event,
int iFirst);
320 double mDip(
const ColourDipolePtr& dip)
const;
325 Vec4 getVProd(
const ColourDipolePtr& dip,
bool anti)
const;
329 Vec4 getVProd(
int iJun,
const ColourDipolePtr& dip,
bool anti)
const;
333 bool checkDist(
const ColourDipolePtr& dip1,
const ColourDipolePtr& dip2)
338 bool findAntiNeighbour(ColourDipolePtr& dip);
342 bool findColNeighbour(ColourDipolePtr& dip);
345 bool checkJunctionTrials();
360 double getLambdaDiff(
const ColourDipolePtr& dip1,
361 const ColourDipolePtr& dip2,
const ColourDipolePtr& dip3,
362 const ColourDipolePtr& dip4,
const int mode)
const;
365 double getLambdaDiff(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
368 void updateDipoleTrials();
371 void updateJunctionTrials();
374 bool checkTimeDilation(
const ColourDipolePtr& dip1 = 0,
375 const ColourDipolePtr& dip2 = 0,
const ColourDipolePtr& dip3 = 0,
376 const ColourDipolePtr& dip4 = 0)
const;
379 bool checkTimeDilation(
const Vec4& p1,
const Vec4& p2,
const double t1,
380 const double t2)
const;
383 Vec4 getDipoleMomentum(
const ColourDipolePtr& dip)
const;
386 void addJunctionIndices(
const int iSinglePar, set<int> &iPar,
387 set<int> &usedJuncs)
const;
390 void setupFormationTimes(
Event & event);
393 double getJunctionMass(
Event & event,
int col);
396 void addJunctionIndices(
const Event & event,
const int iSinglePar,
397 set<int> &iPar, set<int> &usedJuncs)
const;
400 bool reconnectMPIs(
Event& event,
int oldSize);
405 vector<int> iReduceCol, iExpandCol;
409 vector<double> lambdaijMove;
412 double lambda12Move(
int i,
int j) {
413 int iAC = iReduceCol[i];
int jAC = iReduceCol[j];
414 return lambdaijMove[nColMove * min( iAC, jAC) + max( iAC, jAC)];
418 double lambda123Move(
int i,
int j,
int k) {
419 int iAC = iReduceCol[i];
int jAC = iReduceCol[j];
int kAC = iReduceCol[k];
420 return lambdaijMove[nColMove * min( iAC, jAC) + max( iAC, jAC)]
421 + lambdaijMove[nColMove * min( iAC, kAC) + max( iAC, kAC)]
422 - lambdaijMove[nColMove * min( jAC, kAC) + max( jAC, kAC)];
426 bool reconnectMove(
Event& event,
int oldSize);
429 bool reconnectTypeCommon(
Event& event,
int oldSize);
432 map<double,pair<int,int> > reconnectTypeI(
Event& event,
433 vector<vector<ColourDipole> > &dips,
Vec4 decays[2]);
437 map<double,pair<int,int> > reconnectTypeII(
Event& event,
438 vector<vector<ColourDipole> > &dips,
Vec4 decays[2]);
442 double determinant3(vector<vector< double> >& vec);
void list() const
Printing function, mainly intended for debugging.
Definition: ColourReconnection.cc:38
Definition: ColourReconnection.h:36
TrialReconnection class.
Definition: ColourReconnection.h:112
The Event class holds all info on the generated event.
Definition: Event.h:408
Definition: StringFragmentation.h:105
int col
Members.
Definition: ColourReconnection.h:54
Definition: StringInteractions.h:78
ColourDipole(int colIn=0, int iColIn=0, int iAcolIn=0, int colReconnectionIn=0, bool isJunIn=false, bool isAntiJunIn=false, bool isActiveIn=true, bool isRealIn=false)
Constructor.
Definition: ColourReconnection.h:41
StringLength class. It is used to calculate the lambda measure.
Definition: StringLength.h:23
Pythia8::Vec4 dipoleMomentum
Information for caching dipole momentum.
Definition: ColourReconnection.h:61
void reassignBeamPtrs(BeamParticlePtr beamAPtrIn, BeamParticlePtr beamBPtrIn)
New beams possible for handling of hard diffraction.
Definition: ColourReconnection.h:184
ColourParticle class.
Definition: ColourReconnection.h:142
Definition: ColourReconnection.h:80
The ColourReconnection class handles the colour reconnection.
Definition: ColourReconnection.h:167
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
double m(const Vec4 &v1)
Invariant mass and its square.
Definition: Basics.cc:588
bool operator<(const ColourDipolePtr &d1, const ColourDipolePtr &d2)
Comparison operator by index for two dipole pointers.
Definition: ColourReconnection.h:72
ColourReconnection()
Constructor.
Definition: ColourReconnection.h:172