Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-17 09:08:28

0001 // ProcessContainer.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 // This file contains the collected machinery of a process.
0007 // ProcessContainer: contains information on a particular process.
0008 // SetupContainers: administrates the selection/creation of processes.
0009 
0010 #ifndef Pythia8_ProcessContainer_H
0011 #define Pythia8_ProcessContainer_H
0012 
0013 #include "Pythia8/Basics.h"
0014 #include "Pythia8/BeamParticle.h"
0015 #include "Pythia8/Event.h"
0016 #include "Pythia8/Info.h"
0017 #include "Pythia8/ParticleData.h"
0018 #include "Pythia8/PartonDistributions.h"
0019 #include "Pythia8/PhaseSpace.h"
0020 #include "Pythia8/PhysicsBase.h"
0021 #include "Pythia8/PythiaStdlib.h"
0022 #include "Pythia8/ResonanceDecays.h"
0023 #include "Pythia8/Settings.h"
0024 #include "Pythia8/SigmaProcess.h"
0025 #include "Pythia8/SigmaTotal.h"
0026 #include "Pythia8/SigmaOnia.h"
0027 #include "Pythia8/StandardModel.h"
0028 #include "Pythia8/SusyCouplings.h"
0029 #include "Pythia8/SLHAinterface.h"
0030 #include "Pythia8/UserHooks.h"
0031 
0032 namespace Pythia8 {
0033 
0034 //==========================================================================
0035 
0036 // The ProcessContainer class combines pointers to matrix element and
0037 // phase space generator with general generation info.
0038 
0039 class ProcessContainer : public PhysicsBase {
0040 
0041 public:
0042 
0043   // Constructor.
0044   ProcessContainer(SigmaProcessPtr sigmaProcessPtrIn = 0,
0045     PhaseSpacePtr phaseSpacePtrIn = 0) :
0046       sigmaProcessPtr(sigmaProcessPtrIn),
0047       phaseSpacePtr(phaseSpacePtrIn), resDecaysPtr(), gammaKinPtr(),
0048       matchInOut(), idRenameBeams(), setLifetime(), setQuarkMass(),
0049       setLeptonMass(), idNewM(), smearHadronMass(), idSmearHadIn(),
0050       idSmearHadrons(), mRecalculate(), mNewM(), isLHA(), isNonDiff(),
0051       isResolved(), isDiffA(), isDiffB(), isDiffC(), isQCD3body(),
0052       allowNegSig(), isSameSave(), increaseMaximum(), canVetoResDecay(),
0053       lhaStrat(), lhaStratAbs(), processCode(), useStrictLHEFscales(),
0054       isAsymLHA(), betazLHA(), newSigmaMx(), nTry(), nSel(), nAcc(),
0055       nTryStat(), sigmaMx(), sigmaSgn(), sigmaSum(), sigma2Sum(), sigmaNeg(),
0056       sigmaAvg(), sigmaFin(), deltaFin(), weightNow(), wtAccSum(),
0057       beamAhasResGamma(), beamBhasResGamma(), beamHasResGamma(),
0058       beamHasGamma(), beamAgammaMode(), beamBgammaMode(), gammaModeEvent(),
0059       approximatedGammaFlux(), doMerging(), nTryRequested(), nSelRequested(),
0060       nAccRequested(), sigmaTemp(), sigma2Temp(), normVar3() {}
0061 
0062   // Initialize phase space and counters.
0063   bool init(bool isFirst, ResonanceDecays* resDecaysPtrIn,
0064     SLHAinterface* slhaInterfacePtr, GammaKinematics* gammaKinPtrIn);
0065 
0066   // Store or replace Les Houches pointer.
0067   void setLHAPtr( LHAupPtr lhaUpPtrIn,  ParticleData* particleDataPtrIn = 0,
0068     Settings* settingsPtrIn = 0, Rndm* rndmPtrIn = 0)
0069     {lhaUpPtr = lhaUpPtrIn; setLifetime = 0;
0070     if (settingsPtrIn && rndmPtrIn) {
0071       rndmPtr = rndmPtrIn;
0072       setLifetime = settingsPtrIn->mode("LesHouches:setLifetime");
0073     }
0074     if (particleDataPtrIn != 0) particleDataPtr = particleDataPtrIn;
0075     if (sigmaProcessPtr != 0) sigmaProcessPtr->setLHAPtr(lhaUpPtr);
0076     if (phaseSpacePtr != 0) phaseSpacePtr->setLHAPtr(lhaUpPtr);}
0077 
0078   // Switch to new beam particle identities; for similar hadrons only.
0079   void updateBeamIDs() {phaseSpacePtr->updateBeamIDs();}
0080 
0081   // Update the CM energy of the event.
0082   void newECM(double eCM) {phaseSpacePtr->newECM(eCM);}
0083 
0084   // Generate a trial event; accepted or not.
0085   bool trialProcess();
0086 
0087   // Pick flavours and colour flow of process.
0088   bool constructState();
0089 
0090   // Give the hard subprocess (with option for a second hard subprocess).
0091   bool constructProcess( Event& process, bool isHardest = true);
0092 
0093   // Give the Les Houches decay chain for external resonance.
0094   bool constructDecays( Event& process);
0095 
0096   // Do resonance decays.
0097   bool decayResonances( Event& process);
0098 
0099   // Accumulate statistics after user veto.
0100   void accumulate();
0101 
0102   // Reset statistics on events generated so far.
0103   void reset();
0104 
0105   // Set whether (photon) beam is resolved or unresolved.
0106   // Method propagates the choice of photon process type to beam pointers.
0107   void setBeamModes(bool setVMD = false, bool isSampled = true);
0108 
0109   // Process name and code, and the number of final-state particles.
0110   string name()             const {return sigmaProcessPtr->name();}
0111   int    code()             const {return sigmaProcessPtr->code();}
0112   int    nFinal()           const {return sigmaProcessPtr->nFinal();}
0113   bool   isSUSY()           const {return sigmaProcessPtr->isSUSY();}
0114   bool   isNonDiffractive() const {return isNonDiff;}
0115   bool   isSoftQCD()        const {return (code() > 100 && code() < 107);}
0116   bool   isElastic()        const {return (code() == 102);}
0117 
0118   // Member functions for info on generation process.
0119   bool   newSigmaMax() const {return newSigmaMx;}
0120   double sigmaMax()    const {return sigmaMx;}
0121   long   nTried()      const {return nTry;}
0122   long   nSelected()   const {return nSel;}
0123   long   nAccepted()   const {return nAcc;}
0124   double weightSum()   const {return wtAccSum;}
0125   double sigmaSelMC( bool doAccumulate = true)
0126     { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaAvg;}
0127   double sigmaMC(    bool doAccumulate = true)
0128     { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return sigmaFin;}
0129   double deltaMC(    bool doAccumulate = true)
0130     { if (nTry > nTryStat && doAccumulate) sigmaDelta(); return deltaFin;}
0131 
0132   // New value for switched beam identity or energy (for SoftQCD processes).
0133   double sigmaMaxSwitch() {sigmaMx = phaseSpacePtr->sigmaMaxSwitch();
0134     return sigmaMx;}
0135 
0136   // Some kinematics quantities.
0137   int    id1()         const {return sigmaProcessPtr->id(1);}
0138   int    id2()         const {return sigmaProcessPtr->id(2);}
0139   double x1()          const {return phaseSpacePtr->x1();}
0140   double x2()          const {return phaseSpacePtr->x2();}
0141   double Q2Fac()       const {return sigmaProcessPtr->Q2Fac();}
0142   double mHat()        const {return sqrtpos(phaseSpacePtr->sHat());}
0143   double pTHat()       const {return phaseSpacePtr->pTHat();}
0144 
0145   // Tell whether container is for Les Houches events.
0146   bool   isLHAContainer() const {return isLHA;}
0147   int    lhaStrategy()    const {return lhaStrat;}
0148 
0149   // Info on Les Houches events.
0150   int    codeLHASize()       const {return codeLHA.size();}
0151   int    subCodeLHA(int i)   const {return codeLHA[i];}
0152   long   nTriedLHA(int i)    const {return nTryLHA[i];}
0153   long   nSelectedLHA(int i) const {return nSelLHA[i];}
0154   long   nAcceptedLHA(int i) const {return nAccLHA[i];}
0155 
0156   // When two hard processes set or get info whether process is matched.
0157   void   isSame( bool isSameIn) { isSameSave = isSameIn;}
0158   bool   isSame()      const {return isSameSave;}
0159 
0160 private:
0161 
0162   // Constants: could only be changed in the code itself.
0163   static const int N12SAMPLE, N3SAMPLE;
0164 
0165   // Pointer to the subprocess matrix element. Mark if external.
0166   SigmaProcessPtr  sigmaProcessPtr;
0167 
0168   // Pointer to the phase space generator.
0169   PhaseSpacePtr    phaseSpacePtr;
0170 
0171   // Pointer to ResonanceDecays object for sequential resonance decays.
0172   ResonanceDecays* resDecaysPtr;
0173 
0174   // Pointer to LHAup for generating external events.
0175   LHAupPtr         lhaUpPtr;
0176 
0177   // Pointer to the phase space generator of photons from leptons.
0178   GammaKinematics* gammaKinPtr;
0179 
0180   // Possibility to modify Les Houches input.
0181   bool   matchInOut;
0182   int    idRenameBeams, setLifetime, setQuarkMass, setLeptonMass, idNewM[9],
0183          smearHadronMass;
0184   vector<int> idSmearHadIn, idSmearHadrons;
0185   double mRecalculate, mNewM[9];
0186 
0187   // Info on process.
0188   bool   isLHA, isNonDiff, isResolved, isDiffA, isDiffB, isDiffC, isQCD3body,
0189          allowNegSig, isSameSave, increaseMaximum, canVetoResDecay;
0190   int    lhaStrat, lhaStratAbs, processCode;
0191   bool   useStrictLHEFscales;
0192 
0193   // Boost Les Houches events to CM frame (when originally asymmetric).
0194   bool   isAsymLHA;
0195   double betazLHA;
0196 
0197   // Statistics on generation process. (Long integers just in case.)
0198   bool   newSigmaMx;
0199   long   nTry, nSel, nAcc, nTryStat;
0200   double sigmaMx, sigmaSgn, sigmaSum, sigma2Sum, sigmaNeg, sigmaAvg,
0201          sigmaFin, deltaFin, weightNow, wtAccSum;
0202 
0203   // Flags to store whether beam has a (un)resolved photon.
0204   bool   beamAhasResGamma, beamBhasResGamma, beamHasResGamma, beamHasGamma;
0205   int    beamAgammaMode, beamBgammaMode, gammaModeEvent;
0206 
0207   // Use approximated photon flux for process sampling.
0208   bool   approximatedGammaFlux;
0209 
0210   // Check if merging is enabled.
0211   bool   doMerging;
0212 
0213   // Statistics for Les Houches event classification.
0214   vector<int> codeLHA;
0215   vector<long> nTryLHA, nSelLHA, nAccLHA;
0216 
0217   // More fine-grained event counting.
0218   long nTryRequested, nSelRequested, nAccRequested;
0219 
0220   // Temporary summand for handling (weighted) events when vetoes are applied.
0221   double sigmaTemp, sigma2Temp, normVar3;
0222 
0223   // Estimate integrated cross section and its uncertainty.
0224   void sigmaDelta();
0225 
0226 };
0227 
0228 //==========================================================================
0229 
0230 // The SetupContainers class turns the list of user-requested processes
0231 // into a vector of ProcessContainer objects, each with a process.
0232 
0233 class SetupContainers {
0234 
0235 public:
0236 
0237   // Constructor.
0238   SetupContainers() : nVecA(), nVecB() {}
0239 
0240   // Initialization assuming all necessary data already read.
0241   bool init(vector<ProcessContainer*>& containerPtrs, Info* infoPtr);
0242 
0243   // Initialization of a second hard process.
0244   bool init2(vector<ProcessContainer*>& container2Ptrs, Info* infoPtr);
0245 
0246 private:
0247 
0248   // Methods to check that outgoing SUSY particles are allowed ones.
0249   void setupIdVecs( Settings& settings);
0250   bool allowIdVals( int idCheck1, int idCheck2);
0251 
0252   // Arrays of allowed outgoing SUSY particles and their lengths.
0253   vector<int> idVecA, idVecB;
0254   int nVecA, nVecB;
0255 
0256   // Helper class to setup onia production.
0257   SigmaOniaSetup charmonium, bottomonium;
0258 
0259 };
0260 
0261 //==========================================================================
0262 
0263 } // end namespace Pythia8
0264 
0265 #endif // Pythia8_ProcessContainer_H