File indexing completed on 2025-01-18 10:06:19
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #ifndef Pythia8_ColourReconnection_H
0012 #define Pythia8_ColourReconnection_H
0013
0014 #include "Pythia8/Basics.h"
0015 #include "Pythia8/BeamParticle.h"
0016 #include "Pythia8/Event.h"
0017 #include "Pythia8/FragmentationFlavZpT.h"
0018 #include "Pythia8/Info.h"
0019 #include "Pythia8/ParticleData.h"
0020 #include "Pythia8/StringFragmentation.h"
0021 #include "Pythia8/PartonDistributions.h"
0022 #include "Pythia8/PartonSystems.h"
0023 #include "Pythia8/PythiaStdlib.h"
0024 #include "Pythia8/Settings.h"
0025 #include "Pythia8/StringLength.h"
0026 #include "Pythia8/StringInteractions.h"
0027
0028 namespace Pythia8 {
0029
0030
0031
0032
0033
0034
0035
0036 class ColourDipole {
0037
0038 public:
0039
0040
0041 ColourDipole( int colIn = 0, int iColIn = 0, int iAcolIn = 0,
0042 int colReconnectionIn = 0, bool isJunIn = false, bool isAntiJunIn = false,
0043 bool isActiveIn = true, bool isRealIn = false) : col(colIn), iCol(iColIn),
0044 iAcol(iAcolIn), colReconnection(colReconnectionIn), isJun(isJunIn),
0045 isAntiJun(isAntiJunIn),isActive(isActiveIn), isReal(isRealIn)
0046 {iColLeg = 0; iAcolLeg = 0; printed = false; p1p2 = 0.;}
0047
0048 double mDip(Event & event) {
0049 if (isJun || isAntiJun) return 1E9;
0050 else return m(event[iCol].p(),event[iAcol].p());
0051 }
0052
0053
0054 int col, iCol, iAcol, iColLeg, iAcolLeg, colReconnection;
0055 bool isJun, isAntiJun, isActive, isReal, printed;
0056 weak_ptr<ColourDipole> leftDip{}, rightDip{};
0057 vector<weak_ptr<ColourDipole> > colDips, acolDips;
0058 double p1p2;
0059
0060
0061 Pythia8::Vec4 dipoleMomentum;
0062 int ciCol{-1}, ciAcol{-1};
0063 bool pCalculated{false};
0064
0065
0066 void list() const;
0067 long index{0};
0068
0069 };
0070
0071
0072 inline bool operator<(const ColourDipolePtr& d1, const ColourDipolePtr& d2) {
0073 return ( d1 && d2 ? d1->index < d2->index : !d1 );}
0074
0075
0076
0077
0078
0079
0080 class ColourJunction : public Junction {
0081
0082 public:
0083
0084 ColourJunction(const Junction& ju) : Junction(ju), dips(), dipsOrig() { }
0085
0086 ColourJunction(const ColourJunction& ju) : Junction(Junction(ju)), dips(),
0087 dipsOrig() { for(int i = 0;i < 3;++i) {
0088 dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];}
0089 }
0090 ColourJunction& operator=( const ColourJunction& ju) {
0091 this->Junction::operator=(ju);
0092 for(int i = 0;i < 3;++i) {
0093 dips[i] = ju.dips[i]; dipsOrig[i] = ju.dipsOrig[i];
0094 }
0095 return (*this);
0096 }
0097
0098 ColourDipolePtr getColDip(int i) {return dips[i];}
0099 void setColDip(int i, ColourDipolePtr dip) {dips[i] = dip;}
0100 ColourDipolePtr dips[3];
0101 ColourDipolePtr dipsOrig[3];
0102 void list() const;
0103
0104 };
0105
0106
0107
0108
0109
0110
0111
0112 class TrialReconnection {
0113
0114 public:
0115
0116 TrialReconnection(ColourDipolePtr dip1In = 0, ColourDipolePtr dip2In = 0,
0117 ColourDipolePtr dip3In = 0, ColourDipolePtr dip4In = 0, int modeIn = 0,
0118 double lambdaDiffIn = 0) {
0119 dips.push_back(dip1In); dips.push_back(dip2In);
0120 dips.push_back(dip3In); dips.push_back(dip4In);
0121 mode = modeIn; lambdaDiff = lambdaDiffIn;
0122 }
0123
0124 void list() {
0125 cout << "mode: " << mode << " " << "lambdaDiff: " << lambdaDiff << endl;
0126 for (int i = 0;i < int(dips.size()) && dips[i] != 0;++i) {
0127 cout << " "; dips[i]->list(); }
0128 }
0129
0130 vector<ColourDipolePtr> dips;
0131 int mode;
0132 double lambdaDiff;
0133
0134 };
0135
0136
0137
0138
0139
0140
0141
0142 class ColourParticle : public Particle {
0143
0144 public:
0145
0146 ColourParticle(const Particle& ju) : Particle(ju), isJun(), junKind() {}
0147
0148 vector<vector<ColourDipolePtr> > dips;
0149 vector<bool> colEndIncluded, acolEndIncluded;
0150 vector<ColourDipolePtr> activeDips;
0151 bool isJun;
0152 int junKind;
0153
0154
0155 void listParticle();
0156 void listActiveDips();
0157 void listDips();
0158
0159 };
0160
0161
0162
0163
0164
0165
0166
0167 class ColourReconnection : public ColourReconnectionBase {
0168
0169 public:
0170
0171
0172 ColourReconnection() : allowJunctions(), sameNeighbourCol(),
0173 singleReconOnly(), lowerLambdaOnly(), allowDiqJunCR(), nSys(),
0174 nReconCols(), swap1(), swap2(), reconnectMode(), flipMode(),
0175 timeDilationMode(), eCM(), sCM(), pT0(), pT20Rec(), pT0Ref(), ecmRef(),
0176 ecmPow(), reconnectRange(), m0(), mPseudo(), m2Lambda(), fracGluon(),
0177 dLambdaCut(), timeDilationPar(), timeDilationParGeV(), tfrag(), blowR(),
0178 blowT(), rHadron(), kI(), dipMaxDist(), nColMove() {}
0179
0180
0181 bool init();
0182
0183
0184 void reassignBeamPtrs( BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn)
0185 {beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;}
0186
0187
0188 bool next( Event & event, int oldSize);
0189
0190 private:
0191
0192
0193 static const double MINIMUMGAIN, MINIMUMGAINJUN, TINYP1P2;
0194 static const int MAXRECONNECTIONS;
0195
0196
0197 bool allowJunctions, sameNeighbourCol, singleReconOnly, lowerLambdaOnly;
0198 bool allowDiqJunCR;
0199 int nSys, nReconCols, swap1, swap2, reconnectMode, flipMode,
0200 timeDilationMode;
0201 double eCM, sCM, pT0, pT20Rec, pT0Ref, ecmRef, ecmPow, reconnectRange,
0202 m0, mPseudo, m2Lambda, fracGluon, dLambdaCut, timeDilationPar,
0203 timeDilationParGeV, tfrag, blowR, blowT, rHadron, kI;
0204 double dipMaxDist;
0205
0206
0207 vector<ColourDipolePtr> dipoles, usedDipoles;
0208
0209
0210 int dipoleIndex{0};
0211
0212
0213 void addDipole(int colIn = 0, int iColIn = 0, int iAcolIn = 0,
0214 int colReconnectionIn = 0, bool isJunIn = false, bool isAntiJunIn = false,
0215 bool isActiveIn = true, bool isRealIn = false) {
0216 dipoles.push_back(make_shared<ColourDipole>(colIn, iColIn, iAcolIn,
0217 colReconnectionIn, isJunIn, isAntiJunIn, isActiveIn, isRealIn));
0218 dipoles.back()->index = ++dipoleIndex;
0219 }
0220
0221
0222 void addDipole(const ColourDipole& dipole) {
0223 dipoles.push_back(make_shared<ColourDipole>(dipole));
0224 dipoles.back()->index = ++dipoleIndex;
0225 }
0226
0227
0228 vector<ColourJunction> junctions;
0229 vector<ColourParticle> particles;
0230 vector<TrialReconnection> junTrials, dipTrials;
0231 vector<vector<int> > iColJun;
0232 vector<double> formationTimes;
0233
0234
0235 StringFragmentation stringFragmentation;
0236
0237
0238 StringLength stringLength;
0239
0240
0241 bool nextNew( Event & event, int oldSize);
0242
0243
0244 void swapDipoles(ColourDipolePtr& dip1, ColourDipolePtr& dip2,
0245 bool back = false);
0246
0247
0248 void setupDipoles( Event& event, int iFirst = 0);
0249
0250
0251 void makePseudoParticle(ColourDipolePtr& dip, int status,
0252 bool setupDone = false);
0253
0254
0255
0256 bool getJunctionIndices(const ColourDipolePtr& dip, int &iJun, int &i0,
0257 int &i1, int &i2, int &junLeg0, int &junLeg1, int &junLeg2) const;
0258
0259
0260 void makeAllPseudoParticles(Event & event, int iFirst = 0);
0261
0262
0263 void updateEvent( Event& event, int iFirst = 0);
0264
0265 double calculateStringLength( const ColourDipolePtr& dip,
0266 vector<ColourDipolePtr> & dips);
0267
0268
0269 double calculateStringLength( int i, int j) const;
0270
0271
0272
0273 double calculateJunctionLength(const int i, const int j, const int k) const;
0274
0275
0276
0277
0278 double calculateDoubleJunctionLength( const int i, const int j, const int k,
0279 const int l) const;
0280
0281
0282
0283
0284 bool findJunctionParticles( int iJun, vector<int>& iParticles,
0285 vector<bool> &usedJuns, int &nJuns, vector<ColourDipolePtr> &dips);
0286
0287
0288 void singleReconnection( ColourDipolePtr& dip1, ColourDipolePtr& dip2);
0289
0290
0291 void singleJunction(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
0292
0293
0294 void singleJunction(const ColourDipolePtr& dip1, const ColourDipolePtr& dip2,
0295 const ColourDipolePtr& dip3);
0296
0297
0298 void listChain(ColourDipolePtr& dip);
0299
0300
0301 void listAllChains();
0302
0303
0304 void listDipoles( bool onlyActive = false, bool onlyReal = false) const;
0305
0306
0307 void listParticles() const;
0308
0309
0310 void listJunctions() const;
0311
0312
0313 void checkDipoles();
0314
0315
0316 void checkRealDipoles(Event& event, int iFirst);
0317
0318
0319 double mDip(const ColourDipolePtr& dip) const;
0320
0321
0322
0323
0324 Vec4 getVProd(const ColourDipolePtr& dip, bool anti) const;
0325
0326
0327
0328 Vec4 getVProd(int iJun, const ColourDipolePtr& dip, bool anti) const;
0329
0330
0331
0332 bool checkDist(const ColourDipolePtr& dip1, const ColourDipolePtr& dip2)
0333 const;
0334
0335
0336
0337 bool findAntiNeighbour(ColourDipolePtr& dip);
0338
0339
0340
0341 bool findColNeighbour(ColourDipolePtr& dip);
0342
0343
0344 bool checkJunctionTrials();
0345
0346
0347 void storeUsedDips(TrialReconnection& trial);
0348
0349
0350 bool doJunctionTrial(Event& event, TrialReconnection& juncTrial);
0351
0352
0353 void doDipoleTrial(TrialReconnection& trial);
0354
0355
0356 bool doTripleJunctionTrial(Event& event, TrialReconnection& juncTrial);
0357
0358
0359 double getLambdaDiff(const ColourDipolePtr& dip1,
0360 const ColourDipolePtr& dip2, const ColourDipolePtr& dip3,
0361 const ColourDipolePtr& dip4, const int mode) const;
0362
0363
0364 double getLambdaDiff(ColourDipolePtr& dip1, ColourDipolePtr& dip2);
0365
0366
0367 void updateDipoleTrials();
0368
0369
0370 void updateJunctionTrials();
0371
0372
0373 bool checkTimeDilation(const ColourDipolePtr& dip1 = 0,
0374 const ColourDipolePtr& dip2 = 0, const ColourDipolePtr& dip3 = 0,
0375 const ColourDipolePtr& dip4 = 0) const;
0376
0377
0378 bool checkTimeDilation(const Vec4& p1, const Vec4& p2, const double t1,
0379 const double t2) const;
0380
0381
0382 Vec4 getDipoleMomentum(const ColourDipolePtr& dip) const;
0383
0384
0385 void addJunctionIndices(const int iSinglePar, set<int> &iPar,
0386 set<int> &usedJuncs) const;
0387
0388
0389 void setupFormationTimes( Event & event);
0390
0391
0392 double getJunctionMass(Event & event, int col);
0393
0394
0395 void addJunctionIndices(const Event & event, const int iSinglePar,
0396 set<int> &iPar, set<int> &usedJuncs) const;
0397
0398
0399 bool reconnectMPIs( Event& event, int oldSize);
0400
0401
0402
0403
0404 vector<int> iReduceCol, iExpandCol;
0405
0406
0407 int nColMove;
0408 vector<double> lambdaijMove;
0409
0410
0411 double lambda12Move( int i, int j) {
0412 int iAC = iReduceCol[i]; int jAC = iReduceCol[j];
0413 return lambdaijMove[nColMove * min( iAC, jAC) + max( iAC, jAC)];
0414 }
0415
0416
0417 double lambda123Move( int i, int j, int k) {
0418 int iAC = iReduceCol[i]; int jAC = iReduceCol[j]; int kAC = iReduceCol[k];
0419 return lambdaijMove[nColMove * min( iAC, jAC) + max( iAC, jAC)]
0420 + lambdaijMove[nColMove * min( iAC, kAC) + max( iAC, kAC)]
0421 - lambdaijMove[nColMove * min( jAC, kAC) + max( jAC, kAC)];
0422 }
0423
0424
0425 bool reconnectMove( Event& event, int oldSize);
0426
0427
0428 bool reconnectTypeCommon( Event& event, int oldSize);
0429
0430
0431 map<double,pair<int,int> > reconnectTypeI( Event& event,
0432 vector<vector<ColourDipole> > &dips, Vec4 decays[2]);
0433
0434
0435
0436 map<double,pair<int,int> > reconnectTypeII( Event& event,
0437 vector<vector<ColourDipole> > &dips, Vec4 decays[2]);
0438
0439
0440
0441 double determinant3(vector<vector< double> >& vec);
0442 };
0443
0444
0445
0446 }
0447
0448 #endif