Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:28

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