PYTHIA  8.314
ColourReconnection.h
1 // ColourReconnection.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 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 // Header file for the Colour reconnection handling.
7 // Reconnect the colours between the partons before hadronization.
8 // It Contains the following classes:
9 // ColourDipole, ColourParticle, ColourJunction, ColourReconnection.
10 
11 #ifndef Pythia8_ColourReconnection_H
12 #define Pythia8_ColourReconnection_H
13 
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"
27 
28 namespace Pythia8 {
29 
30 //==========================================================================
31 
32 // Contain a single colour chain. It always start from a quark and goes to
33 // an anti quark or from an anti-junction to at junction
34 // (or possible combinations).
35 
36 class ColourDipole {
37 
38 public:
39 
40  // Constructor.
41  ColourDipole( int colIn = 0, int iColIn = 0, int iAcolIn = 0,
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.;}
47 
48  double mDip(Event & event) {
49  if (isJun || isAntiJun) return 1E9;
50  else return m(event[iCol].p(),event[iAcol].p());
51  }
52 
53  // Members.
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;
58  double p1p2;
59 
60  // Information for caching dipole momentum.
62  int ciCol{-1}, ciAcol{-1};
63  bool pCalculated{false};
64 
65  // Printing function, mainly intended for debugging.
66  void list() const;
67  long index{0};
68 
69 };
70 
71 // Comparison operator by index for two dipole pointers.
72 inline bool operator<(const ColourDipolePtr& d1, const ColourDipolePtr& d2) {
73  return ( d1 && d2 ? d1->index < d2->index : !d1 );}
74 
75 //==========================================================================
76 
77 // Junction class. In addition to the normal junction class, also contains a
78 // list of dipoles connected to it.
79 
80 class ColourJunction : public Junction {
81 
82 public:
83 
84  ColourJunction(const Junction& ju) : Junction(ju), dips(), dipsOrig() { }
85 
86  ColourJunction(const ColourJunction& ju) : Junction(Junction(ju)), dips(),
87  dipsOrig() { for(int i = 0;i < 3;++i) {
88  dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];}
89  }
90  ColourJunction& operator=( const ColourJunction& ju) {
91  this->Junction::operator=(ju);
92  for(int i = 0;i < 3;++i) {
93  dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];
94  }
95  return (*this);
96  }
97 
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];
102  void list() const;
103 
104 };
105 
106 //==========================================================================
107 
108 // TrialReconnection class.
109 
110 //--------------------------------------------------------------------------
111 
113 
114 public:
115 
116  TrialReconnection(ColourDipolePtr dip1In = 0, ColourDipolePtr dip2In = 0,
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;
122  }
123 
124  void list() {
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(); }
128  }
129 
130  vector<ColourDipolePtr> dips;
131  int mode;
132  double lambdaDiff;
133 
134 };
135 
136 //==========================================================================
137 
138 // ColourParticle class.
139 
140 //--------------------------------------------------------------------------
141 
142 class ColourParticle : public Particle {
143 
144 public:
145 
146  ColourParticle(const Particle& ju) : Particle(ju), isJun(), junKind() {}
147 
148  vector<vector<ColourDipolePtr> > dips;
149  vector<bool> colEndIncluded, acolEndIncluded;
150  vector<ColourDipolePtr> activeDips;
151  bool isJun;
152  int junKind;
153 
154  // Printing functions, intended for debugging.
155  void listParticle();
156  void listActiveDips();
157  void listDips();
158 
159 };
160 
161 //==========================================================================
162 
163 // The ColourReconnection class handles the colour reconnection.
164 
165 //--------------------------------------------------------------------------
166 
168 
169 public:
170 
171  // Constructor
172  ColourReconnection() : allowJunctions(), sameNeighbourCol(),
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() {}
179 
180  // Initialization.
181  bool init();
182 
183  // New beams possible for handling of hard diffraction.
184  void reassignBeamPtrs( BeamParticlePtr beamAPtrIn,
185  BeamParticlePtr beamBPtrIn) {
186  beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;}
187 
188  // Do colour reconnection for current event.
189  bool next( Event & event, int oldSize);
190 
191 private:
192 
193  // Constants: could only be changed in the code itself.
194  static const double MINIMUMGAIN, MINIMUMGAINJUN, TINYP1P2;
195  static const int MAXRECONNECTIONS;
196 
197  // Variables needed.
198  bool allowJunctions, sameNeighbourCol, singleReconOnly, lowerLambdaOnly;
199  bool allowDiqJunCR;
200  int nSys, nReconCols, swap1, swap2, reconnectMode, flipMode,
201  timeDilationMode;
202  double eCM, sCM, pT0, pT20Rec, pT0Ref, ecmRef, ecmPow, reconnectRange,
203  m0, mPseudo, m2Lambda, fracGluon, dLambdaCut, timeDilationPar,
204  timeDilationParGeV, tfrag, blowR, blowT, rHadron, kI;
205  double dipMaxDist;
206 
207  // List of current dipoles.
208  vector<ColourDipolePtr> dipoles, usedDipoles;
209 
210  // Last used dipole index, used for sorting.
211  int dipoleIndex{0};
212 
213  // Add a dipole and increment index.
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;
220  }
221 
222  // Add a dipole and increment index.
223  void addDipole(const ColourDipole& dipole) {
224  dipoles.push_back(make_shared<ColourDipole>(dipole));
225  dipoles.back()->index = ++dipoleIndex;
226  }
227 
228  // Lists of particles, junctions and trials.
229  vector<ColourJunction> junctions;
230  vector<ColourParticle> particles;
231  vector<TrialReconnection> junTrials, dipTrials;
232  vector<vector<int> > iColJun;
233  vector<double> formationTimes;
234 
235  // This is only to access the function call junctionRestFrame.
236  StringFragmentation stringFragmentation;
237 
238  // This class is used to calculate the string length.
239  StringLength stringLength;
240 
241  // Do colour reconnection for the event using the new model.
242  bool nextNew( Event & event, int oldSize);
243 
244  // Simple test swap between two dipoles.
245  void swapDipoles(ColourDipolePtr& dip1, ColourDipolePtr& dip2,
246  bool back = false);
247 
248  // Setup the dipoles.
249  void setupDipoles( Event& event, int iFirst = 0);
250 
251  // Form pseuparticle of a given dipole (or junction system).
252  void makePseudoParticle(ColourDipolePtr& dip, int status,
253  bool setupDone = false);
254 
255  // Find the indices in the particle list of the junction and also their
256  // respectively leg numbers.
257  bool getJunctionIndices(const ColourDipolePtr& dip, int &iJun, int &i0,
258  int &i1, int &i2, int &junLeg0, int &junLeg1, int &junLeg2) const;
259 
260  // Form all possible pseudoparticles.
261  void makeAllPseudoParticles(Event & event, int iFirst = 0);
262 
263  // Update all colours in the event.
264  void updateEvent( Event& event, int iFirst = 0);
265 
266  double calculateStringLength( const ColourDipolePtr& dip,
267  vector<ColourDipolePtr> & dips);
268 
269  // Calculate the string length for two event indices.
270  double calculateStringLength( int i, int j) const;
271 
272  // Calculate the length of a single junction
273  // given the 3 entries in the particle list.
274  double calculateJunctionLength(const int i, const int j, const int k) const;
275 
276  // Calculate the length of a double junction,
277  // given the 4 entries in the particle record.
278  // First two are expected to be the quarks and second two the anti quarks.
279  double calculateDoubleJunctionLength( const int i, const int j, const int k,
280  const int l) const;
281 
282  // Find all the particles connected in the junction.
283  // If a single junction, the size of iParticles should be 3.
284  // For multiple junction structures, the size will increase.
285  bool findJunctionParticles( int iJun, vector<int>& iParticles,
286  vector<bool> &usedJuns, int &nJuns, vector<ColourDipolePtr> &dips);
287 
288  // Do a single trial reconnection.
289  void singleReconnection( ColourDipolePtr& dip1, ColourDipolePtr& dip2);
290 
291  // Do a single trial reconnection to form a junction.
292  void singleJunction(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
293 
294  // Do a single trial reconnection to form a junction.
295  void singleJunction(const ColourDipolePtr& dip1, const ColourDipolePtr& dip2,
296  const ColourDipolePtr& dip3);
297 
298  // Print the chain containing the dipole.
299  void listChain(ColourDipolePtr& dip);
300 
301  // Print all the chains.
302  void listAllChains();
303 
304  // Print dipoles, intended for debuggning purposes.
305  void listDipoles( bool onlyActive = false, bool onlyReal = false) const;
306 
307  // Print particles, intended for debugging purposes.
308  void listParticles() const;
309 
310  // Print junctions, intended for debugging purposes.
311  void listJunctions() const;
312 
313  // Check that the current dipole setup is consistent. Debug purpose only.
314  void checkDipoles();
315 
316  // Check that the current dipole setup is consistent. Debug purpose only.
317  void checkRealDipoles(Event& event, int iFirst);
318 
319  // Calculate the invariant mass of a dipole.
320  double mDip(const ColourDipolePtr& dip) const;
321 
322  // Find a vertex of the (anti)-colour side of the dipole. If
323  // connected to a junction, recurse using the other junction
324  // dipoles.
325  Vec4 getVProd(const ColourDipolePtr& dip, bool anti) const;
326 
327  // Find an average vertex of the (anti)-colour sides of the dipoles
328  // connected to the given junction (not incuding the given dipole).
329  Vec4 getVProd(int iJun, const ColourDipolePtr& dip, bool anti) const;
330 
331  // Check that the distance between the impact parameter centers of
332  // the dipoles are within the allowed range of dipMaxDist.
333  bool checkDist(const ColourDipolePtr& dip1, const ColourDipolePtr& dip2)
334  const;
335 
336  // Find the neighbour to anti colour side. Return false if the dipole is
337  // connected to a junction or the new particle has a junction inside of it.
338  bool findAntiNeighbour(ColourDipolePtr& dip);
339 
340  // Find the neighbour to colour side. Return false if the dipole is
341  // connected to a junction or the new particle has a junction inside of it.
342  bool findColNeighbour(ColourDipolePtr& dip);
343 
344  // Check that trials do not contain junctions / unusable pseudoparticles.
345  bool checkJunctionTrials();
346 
347  // Store all dipoles connected to the ones used in the junction connection.
348  void storeUsedDips(TrialReconnection& trial);
349 
350  // Change colour structure to describe the reconnection in juncTrial.
351  bool doJunctionTrial(Event& event, TrialReconnection& juncTrial);
352 
353  // Change colour structure to describe the reconnection in juncTrial.
354  void doDipoleTrial(TrialReconnection& trial);
355 
356  // Change colour structure if it is three dipoles forming a junction system.
357  bool doTripleJunctionTrial(Event& event, TrialReconnection& juncTrial);
358 
359  // Calculate the difference between the old and new lambda.
360  double getLambdaDiff(const ColourDipolePtr& dip1,
361  const ColourDipolePtr& dip2, const ColourDipolePtr& dip3,
362  const ColourDipolePtr& dip4, const int mode) const;
363 
364  // Calculate the difference between the old and new lambda (dipole swap).
365  double getLambdaDiff(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
366 
367  // Update the list of dipole trial swaps to account for latest swap.
368  void updateDipoleTrials();
369 
370  // Update the list of dipole trial swaps to account for latest swap.
371  void updateJunctionTrials();
372 
373  // Check whether up to four dipoles are 'causally' connected.
374  bool checkTimeDilation(const ColourDipolePtr& dip1 = 0,
375  const ColourDipolePtr& dip2 = 0, const ColourDipolePtr& dip3 = 0,
376  const ColourDipolePtr& dip4 = 0) const;
377 
378  // Check whether two four momenta are casually connected.
379  bool checkTimeDilation(const Vec4& p1, const Vec4& p2, const double t1,
380  const double t2) const;
381 
382  // Find the momentum of the dipole.
383  Vec4 getDipoleMomentum(const ColourDipolePtr& dip) const;
384 
385  // Find all particles connected to a junction system (particle list).
386  void addJunctionIndices(const int iSinglePar, set<int> &iPar,
387  set<int> &usedJuncs) const;
388 
389  // Find all the formation times.
390  void setupFormationTimes( Event & event);
391 
392  // Get the mass of all partons connected to a junction system (event list).
393  double getJunctionMass(Event & event, int col);
394 
395  // Find all particles connected to a junction system (event list).
396  void addJunctionIndices(const Event & event, const int iSinglePar,
397  set<int> &iPar, set<int> &usedJuncs) const;
398 
399  // The old MPI-based scheme.
400  bool reconnectMPIs( Event& event, int oldSize);
401 
402  // Vectors and methods needed for the new gluon-move model.
403 
404  // Array of (indices of) all final coloured particles.
405  vector<int> iReduceCol, iExpandCol;
406 
407  // Array of all lambda distances between coloured partons.
408  int nColMove;
409  vector<double> lambdaijMove;
410 
411  // Function to return lambda value from array.
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)];
415  }
416 
417  // Function to return lambda(i,j) + lambda(i,k) - lambda(j,k).
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)];
423  }
424 
425  // The new gluon-move scheme.
426  bool reconnectMove( Event& event, int oldSize);
427 
428  // The common part for both Type I and II reconnections in e+e..
429  bool reconnectTypeCommon( Event& event, int oldSize);
430 
431  // The e+e- type I CR model.
432  map<double,pair<int,int> > reconnectTypeI( Event& event,
433  vector<vector<ColourDipole> > &dips, Vec4 decays[2]);
434  // bool reconnectTypeI( Event& event, int oldSize);
435 
436  // The e+e- type II CR model.
437  map<double,pair<int,int> > reconnectTypeII( Event& event,
438  vector<vector<ColourDipole> > &dips, Vec4 decays[2]);
439  // bool reconnectTypeII( Event& event, int oldSize);
440 
441  // calculate the determinant of 3 * 3 matrix.
442  double determinant3(vector<vector< double> >& vec);
443 };
444 
445 //==========================================================================
446 
447 } // end namespace Pythia8
448 
449 #endif // Pythia8_ColourReconnection_H
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: Event.h:338
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: Event.h:32
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
Definition: Basics.h:32
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