Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-16 09:03:53

0001 // VinciaHistory.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2025 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 
0006 // Authors: Helen Brooks, Christian T Preuss
0007 // This file contains the VinciaHistory class.
0008 
0009 #ifndef VINCIA_History_H
0010 #define VINCIA_History_H
0011 
0012 #include "Pythia8/History.h"
0013 #include "Pythia8/Info.h"
0014 #include "Pythia8/Vincia.h"
0015 #include "Pythia8/VinciaAntennaFunctions.h"
0016 #include "Pythia8/VinciaMergingHooks.h"
0017 
0018 namespace Pythia8 {
0019 
0020 //==========================================================================
0021 
0022 // Convenient shorthand for storing ordered list of chains.
0023 
0024 struct PseudoChain {
0025   // Ordered list of concatenated chains.
0026   vector<int> chainlist;
0027   // Index unique up to chain content (not ordering).
0028   int index;
0029   // Charge Index.
0030   int cindex;
0031   // Does any of the chains contain an initial state parton.
0032   bool hasInitial;
0033   // Flavour at start of first and end of last.
0034   int flavStart;
0035   int flavEnd;
0036   // Charge.
0037   int charge;
0038 };
0039 
0040 //==========================================================================
0041 
0042 // A class for associating disconnected colour chains with either
0043 // resonances or the hard QCD scattering. Once fixed, the sector
0044 // history is deterministic.
0045 
0046 class ColourFlow {
0047 
0048  public:
0049 
0050   // Constructor.
0051   ColourFlow() : nChains(0), nBeamChainsMin(0), nBeamChainsMax(0), nRes(0) {
0052     for(int i=0; i<4; i++) {countChainsByChargeIndex[i]=0;
0053       countResByChargeIndex[i]=0;}}
0054 
0055   // Methods.
0056 
0057   // Add or assign chains.
0058   void addChain(int charge, int flavStart, int flavEnd, bool hasInitialIn);
0059   void selectResChains(int index,int iorder, int id);
0060   void selectBeamChains(int index,int iorder);
0061 
0062   // Initialise from hard process information.
0063   bool initHard(map<int, map<int,int> >& countRes,
0064     shared_ptr<VinciaMergingHooks> vinMergingHooksPtr);
0065 
0066   // Get information about current status of (pseudo)chains.
0067   bool checkChains(int cIndex);
0068   bool checkChains();
0069   int getNChainsLeft();
0070   int maxLength();
0071   int minLength();
0072 
0073   // Print information about (pseudo)chains.
0074   void print(bool printpsch = false);
0075 
0076   // Members.
0077 
0078   // Chains that arise from the decay of a resonance.
0079   map<int, vector<PseudoChain> > resChains;
0080 
0081   // Remaining list of ordered chains after resonances stripped off.
0082   vector<PseudoChain> beamChains;
0083 
0084   // Maps to all viable combinations of chains (of each charge)
0085   // from an index which depends on the identities of chains in it
0086   // - will be empty when resChains and beamChains set.
0087   map< int , vector< PseudoChain > > pseudochains;
0088 
0089   // Easy way to look up all pseudochains which contain a given chain.
0090   // Map from chain number to pseudochain index.
0091   map<int,vector<int> > chainToIndices;
0092 
0093   // Map from chain to the flavour at the start or end of that chain.
0094   map<int,int> chainStartToFlav;
0095   map<int,int> chainEndToFlav;
0096   // Keep track of chains which contain an initial state parton.
0097   map<int,bool> hasInitial;
0098   // Keep track of charge information of chains.
0099   map<int,int> chainToCharge;
0100 
0101  private:
0102 
0103   // Add or assign chains.
0104   void addChain(int oldIndex, int chainsIndex, int iChain,
0105     vector<int> & newIndices);
0106   void selectPseudochain(vector<int> & psch);
0107   void selectChain(int iChain);
0108 
0109   // Create empty vectors in resChains, and count each type.
0110   void addResonances(vector<int>& idsIn, map<int, map<int,int> >& idCounter,
0111     int charge, bool fc);
0112 
0113   // Convert charge information to an index.
0114   int getChargeIndex(int charge, bool fc);
0115 
0116   // Convert pseudochain to unique ID.
0117   int getID(PseudoChain& psc) {
0118     int id = 0;
0119     int iChains = psc.chainlist.size()-1;
0120     for (int iPwr(iChains); iPwr>=0; --iPwr)
0121       id += (psc.chainlist.at(iPwr)+1)*pow(10,iChains-iPwr);
0122     return id;
0123   }
0124 
0125   // List of found pseudochains to not double count.
0126   vector<int> pseudochainIDs;
0127 
0128   // Counters.
0129   int nChains;
0130   int nBeamChainsMin, nBeamChainsMax, nRes;
0131   map<int,int> countChainsByChargeIndex;
0132   map<int,int> countResByChargeIndex;
0133 
0134 };
0135 
0136 //==========================================================================
0137 
0138 // Class for a single step in the history of a process.
0139 
0140 class HistoryNode {
0141 
0142  public:
0143 
0144   // Constructors.
0145   HistoryNode() {};
0146   HistoryNode(Event& stateIn, vector< vector<int> > chainsIn,
0147     double scaleIn) : HistoryNode() {
0148     state = stateIn;
0149     clusterableChains = chainsIn;
0150     qEvolNow = scaleIn;
0151     hasRes = false;
0152     iRes = 0;
0153     idRes = 0;
0154     isInitPtr=false;
0155   };
0156 
0157   // Methods.
0158   void initPtr(VinciaCommon* vinComPtrIn, Resolution* resPtrIn,
0159     AntennaSetFSR* antSetPtrIn) {
0160     resPtr = resPtrIn;
0161     vinComPtr = vinComPtrIn;
0162     antSetFSRptr = antSetPtrIn;
0163     isInitPtr=true;
0164     nMinQQbar = 0;
0165   }
0166 
0167   // Set clusterList.
0168   int getNClusterings(shared_ptr<VinciaMergingHooks> vinMergingHooksPtr,
0169     Logger* loggerPtr, int verboseIn);
0170   void setClusterList(shared_ptr<VinciaMergingHooks> vinMergingHooksPtr,
0171     Logger* loggerPtr, int verboseIn);
0172 
0173   // Perform the clusterings according to the resolution criterion.
0174   bool cluster(HistoryNode& nodeClus,
0175     Logger* loggerPtr, int verboseIn);
0176 
0177   // Get energy fractions (used in PDF ratios).
0178   double xA() const {return 2. * state[3].e() / state[0].e();}
0179   double xB() const {return 2. * state[4].e() / state[0].e();}
0180 
0181   // Get flavours of beams (used in PDF ratios).
0182   int idA() const {return state[3].id();}
0183   int idB() const {return state[4].id();}
0184 
0185   // Get colour types of beams (used in PDF ratios).
0186   int colTypeA() const {return state[3].colType();}
0187   int colTypeB() const {return state[4].colType();}
0188 
0189   // Get evolution scale (used in trial shower).
0190   double getEvolNow() const {return qEvolNow;}
0191 
0192   // Setter methods.
0193   void setEvolScale(double scaleIn) {qEvolNow = scaleIn;}
0194 
0195   // Current state.
0196   Event state;
0197 
0198   // Resonance info.
0199   bool hasRes;
0200   int  iRes;
0201   int  idRes;
0202 
0203   // Minimal number of qqbar pairs.
0204   int nMinQQbar;
0205 
0206   // List of unclusterable lists of colour-connected partons.
0207   vector<vector<int>> clusterableChains;
0208 
0209   // Information corresponding to the last clustering.
0210   VinciaClustering lastClustering;
0211 
0212  private:
0213 
0214   // Perform a clustering.
0215   bool doClustering(VinciaClustering& clus, Event& clusEvent,
0216     vector<vector<int>>& clusChains, Logger* loggerPtr, int verboseIn);
0217 
0218   // Methods to calculate resolution and evolution scales.
0219   double calcResolution(VinciaClustering& clusIn) {
0220     return resPtr->q2sector(clusIn);}
0221   double calcEvolScale(VinciaClustering& clusIn) {
0222     return resPtr->q2evol(clusIn);}
0223 
0224   // Members.
0225 
0226   // Vincia pointers.
0227   Resolution*    resPtr;
0228   VinciaCommon*  vinComPtr;
0229   AntennaSetFSR* antSetFSRptr;
0230 
0231   bool isInitPtr;
0232 
0233   // The value of the evolution scale.
0234   double qEvolNow;
0235 
0236   // List of next possible clusterings.
0237   // Map is from corresponding resolution criterion.
0238   map<double, VinciaClustering> clusterList;
0239 
0240 };
0241 typedef map<int, vector<HistoryNode> > HistoryNodes;
0242 
0243 //==========================================================================
0244 
0245 // History class for the Vincia shower.
0246 
0247 class VinciaHistory {
0248 
0249 public:
0250 
0251   // Constructor.
0252   VinciaHistory(Event &stateIn,
0253     BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
0254     MergingHooksPtr mergingHooksPtrIn,
0255     PartonLevel* trialPartonLevelPtrIn,
0256     ParticleData* particleDataPtrIn,
0257     Info* infoPtrIn);
0258 
0259   // Check if a valid history was found.
0260   bool isValid() {return foundValidHistory;}
0261 
0262   // Check if history failed merging scale cut.
0263   bool isBelowMS() {return failedMSCut;}
0264 
0265   // Perform a trial shower and get the CKKW-L weight
0266   // (excluding last no-emission probability that is
0267   // calculated in the main shower).
0268   double getWeightCKKWL();
0269 
0270   // Find the first clustered state above the merging scale.
0271   Event getFirstClusteredEventAboveTMS();
0272 
0273   // What was the multiplicity of this event?
0274   int getNClusterSteps();
0275 
0276   // Get the shower starting scale.
0277   double getRestartScale();
0278 
0279   // Should we overwrite the original event?
0280   // E.g. if an MPI was generated.
0281   bool  hasNewProcess() {return hasNewProcessSav;}
0282   Event getNewProcess() {return newProcess;}
0283   bool  doAbort() {return aborted;}
0284 
0285 private:
0286 
0287   // Loop over colPerms and set historyBest.
0288   void findBestHistory();
0289 
0290   // Find all viable colour orderings and return number.
0291   unsigned int countPerms();
0292 
0293   // Fetch the colour chains from the hard event.
0294   bool getColChains();
0295 
0296   // Find all selections of chains for colour singlets.
0297   bool assignResChains(map<int, map<int,int>>& idCounter,
0298     vector<ColourFlow>& flowsSoFar);
0299 
0300   // Find all selections of chains for colour singlets.
0301   bool assignResFromEvent(map<int, map<int,int>>& idCounter,
0302     vector<ColourFlow>& flowsSoFar);
0303 
0304   // Find all selections of chains for everything that is
0305   // not a colour singlet resonance.
0306   bool assignBeamChains(vector<ColourFlow>& flowsSoFar);
0307 
0308   // Make a single selection in all possible ways.
0309   bool assignNext(vector<ColourFlow>& flowsSoFar, bool isRes = false,
0310     int id = 0, int cIndex = 0);
0311 
0312   // Make a specific selection for a resonance.
0313   bool assignThis(vector<ColourFlow>& flowsSoFar, int id, int cIndex,
0314     vector<int>& chains);
0315 
0316   // Is this (completed flow) compatible with the number of
0317   // resonances of each type?
0318   bool check(ColourFlow& flow);
0319 
0320   // Construct history for a given colour permutation.
0321   tuple<bool, double, HistoryNodes> findHistoryPerm(ColourFlow& flow);
0322 
0323   // Check if history failed merging scale cut.
0324   bool checkMergingCut(HistoryNodes& history);
0325 
0326   // Initialise history nodes for each system.
0327   HistoryNodes initHistoryNodes(ColourFlow& flow );
0328 
0329   // Translate abstract book-keeping of colour ordering into
0330   // systems of particles.
0331   map<int, vector<vector<int>>> getSystems(ColourFlow& flow,
0332     map<int, int>& sysToRes);
0333 
0334   // Decide if state is the Born topology.
0335   bool isBorn(const HistoryNode& nodeIn, bool isRes);
0336 
0337   // Set up beams for given history node.
0338   bool setupBeams(const HistoryNode* node, double scale2);
0339 
0340   // Calculate criterion for testing whether to keep history.
0341   double calcME2guess(vector<HistoryNode>& history,bool isRes);
0342   // Calculate ME for Born-level event.
0343   double calcME2Born(const HistoryNode& bornNode, bool isRes);
0344 
0345   // Calculate the antenna function for a given clustering.
0346   double calcAntFun(const VinciaClustering& clusNow);
0347 
0348   // Calculate PDF ratio to multiply CKKW-L weight.
0349   double calcPDFRatio(const HistoryNode* nodeNow,
0350     double pT2now, double pT2next);
0351 
0352   // Calculate alphaS ratio to multiply CKKW-L weight.
0353   double calcAlphaSRatio(const HistoryNode& node);
0354 
0355   // Return the kinematic maximum for this event.
0356   double getStartScale(Event& event, bool isRes);
0357 
0358   // Perform a trial branching and return scale.
0359   double qNextTrial(double qStart, Event& evtIn);
0360 
0361   // Print colour chains.
0362   void printChains();
0363 
0364   // Verbosity.
0365   int verbose;
0366 
0367   // Beams -- for PDFs.
0368   BeamParticle beamA, beamB;
0369 
0370   // MergingHooks.
0371   shared_ptr<VinciaMergingHooks> vinMergingHooksPtr{};
0372 
0373   // PartonLevel pointer.
0374   PartonLevel* trialPartonLevel{};
0375 
0376   // Particle data pointer.
0377   ParticleData* particleDataPtr{};
0378 
0379   // Other Pythia pointers.
0380   Info* infoPtr{};
0381   Logger* loggerPtr{};
0382 
0383   // Vincia pointers.
0384   shared_ptr<VinciaFSR> fsrShowerPtr{};
0385   shared_ptr<VinciaISR> isrShowerPtr{};
0386   MECs*                 mecsPtr{};
0387   VinciaCommon*         vinComPtr{};
0388   Resolution*           resPtr{};
0389   AntennaSetFSR*        antSetFSRptr{};
0390 
0391   // Check we found a valid history at all.
0392   bool foundValidHistory;
0393 
0394   // Check we vetoed due failing the merging scale cut.
0395   bool failedMSCut;
0396 
0397   // This is the best PS history so far.
0398   HistoryNodes historyBest;
0399 
0400   // Criterion to minimise (best so far).
0401   // Call calcME2guess(..) to evaluate.
0402   double ME2guessBest;
0403 
0404   // Initial colour permutation, stored as
0405   // list of sets of colour connected partons.
0406   vector<vector<int>> colChainsSav;
0407 
0408   // Track if chain contains initial state quarks.
0409   map<int, bool> chainHasInitial;
0410 
0411   // Keep track of resonances in the event record.
0412   map<int, vector<int>> resIDToIndices;
0413   map<int, vector<int>> resIndexToChains;
0414 
0415   // All colour flows compatible with Born process.
0416   vector<ColourFlow> colPerms;
0417 
0418   // ME generated event.
0419   Event state;
0420 
0421   // Number of quarks  and gluon loops in event.
0422   int nQSave, nGluonLoopsSave;
0423 
0424   // The merging scale and whether it is the evolution variable.
0425   double qms;
0426   bool   msIsEvolVar;
0427 
0428   // The maximum multiplicity of our ME generator.
0429   int nMax, nMaxRes;
0430 
0431   // Possible new hard process info (if MPI was generated).
0432   bool   hasNewProcessSav;
0433   Event  newProcess;
0434   double newProcessScale;
0435 
0436   // Flag to signal if something went wrong.
0437   bool aborted;
0438 
0439 };
0440 
0441 //==========================================================================
0442 
0443 } // end namespace Pythia8
0444 
0445 #endif // Pythia8_VinciaHistory_H