PYTHIA  8.314
BeamSetup.h
1 // BeamSetup.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 // This file contains a helper class for the setup of beam flavour,
7 // kinematics and PDFs.
8 
9 #ifndef Pythia8_BeamSetup_H
10 #define Pythia8_BeamSetup_H
11 
12 #include "Pythia8/Basics.h"
13 #include "Pythia8/BeamParticle.h"
14 #include "Pythia8/BeamShape.h"
15 #include "Pythia8/HadronLevel.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/LesHouches.h"
18 #include "Pythia8/ParticleData.h"
19 #include "Pythia8/PartonDistributions.h"
20 #include "Pythia8/PartonLevel.h"
21 #include "Pythia8/PhysicsBase.h"
22 #include "Pythia8/ProcessLevel.h"
23 #include "Pythia8/PythiaStdlib.h"
24 #include "Pythia8/Settings.h"
25 
26 namespace Pythia8 {
27 
28 //==========================================================================
29 
30 // The BeamSetup class contains a number of routines auxiliary to Pythia
31 // to set up beam flavour, kinematics and PDFs.
32 
33 class BeamSetup : public PhysicsBase {
34 
35 public:
36 
37  // Constructor.
39  beamAPtr = make_shared<BeamParticle>();
40  beamBPtr = make_shared<BeamParticle>();
41  beamPomAPtr = make_shared<BeamParticle>();
42  beamPomBPtr = make_shared<BeamParticle>();
43  beamGamAPtr = make_shared<BeamParticle>();
44  beamGamBPtr = make_shared<BeamParticle>();
45  beamVMDAPtr = make_shared<BeamParticle>();
46  beamVMDBPtr = make_shared<BeamParticle>();
47  }
48 
49  // Possibility to pass in pointers to PDF's.
50  bool setPDFPtr( PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn,
51  PDFPtr pdfHardAPtrIn = nullptr, PDFPtr pdfHardBPtrIn = nullptr,
52  PDFPtr pdfPomAPtrIn = nullptr, PDFPtr pdfPomBPtrIn = nullptr,
53  PDFPtr pdfGamAPtrIn = nullptr, PDFPtr pdfGamBPtrIn = nullptr,
54  PDFPtr pdfHardGamAPtrIn = nullptr, PDFPtr pdfHardGamBPtrIn = nullptr,
55  PDFPtr pdfUnresAPtrIn = nullptr, PDFPtr pdfUnresBPtrIn = nullptr,
56  PDFPtr pdfUnresGamAPtrIn = nullptr, PDFPtr pdfUnresGamBPtrIn = nullptr,
57  PDFPtr pdfVMDAPtrIn = nullptr, PDFPtr pdfVMDBPtrIn = nullptr);
58  bool setPDFAPtr( PDFPtr pdfAPtrIn );
59  bool setPDFBPtr( PDFPtr pdfBPtrIn );
60 
61  // Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
62  bool setPhotonFluxPtr( PDFPtr photonFluxAIn, PDFPtr photonFluxBIn) {
63  if ( photonFluxAIn ) pdfGamFluxAPtr = photonFluxAIn;
64  if ( photonFluxBIn ) pdfGamFluxBPtr = photonFluxBIn;
65  return true;}
66 
67  // Possibility to pass in pointer to external LHA-interfaced generator.
68  bool setLHAupPtr( LHAupPtr lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;
69  useNewLHA = false; return true;}
70 
71  // For a given particle id, get a particle that represents its properties,
72  // i.e. a particle with the same PDF shape and parameters.
73  int represent(int idIn) const;
74 
75  // Switch to new beam particle identities; for similar hadrons only.
76  bool setBeamIDs( int idAin, int idBin = 0);
77 
78  // Switch beam kinematics.
79  bool setKinematics(double eCMIn);
80  bool setKinematics(double eAIn, double eBIn);
81  bool setKinematics(double pxAIn, double pyAIn, double pzAIn,
82  double pxBIn, double pyBIn, double pzBIn);
83  bool setKinematics(Vec4 pAIn, Vec4 pBIn);
84 
85  // Possibility to pass in pointer for beam shape.
86  bool setBeamShapePtr( BeamShapePtr beamShapePtrIn) {
87  beamShapePtr = beamShapePtrIn; return true;}
88 
89  // Possibility to access the pointer to the BeamShape object.
90  BeamShapePtr getBeamShapePtr() { return beamShapePtr; }
91 
92  // Return a parton density set among list of possibilities.
93  PDFPtr getPDFPtr(int idIn, int sequence = 1, string beam = "A",
94  bool resolved = true);
95 
96  // Return a map of the PDF pointers.
97  map<string, PDFPtr> getPDFPtr();
98 
99  // Set up frame of beams, Les Houches input, and switches for beam handling.
100  bool initFrame();
101 
102  // Initialize kinematics and PDFs of beams.
103  bool initBeams(bool doNonPertIn, StringFlav* flavSelPtr);
104 
105  // Return whether VMD states sampled.
106  bool getVMDsideA() { return doVMDsideA; }
107  bool getVMDsideB() { return doVMDsideB; }
108 
109  // Clear all beams.
110  void clear();
111 
112  // Pick new beam valence flavours (for pi0, eta, K0S, Pomeron, etc.).
113  void newValenceContent();
114 
115  // Recalculate kinematics for each event when beam momentum has a spread.
116  void nextKinematics();
117 
118  // Boost from CM frame to lab frame, or inverse. Set production vertex.
119  void boostAndVertex( Event& process, Event& event, bool toLab,
120  bool setVertex);
121 
122  // Print parton lists for the main beams. For debug mainly.
123  void list() const { beamAPtr->list(); beamBPtr->list(); }
124 
125  // Some data values are kept public so that the Pythia class can access them.
126  bool doLHA = false, useNewLHA = false, skipInit = false,
127  doMomentumSpread = {}, doVertexSpread = {}, doVarEcm = {},
128  allowIDAswitch = {}, hasSwitchedIDs = {}, beamA2gamma = {},
129  beamB2gamma = {};
130  int idA = {}, idB = {}, frameType = {}, boostType = {}, iPDFAsave = {},
131  gammaMode = {};
132  double mA = {}, mB = {}, pxA = {}, pxB = {}, pyA = {}, pyB = {}, pzA = {},
133  pzB = {}, eA = {}, eB = {}, pzAcm = {}, pzBcm = {}, eCM = {},
134  betaZ = {}, gammaZ = {};
135  Vec4 pAinit = {}, pBinit = {}, pAnow = {}, pBnow = {};
136  RotBstMatrix MfromCM = {}, MtoCM = {};
137  LHAupPtr lhaUpPtr = {};
138 
139  // The two incoming beams.
140  BeamParticlePtr beamAPtr = {};
141  BeamParticlePtr beamBPtr = {};
142 
143  // Alternative Pomeron beam-inside-beam.
144  BeamParticlePtr beamPomAPtr = {};
145  BeamParticlePtr beamPomBPtr = {};
146 
147  // Alternative photon beam-inside-beam.
148  BeamParticlePtr beamGamAPtr = {};
149  BeamParticlePtr beamGamBPtr = {};
150 
151  // Alternative VMD beam-inside-beam.
152  BeamParticlePtr beamVMDAPtr = {};
153  BeamParticlePtr beamVMDBPtr = {};
154 
155  // Hadron types for rapid switching.
156  vector<int> idAList = { 2212, 211, 311, 221,
157  331, 333, 411, 431, 443, 511, 531, 541, 553, 3212, 3312, 3334,
158  4112, 4312, 4332, 5112, 5312, 5332};
159 
160 protected:
161 
162  void onInitInfoPtr() override {
163  registerSubObject(*beamAPtr);
164  registerSubObject(*beamBPtr);
165  registerSubObject(*beamPomAPtr);
166  registerSubObject(*beamPomBPtr);
167  registerSubObject(*beamGamAPtr);
168  registerSubObject(*beamGamBPtr);
169  registerSubObject(*beamVMDAPtr);
170  registerSubObject(*beamVMDBPtr);
171  }
172 
173 private:
174 
175  // Initialization data, plus some event-specific.
176  bool doNonPert = {}, doDiffraction = {}, doSoftQCD = {},
177  doHardDiff = {}, doProcessLevel = {}, doPartonVertex = {},
178  doVertexPlane = {}, isUnresolvedA = {}, isUnresolvedB = {},
179  doVMDsideA = {}, doVMDsideB = {}, beamAResGamma = {},
180  beamBResGamma = {}, beamAUnresGamma = {}, beamBUnresGamma = {};
181 
182  // Pointers to the PDFs of beams, with several alternatives.
183  PDFPtr pdfAPtr = {}, pdfBPtr = {}, pdfHardAPtr = {}, pdfHardBPtr = {},
184  pdfPomAPtr = {}, pdfPomBPtr = {}, pdfGamAPtr = {}, pdfGamBPtr = {},
185  pdfHardGamAPtr = {}, pdfHardGamBPtr = {}, pdfUnresAPtr = {},
186  pdfUnresBPtr = {}, pdfUnresGamAPtr = {}, pdfUnresGamBPtr = {},
187  pdfGamFluxAPtr = {}, pdfGamFluxBPtr = {}, pdfVMDAPtr = {},
188  pdfVMDBPtr = {};
189 
190  // Array of PDFs to be used when idA can be changed between events.
191  vector<PDFPtr> pdfASavePtrs = {};
192 
193  // Pointer to BeamShape object for beam momentum and interaction vertex.
194  BeamShapePtr beamShapePtr = {};
195 
196  // Check that beams and beam combination can be handled.
197  bool checkBeams();
198 
199  // Calculate kinematics at initialization.
200  bool initKinematics();
201 
202  // Set up pointers to PDFs.
203  bool initPDFs();
204 
205  // Create an LHAPDF plugin PDF.
206  PDFPtr initLHAPDF(int idIn, string cfg);
207 
208 };
209 
210 //==========================================================================
211 
212 } // end namespace Pythia8
213 
214 #endif // Pythia8_BeamSetup_H
bool setPDFAPtr(PDFPtr pdfAPtrIn)
Routine to pass in pointers to PDF&#39;s. Usage optional.
Definition: BeamSetup.cc:105
Definition: PhysicsBase.h:27
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:56
Definition: BeamSetup.h:33
void newValenceContent()
Pick new beam valence flavours (for pi0, eta, K0S, Pomeron, etc.).
Definition: BeamSetup.cc:618
The Event class holds all info on the generated event.
Definition: Event.h:408
bool setKinematics(double eCMIn)
Switch beam kinematics.
Definition: BeamSetup.cc:258
bool setPhotonFluxPtr(PDFPtr photonFluxAIn, PDFPtr photonFluxBIn)
Set photon fluxes externally. Used with option "PDF:lepton2gammaSet = 2".
Definition: BeamSetup.h:62
void clear()
Clear all beams.
Definition: BeamSetup.cc:603
BeamParticlePtr beamAPtr
The two incoming beams.
Definition: BeamSetup.h:140
void onInitInfoPtr() override
Definition: BeamSetup.h:162
bool setPDFPtr(PDFPtr pdfAPtrIn, PDFPtr pdfBPtrIn, PDFPtr pdfHardAPtrIn=nullptr, PDFPtr pdfHardBPtrIn=nullptr, PDFPtr pdfPomAPtrIn=nullptr, PDFPtr pdfPomBPtrIn=nullptr, PDFPtr pdfGamAPtrIn=nullptr, PDFPtr pdfGamBPtrIn=nullptr, PDFPtr pdfHardGamAPtrIn=nullptr, PDFPtr pdfHardGamBPtrIn=nullptr, PDFPtr pdfUnresAPtrIn=nullptr, PDFPtr pdfUnresBPtrIn=nullptr, PDFPtr pdfUnresGamAPtrIn=nullptr, PDFPtr pdfUnresGamBPtrIn=nullptr, PDFPtr pdfVMDAPtrIn=nullptr, PDFPtr pdfVMDBPtrIn=nullptr)
Possibility to pass in pointers to PDF&#39;s.
Definition: BeamSetup.cc:21
bool initFrame()
Set up frame of beams, Les Houches input, and switches for beam handling.
Definition: BeamSetup.cc:342
bool doLHA
Some data values are kept public so that the Pythia class can access them.
Definition: BeamSetup.h:126
BeamSetup()
Constructor.
Definition: BeamSetup.h:38
map< string, PDFPtr > getPDFPtr()
Return a map of the PDF pointers.
Definition: BeamSetup.cc:1600
bool setBeamIDs(int idAin, int idBin=0)
Switch to new beam particle identities; for similar hadrons only.
Definition: BeamSetup.cc:180
bool setLHAupPtr(LHAupPtr lhaUpPtrIn)
Possibility to pass in pointer to external LHA-interfaced generator.
Definition: BeamSetup.h:68
bool getVMDsideA()
Return whether VMD states sampled.
Definition: BeamSetup.h:106
void list() const
Print parton lists for the main beams. For debug mainly.
Definition: BeamSetup.h:123
bool initBeams(bool doNonPertIn, StringFlav *flavSelPtr)
Initialize kinematics and PDFs of beams.
Definition: BeamSetup.cc:506
BeamParticlePtr beamVMDAPtr
Alternative VMD beam-inside-beam.
Definition: BeamSetup.h:152
Definition: Basics.h:254
bool setPDFBPtr(PDFPtr pdfBPtrIn)
Routine to pass in pointers to PDF&#39;s. Usage optional.
Definition: BeamSetup.cc:129
int represent(int idIn) const
Definition: BeamSetup.cc:151
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:87
BeamShapePtr getBeamShapePtr()
Possibility to access the pointer to the BeamShape object.
Definition: BeamSetup.h:90
vector< int > idAList
Hadron types for rapid switching.
Definition: BeamSetup.h:156
void boostAndVertex(Event &process, Event &event, bool toLab, bool setVertex)
Boost from CM frame to lab frame, or inverse. Set production vertex.
Definition: BeamSetup.cc:702
BeamParticlePtr beamGamAPtr
Alternative photon beam-inside-beam.
Definition: BeamSetup.h:148
bool setBeamShapePtr(BeamShapePtr beamShapePtrIn)
Possibility to pass in pointer for beam shape.
Definition: BeamSetup.h:86
Definition: Basics.h:32
BeamParticlePtr beamPomAPtr
Alternative Pomeron beam-inside-beam.
Definition: BeamSetup.h:144
void nextKinematics()
Recalculate kinematics for each event when beam momentum has a spread.
Definition: BeamSetup.cc:633