PYTHIA  8.314
BeamParticle.h
1 // BeamParticle.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 // Header file for information on incoming beams.
7 // ResolvedParton: an initiator or remnant in beam.
8 // BeamParticle: contains partons, parton densities, etc.
9 
10 #ifndef Pythia8_BeamParticle_H
11 #define Pythia8_BeamParticle_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/Event.h"
15 #include "Pythia8/FragmentationFlavZpT.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/ParticleData.h"
18 #include "Pythia8/PartonDistributions.h"
19 #include "Pythia8/PhysicsBase.h"
20 #include "Pythia8/PythiaStdlib.h"
21 #include "Pythia8/Settings.h"
22 
23 namespace Pythia8 {
24 
25 //==========================================================================
26 
27 // This class holds info on a parton resolved inside the incoming beam,
28 // i.e. either an initiator (part of a hard or a multiparton interaction)
29 // or a remnant (part of the beam remnant treatment).
30 
31 // The companion code is -1 from onset and for g, is -2 for an unmatched
32 // sea quark, is >= 0 for a matched sea quark, with the number giving the
33 // companion position, and is -3 for a valence quark.
34 
35 // Rescattering partons properly do not belong here, but bookkeeping is
36 // simpler with them, so they are stored with companion code -10.
37 
39 
40 public:
41 
42  // Constructor.
43  ResolvedParton( int iPosIn = 0, int idIn = 0, double xIn = 0.,
44  int companionIn = -1) : iPosRes(iPosIn), idRes(idIn), xRes(xIn),
45  companionRes(companionIn), xqCompRes(0.), mRes(0.), factorRes(1.),
46  colRes(0), acolRes(0) { }
47 
48  // Set info on initiator or remnant parton.
49  void iPos( int iPosIn) {iPosRes = iPosIn;}
50  void id( int idIn) {idRes = idIn;}
51  void x( double xIn) {xRes = xIn;}
52  void update( int iPosIn, int idIn, double xIn) {iPosRes = iPosIn;
53  idRes = idIn; xRes = xIn;}
54  void companion( int companionIn) {companionRes = companionIn;}
55  void xqCompanion( double xqCompIn) {xqCompRes = xqCompIn;}
56  void p(Vec4 pIn) {pRes = pIn;}
57  void px(double pxIn) {pRes.px(pxIn);}
58  void py(double pyIn) {pRes.py(pyIn);}
59  void pz(double pzIn) {pRes.pz(pzIn);}
60  void e(double eIn) {pRes.e(eIn);}
61  void m(double mIn) {mRes = mIn;}
62  void col(int colIn) {colRes = colIn;}
63  void acol(int acolIn) {acolRes = acolIn;}
64  void cols(int colIn = 0,int acolIn = 0)
65  {colRes = colIn; acolRes = acolIn;}
66  void scalePT( double factorIn) {pRes.px(factorIn * pRes.px());
67  pRes.py(factorIn * pRes.py()); factorRes *= factorIn;}
68  void scaleX( double factorIn) {xRes *= factorIn;}
69 
70  // Get info on initiator or remnant parton.
71  int iPos() const {return iPosRes;}
72  int id() const {return idRes;}
73  double x() const {return xRes;}
74  int companion() const {return companionRes;}
75  bool isValence() const {return (companionRes == -3);}
76  bool isUnmatched() const {return (companionRes == -2);}
77  bool isCompanion() const {return (companionRes >= 0);}
78  bool isFromBeam() const {return (companionRes > -10);}
79  double xqCompanion() const {return xqCompRes;}
80  Vec4 p() const {return pRes;}
81  double px() const {return pRes.px();}
82  double py() const {return pRes.py();}
83  double pz() const {return pRes.pz();}
84  double e() const {return pRes.e();}
85  double m() const {return mRes;}
86  double pT() const {return pRes.pT();}
87  double mT2() const {return (mRes >= 0.)
88  ? mRes*mRes + pRes.pT2() : - mRes*mRes + pRes.pT2();}
89  double pPos() const {return pRes.e() + pRes.pz();}
90  double pNeg() const {return pRes.e() - pRes.pz();}
91  int col() const {return colRes;}
92  int acol() const {return acolRes;}
93  double pTfactor() const {return factorRes;}
94  bool hasCol() const {return (idRes == 21 || (idRes > 0 && idRes < 9)
95  || (-idRes > 1000 && -idRes < 10000 && (-idRes/10)%10 == 0));}
96  bool hasAcol() const {return (idRes == 21 || (-idRes > 0 && -idRes < 9)
97  || (idRes > 1000 && idRes < 10000 && (idRes/10)%10 == 0));}
98 
99 private:
100 
101  // Properties of a resolved parton.
102  int iPosRes, idRes;
103  double xRes;
104  // Companion code and distribution value, if any.
105  int companionRes;
106  double xqCompRes;
107  // Four-momentum and mass; for remnant kinematics construction.
108  Vec4 pRes;
109  double mRes, factorRes;
110  // Colour codes.
111  int colRes, acolRes;
112 
113 };
114 
115 //==========================================================================
116 
117 // The xfModPrepData struct saves information on the amount of x already used,
118 // to speed up PDF evaluation.
119 
120 typedef struct {
121  double xValTot;
122  double xValLeft;
123  double xLeft;
124  double xCompAdded;
125  double rescaleGS;
126 } xfModPrepData;
127 
128 //==========================================================================
129 
130 // This class holds info on a beam particle in the evolution of
131 // initial-state radiation and multiparton interactions.
132 
133 class BeamParticle : public PhysicsBase {
134 
135 public:
136 
137  // Constructor.
138  BeamParticle() : pdfBeamPtr(), pdfHardBeamPtr(), pdfUnresBeamPtr(),
139  pdfBeamPtrSave(), pdfHardBeamPtrSave(), pdfSavePtrs(),
140  pdfSaveIdx(-1), flavSelPtr(), allowJunction(), beamJunction(),
141  maxValQuark(), companionPower(), valencePowerMeson(), valencePowerUinP(),
142  valencePowerDinP(), valenceDiqEnhance(), pickQuarkNorm(), pickQuarkPower(),
143  diffPrimKTwidth(), diffLargeMassSuppress(), beamSat(), gluonPower(),
144  xGluonCutoff(), heavyQuarkEnhance(), idBeam(), idBeamAbs(), idVMDBeam(),
145  mBeam(), mVMDBeam(), scaleVMDBeam(), isUnresolvedBeam(), isLeptonBeam(),
146  isHadronBeam(), isMesonBeam(), isBaryonBeam(), isGammaBeam(), nValKinds(),
147  idVal(), nVal(), idSave(), iSkipSave(), nValLeft(), xqgTot(), xqVal(),
148  xqgSea(), xqCompSum(), doISR(), doMPI(), doND(), isResolvedGamma(),
149  hasResGammaInBeam(), isResUnres(), hasVMDstateInBeam(), initGammaBeam(),
150  pTminISR(), pTminMPI(), pT2gm2qqbar(), iGamVal(), iPosVal(), gammaMode(),
151  xGm(), Q2gm(), kTgamma(), phiGamma(), cPowerCache(-100), xsCache(-1),
152  resCache(), resolved(), nInit(0), hasJunctionBeam(), junCol(), nJuncs(),
153  nAjuncs(), nDiffJuncs(), allowBeamJunctions(), Q2ValFracSav(-1.),
154  uValInt(), dValInt(), idVal1(), idVal2(), idVal3(), zRel(), pxRel(),
155  pyRel() { }
156 
157  // Initialize data on a beam particle and save pointers.
158  void init( int idIn, double pzIn, double eIn, double mIn,
159  PDFPtr pdfInPtr, PDFPtr pdfHardInPtr, bool isUnresolvedIn,
160  StringFlav* flavSelPtrIn);
161 
162  // Initialize only the id.
163  void initID( int idIn) { idBeam = idIn; initBeamKind();}
164 
165  // Initialize only the two pdf pointers.
166  void initPDFPtr(PDFPtr pdfInPtr, PDFPtr pdfHardInPtr) {
167  pdfBeamPtr = pdfInPtr; pdfHardBeamPtr = pdfHardInPtr; }
168 
169  // Initialize additional PDF pointer for unresolved beam.
170  void initUnres(PDFPtr pdfUnresInPtr);
171 
172  // Initialize array of PDFs for switching between them.
173  void initSwitchID( const vector<PDFPtr>& pdfSavePtrsIn) {
174  pdfSavePtrs = pdfSavePtrsIn; }
175 
176  // For mesons like pi0 valence content varies from event to event.
177  void newValenceContent();
178  void setValenceContent(int idq1, int idq2 = 0, int idq3 = 0);
179 
180  // Switch to new beam particle identity; for similar hadrons only.
181  void setBeamID( int idIn, int iPDFin = -1) {idBeam = idIn;
182  if ( iPDFin >= 0 && iPDFin < int(pdfSavePtrs.size())
183  && iPDFin != pdfSaveIdx ) {
184  pdfBeamPtr = pdfSavePtrs[iPDFin]; pdfHardBeamPtr = pdfBeamPtr;
185  pdfSaveIdx = iPDFin; }
186  mBeam = particleDataPtr->m0(idIn); pdfBeamPtr->setBeamID(idIn); }
187 
188  // Set new pZ and E, but keep the rest the same.
189  void newPzE( double pzIn, double eIn) {pBeam = Vec4( 0., 0., pzIn, eIn);}
190 
191  // Set new mass. Used with photons when virtuality is sampled.
192  void newM( double mIn) { mBeam = mIn; }
193 
194  // Member functions for output.
195  int id() const {return idBeam;}
196  int idVMD() const {return idVMDBeam;}
197  Vec4 p() const {return pBeam;}
198  double px() const {return pBeam.px();}
199  double py() const {return pBeam.py();}
200  double pz() const {return pBeam.pz();}
201  double e() const {return pBeam.e();}
202  double m() const {return mBeam;}
203  double mVMD() const {return mVMDBeam;}
204  double scaleVMD() const {return scaleVMDBeam;}
205  bool isLepton() const {return isLeptonBeam;}
206  bool isUnresolved() const {return isUnresolvedBeam;}
207  // As hadrons here we only count those we know how to handle remnants for.
208  bool isHadron() const {return isHadronBeam;}
209  bool isMeson() const {return isMesonBeam;}
210  bool isBaryon() const {return isBaryonBeam;}
211  bool isGamma() const {return isGammaBeam;}
212  bool hasResGamma() const {return hasResGammaInBeam;}
213  bool hasVMDstate() const {return hasVMDstateInBeam;}
214 
215  // Maximum x remaining after previous MPI and ISR, plus safety margin.
216  double xMax(int iSkip = -1);
217 
218  // Special hard-process parton distributions (can agree with standard ones).
219  double xfHard(int idIn, double x, double Q2)
220  {return pdfHardBeamPtr->xf(idIn, x, Q2);}
221 
222  // Overestimate for PDFs. Same as normal except photons inside leptons.
223  double xfMax(int idIn, double x, double Q2)
224  {return pdfHardBeamPtr->xfMax(idIn, x, Q2);}
225 
226  // Accurate and approximated photon flux and PDFs.
227  double xfFlux(int idIn, double x, double Q2)
228  {return pdfHardBeamPtr->xfFlux(idIn, x, Q2);}
229  double xfApprox(int idIn, double x, double Q2)
230  {return pdfHardBeamPtr->xfApprox(idIn, x, Q2);}
231  double xfGamma(int idIn, double x, double Q2)
232  {return pdfHardBeamPtr->xfGamma(idIn, x, Q2);}
233 
234  // Do not sample the x_gamma value to get correct cross section with
235  // possible second call.
236  double xfSame(int idIn, double x, double Q2)
237  {return pdfHardBeamPtr->xfSame(idIn, x, Q2);}
238 
239  // Standard parton distributions.
240  double xf(int idIn, double x, double Q2)
241  {return pdfBeamPtr->xf(idIn, x, Q2);}
242 
243  // Ditto, split into valence and sea parts (where gluon counts as sea).
244  double xfVal(int idIn, double x, double Q2)
245  {return pdfBeamPtr->xfVal(idIn, x, Q2);}
246  double xfSea(int idIn, double x, double Q2)
247  {return pdfBeamPtr->xfSea(idIn, x, Q2);}
248 
249  // Rescaled parton distributions, as needed for MPI and ISR.
250  // For ISR also allow split valence/sea, and only return relevant part.
251  double xfMPI(int idIn, double x, double Q2, xfModPrepData& xfData)
252  {return xfISR(-1, idIn, x, Q2, xfData);}
253  double xfMPI(int idIn, double x, double Q2)
254  {xfModPrepData xfData = xfModPrep(-1, Q2);
255  return xfISR(-1, idIn, x, Q2, xfData);}
256  double xfISR(int indexMPI, int idIn, double x, double Q2,
257  xfModPrepData& xfData) {return xfModified( indexMPI, idIn, x, Q2, xfData);}
258  double xfISR(int indexMPI, int idIn, double x, double Q2)
259  {xfModPrepData xfData= xfModPrep(indexMPI, Q2);
260  return xfModified( indexMPI, idIn, x, Q2, xfData);}
261 
262  // Check whether x and Q2 values fall inside the fit bounds (LHAPDF6 only).
263  bool insideBounds(double x, double Q2)
264  {return pdfBeamPtr->insideBounds(x,Q2);}
265 
266  // Access the running alpha_s of a PDF set (LHAPDF6 only).
267  double alphaS(double Q2) {return pdfBeamPtr->alphaS(Q2);}
268 
269  // Return quark masses used in the PDF fit (LHAPDF6 only).
270  double mQuarkPDF(int idIn) {return pdfBeamPtr->mQuarkPDF(idIn);}
271 
272  // Return number of members in PDF family (LHAPDF6 only).
273  int nMembers() {return pdfBeamPtr->nMembers();}
274 
275  // Calculate envelope of PDF predictions
276  void calcPDFEnvelope(int idNow, double xNow, double Q2Now, int valSea) {
277  pdfBeamPtr->calcPDFEnvelope(idNow,xNow,Q2Now,valSea);}
278  void calcPDFEnvelope(pair<int,int> idNows, pair<double,double> xNows,
279  double Q2Now, int valSea) {
280  pdfBeamPtr->calcPDFEnvelope(idNows,xNows,Q2Now,valSea);}
281  PDF::PDFEnvelope getPDFEnvelope() { return pdfBeamPtr->getPDFEnvelope(); }
282 
283  // Decide whether chosen quark is valence, sea or companion.
284  int pickValSeaComp();
285 
286  // Initialize kind of incoming beam particle.
287  void initBeamKind();
288 
289  // Overload index operator to access a resolved parton from the list.
290  ResolvedParton& operator[](int i) {return resolved[i];}
291  const ResolvedParton& operator[](int i) const {return resolved[i];}
292  ResolvedParton& at(int i) {return resolved.at(i);}
293  const ResolvedParton& at(int i) const {return resolved.at(i);}
294 
295  // Total number of partons extracted from beam, and initiators only.
296  int size() const {return resolved.size();}
297  int sizeInit() const {return nInit;}
298 
299  // Clear list of resolved partons.
300  void clear() {resolved.resize(0); nInit = 0;}
301 
302  // Reset variables related to photon beam.
303  void resetGamma() {iGamVal = -1; iPosVal = -1; pT2gm2qqbar = 0.;
304  isResolvedGamma = (gammaMode == 1) ? true : false;}
305 
306  // Reset variables related to photon beam inside a lepton.
307  void resetGammaInLepton() {xGm = 1.; kTgamma = 0.; phiGamma = 0.;}
308 
309  // Add a resolved parton to list.
310  int append( int iPos, int idIn, double x, int companion = -1)
311  {resolved.push_back( ResolvedParton( iPos, idIn, x, companion) );
312  return resolved.size() - 1;}
313 
314  // Remove the last particle from the beam. Reset companion code if needed.
315  void popBack() { int iComp = resolved.back().companion();
316  resolved.pop_back(); if ( iComp >= 0 ) { iSkipSave = iComp;
317  idSave = resolved[iComp].id(); pickValSeaComp(); } }
318 
319  // Print extracted parton list; for debug mainly.
320  void list() const;
321 
322  // How many different flavours, and how many quarks of given flavour.
323  int nValenceKinds() const {return nValKinds;}
324  int nValence(int idIn) const {for (int i = 0; i < nValKinds; ++i)
325  if (idIn == idVal[i]) return nVal[i];
326  return 0;}
327 
328  // Test whether a lepton is to be considered as unresolved.
329  bool isUnresolvedLepton();
330 
331  // Add extra remnant flavours to make valence and sea come out right.
332  bool remnantFlavours(Event& event, bool isDIS = false);
333 
334  // Correlate all initiators and remnants to make a colour singlet.
335  bool remnantColours(Event& event, vector<int>& colFrom,
336  vector<int>& colTo);
337 
338  // Pick unrescaled x of remnant parton (valence or sea).
339  double xRemnant(int i);
340 
341  // Tell whether a junction has been resolved, and its junction colours.
342  bool hasJunction() const {return hasJunctionBeam;}
343  int junctionCol(int i) const {return junCol[i];}
344  void junctionCol(int i, int col) {junCol[i] = col;}
345 
346  // For a diffractive system, decide whether to kick out gluon or quark.
347  bool pickGluon(double mDiff);
348 
349  // Pick a valence quark at random, and provide the remaining flavour.
350  int pickValence();
351  int pickRemnant() const {return idVal2;}
352 
353  // Share lightcone momentum between two remnants in a diffractive system.
354  // At the same time generate a relative pT for the two.
355  double zShare( double mDiff, double m1, double m2);
356  double pxShare() const {return pxRel;}
357  double pyShare() const {return pyRel;}
358 
359  // Add extra remnant flavours to make valence and sea come out right.
360  bool remnantFlavoursNew(Event& event);
361 
362  // Find the colour setup of the removed partons from the scatterings.
363  void findColSetup(Event& event);
364 
365  // Set initial colours.
366  void setInitialCol(Event & event);
367 
368  // Update colours.
369  void updateCol(vector<pair<int,int> > colourChanges);
370 
371  vector<pair <int,int> > getColUpdates() {return colUpdates;}
372 
373  // Set valence content for photon beams and position of first valence quark.
374  bool gammaInitiatorIsVal(int iResolved, int id, double x, double Q2);
375  bool gammaInitiatorIsVal(int iResolved, double Q2);
376  int getGammaValFlavour() { return abs(idVal[0]); }
377  int gammaValSeaComp(int iResolved);
378  void posVal(int iPosValIn) { iPosVal = iPosValIn; }
379  void gamVal(int iGamValIn) { iGamVal = iGamValIn; }
380  int gamVal() { return iGamVal; }
381 
382  // Set and get the state (resolved and/or unresolved) of photon beam.
383  void resolvedGamma(bool isResolved) { isResolvedGamma = isResolved; }
384  bool resolvedGamma() const { return isResolvedGamma; }
385  void setGammaMode(int gammaModeIn);
386  int getGammaMode() const { return gammaMode; }
387  bool isResolvedUnresolved() const { return isResUnres; }
388  void initGammaInBeam() { initGammaBeam = true; }
389  bool gammaInBeam() const { return initGammaBeam; }
390 
391  // Set state of VMD inside gamma.
392  void setVMDstate(bool isVMDIn, int idIn, double mIn, double scaleIn,
393  bool reassignState = false) {
394  hasVMDstateInBeam = isVMDIn;
395  idVMDBeam = idIn;
396  mVMDBeam = mIn;
397  scaleVMDBeam = scaleIn;
398  if (reassignState) {
399  idBeam = idVMDBeam;
400  mBeam = mVMDBeam;
401  pdfBeamPtr->setVMDscale(scaleVMDBeam);
402  }
403  }
404 
405  // Store the pT2 value of gamma->qqbar splitting.
406  void pT2gamma2qqbar(double pT2in) { pT2gm2qqbar = pT2in; }
407  double pT2gamma2qqbar() { return pT2gm2qqbar; }
408 
409  // Store the pT value for the latest MPI.
410  void pTMPI(double pTminMPIin) { pTminMPI = pTminMPIin; }
411 
412  // Check whether room for beam remnants.
413  bool roomFor1Remnant(double eCM);
414  bool roomFor1Remnant(int id1, double x1, double eCM);
415  bool roomFor2Remnants(int id1, double x1, double eCM);
416  bool roomForRemnants(BeamParticle beamOther);
417 
418  // Evaluate the remnant mass with initiator idIn.
419  double remnantMass(int idIn);
420 
421  // Functions to approximate pdfs for ISR.
422  double gammaPDFxDependence(int flavour, double x)
423  { return pdfBeamPtr->gammaPDFxDependence(flavour, x); }
424  double gammaPDFRefScale(int flavour)
425  { return pdfBeamPtr->gammaPDFRefScale(flavour); }
426  double xIntegratedPDFs(double Q2)
427  { return pdfBeamPtr->xfIntegratedTotal(Q2); }
428 
429  // Save the x_gamma value after latest PDF call or set it later if ND.
430  void xGammaPDF() { xGm = pdfHardBeamPtr->xGamma(); }
431  void xGamma(double xGmIn) { xGm = xGmIn; }
432  void Q2Gamma(double Q2GmIn) { Q2gm = Q2GmIn; }
433  void newGammaKTPhi(double kTIn, double phiIn)
434  { kTgamma = kTIn; phiGamma = phiIn; }
435 
436  // Get the kinematic limits for photons emitted by the beam.
437  double xGammaMin() { return pdfHardBeamPtr->getXmin(); }
438  double xGammaHadr() { return pdfHardBeamPtr->getXhadr(); }
439  double gammaFluxIntApprox() { return pdfHardBeamPtr->intFluxApprox(); }
440 
441  // Do photon flux use an approximation for sampling.
442  bool hasApproxGammaFlux() { return pdfHardBeamPtr->hasApproxGammaFlux(); }
443 
444  // Get the kinematics related photons form lepton beams.
445  double xGamma() const { return xGm; }
446  double Q2Gamma() const { return Q2gm; }
447  double gammaKTx() const { return kTgamma*cos(phiGamma); }
448  double gammaKTy() const { return kTgamma*sin(phiGamma); }
449  double gammaKT() const { return kTgamma; }
450  double gammaPhi() const { return phiGamma; }
451 
452  // Keep track of pomeron momentum fraction.
453  void xPom(double xpom = -1.0)
454  { if ( pdfBeamPtr ) pdfBeamPtr->xPom(xpom); }
455 
456  // Sample x and Q2 for emitted photons according to flux.
457  double sampleXgamma(double xMinIn)
458  { xGm = pdfHardBeamPtr->sampleXgamma(xMinIn); return xGm; }
459  double sampleQ2gamma(double Q2min)
460  { Q2gm = pdfHardBeamPtr->sampleQ2gamma(Q2min); return Q2gm;}
461 
462  // Prepare data on how much x has been used, for speedup of PDF evaluation.
463  xfModPrepData xfModPrep( int iSkip, double Q2);
464 
465 private:
466 
467  // After initInfoPtr, remove circular dependencies.
468  virtual void onInitInfoPtr() override {
469  beamAPtr = nullptr;
470  beamBPtr = nullptr;
471  beamPomAPtr = nullptr;
472  beamPomBPtr = nullptr;
473  beamGamAPtr = nullptr;
474  beamGamBPtr = nullptr;
475  beamVMDAPtr = nullptr;
476  beamVMDBPtr = nullptr;
477  userHooksPtr = nullptr;
478  }
479 
480  // Constants: could only be changed in the code itself.
481  static const double XMINUNRESOLVED, POMERONMASS, XMAXCOMPANION, TINYZREL;
482  static const int NMAX, NRANDOMTRIES;
483 
484  // Pointers to PDF sets.
485  PDFPtr pdfBeamPtr;
486  PDFPtr pdfHardBeamPtr;
487 
488  // Pointer to unresolved PDF and two others to save the resolved ptrs.
489  PDFPtr pdfUnresBeamPtr;
490  PDFPtr pdfBeamPtrSave;
491  PDFPtr pdfHardBeamPtrSave;
492 
493  // Array of PDFs to be used when idA can be changed between events.
494  vector<PDFPtr> pdfSavePtrs;
495  int pdfSaveIdx;
496 
497  // Pointer to class for flavour generation.
498  StringFlav* flavSelPtr;
499 
500  // Initialization data, normally only set once.
501  bool allowJunction, beamJunction;
502  int maxValQuark, companionPower;
503  double valencePowerMeson, valencePowerUinP, valencePowerDinP,
504  valenceDiqEnhance, pickQuarkNorm, pickQuarkPower,
505  diffPrimKTwidth, diffLargeMassSuppress, beamSat, gluonPower,
506  xGluonCutoff, heavyQuarkEnhance[6];
507 
508  // Basic properties of a beam particle.
509  int idBeam, idBeamAbs, idVMDBeam;
510  Vec4 pBeam;
511  double mBeam, mVMDBeam, scaleVMDBeam;
512 
513  // Beam kind. Valence flavour content for hadrons.
514  bool isUnresolvedBeam, isLeptonBeam, isHadronBeam, isMesonBeam,
515  isBaryonBeam, isGammaBeam;
516  int nValKinds, idVal[3], nVal[3];
517 
518  // Current parton density, by valence, sea and companion.
519  int idSave, iSkipSave, nValLeft[3];
520  double xqgTot, xqVal, xqgSea, xqCompSum;
521 
522  // Variables related to photon beams (also inside lepton).
523  bool doISR, doMPI, doND, isResolvedGamma, hasResGammaInBeam,
524  isResUnres, hasVMDstateInBeam, initGammaBeam;
525  double pTminISR, pTminMPI, pT2gm2qqbar;
526  int iGamVal, iPosVal, gammaMode;
527 
528  // Variables for photon from lepton.
529  double xGm, Q2gm, kTgamma, phiGamma;
530 
531  // Cache values for xCompFrac method.
532  int cPowerCache;
533  double xsCache, resCache;
534 
535  // The list of resolved partons.
536  vector<ResolvedParton> resolved;
537 
538  // Status after all initiators have been accounted for. Junction content.
539  int nInit;
540  bool hasJunctionBeam;
541  int junCol[3];
542 
543  // Variables for new colour reconnection;
544  pair <int,int> colSetup;
545  vector<int> acols, cols;
546  vector<bool> usedCol,usedAcol;
547  vector< pair<int,int> > colUpdates;
548  int nJuncs, nAjuncs, nDiffJuncs;
549  bool allowBeamJunctions;
550 
551  // Routine to calculate PDFs given previous interactions.
552  double xfModified( int iSkip, int idIn, double x, double Q2) {
553  xfModPrepData xfData = xfModPrep(iSkip, Q2);
554  return xfModified(iSkip, idIn, x, Q2, xfData);}
555  double xfModified( int iSkip, int idIn, double x, double Q2,
556  xfModPrepData& data);
557  // In case of resolved.size() === 0.
558  double xfModified0( int iSkip, int idIn, double x, double Q2);
559 
560  // Fraction of hadron momentum sitting in a valence quark distribution.
561  double xValFrac(int j, double Q2);
562  double Q2ValFracSav, uValInt, dValInt;
563 
564  // Fraction of hadron momentum sitting in a companion quark distribution.
565  double xCompFrac(double xs);
566 
567  // Value of companion quark PDF, also given the sea quark x.
568  double xCompDist(double xc, double xs);
569 
570  // Valence quark subdivision for diffractive systems.
571  int idVal1, idVal2, idVal3;
572  double zRel, pxRel, pyRel;
573 
574  // Update a single (anti-) colour of the event.
575  void updateSingleCol(int oldCol, int newCol);
576 
577  // Find a single (anti-) colour in the beam,
578  // if a beam remnant is set the new colour.
579  int findSingleCol(Event& event, bool isAcol, bool useHardScatters);
580 
581 };
582 
583 //==========================================================================
584 
585 } // end namespace Pythia8
586 
587 #endif // Pythia8_BeamParticle_H
ResolvedParton(int iPosIn=0, int idIn=0, double xIn=0., int companionIn=-1)
Constructor.
Definition: BeamParticle.h:43
double xfMax(int idIn, double x, double Q2)
Overestimate for PDFs. Same as normal except photons inside leptons.
Definition: BeamParticle.h:223
void setBeamID(int idIn, int iPDFin=-1)
Switch to new beam particle identity; for similar hadrons only.
Definition: BeamParticle.h:181
Definition: BeamParticle.h:120
void resolvedGamma(bool isResolved)
Set and get the state (resolved and/or unresolved) of photon beam.
Definition: BeamParticle.h:383
double xfFlux(int idIn, double x, double Q2)
Accurate and approximated photon flux and PDFs.
Definition: BeamParticle.h:227
Definition: PhysicsBase.h:27
The Event class holds all info on the generated event.
Definition: Event.h:408
void newM(double mIn)
Set new mass. Used with photons when virtuality is sampled.
Definition: BeamParticle.h:192
bool hasJunction() const
Tell whether a junction has been resolved, and its junction colours.
Definition: BeamParticle.h:342
Definition: BeamParticle.h:133
double xfHard(int idIn, double x, double Q2)
Special hard-process parton distributions (can agree with standard ones).
Definition: BeamParticle.h:219
void clear()
Clear list of resolved partons.
Definition: BeamParticle.h:300
BeamParticle()
Constructor.
Definition: BeamParticle.h:138
void iPos(int iPosIn)
Set info on initiator or remnant parton.
Definition: BeamParticle.h:49
void setVMDstate(bool isVMDIn, int idIn, double mIn, double scaleIn, bool reassignState=false)
Set state of VMD inside gamma.
Definition: BeamParticle.h:392
void newPzE(double pzIn, double eIn)
Set new pZ and E, but keep the rest the same.
Definition: BeamParticle.h:189
void initPDFPtr(PDFPtr pdfInPtr, PDFPtr pdfHardInPtr)
Initialize only the two pdf pointers.
Definition: BeamParticle.h:166
int nValenceKinds() const
How many different flavours, and how many quarks of given flavour.
Definition: BeamParticle.h:323
bool hasApproxGammaFlux()
Do photon flux use an approximation for sampling.
Definition: BeamParticle.h:442
void resetGammaInLepton()
Reset variables related to photon beam inside a lepton.
Definition: BeamParticle.h:307
double xGamma() const
Get the kinematics related photons form lepton beams.
Definition: BeamParticle.h:445
int size() const
Total number of partons extracted from beam, and initiators only.
Definition: BeamParticle.h:296
void xGammaPDF()
Save the x_gamma value after latest PDF call or set it later if ND.
Definition: BeamParticle.h:430
void pTMPI(double pTminMPIin)
Store the pT value for the latest MPI.
Definition: BeamParticle.h:410
Error envelope from PDF uncertainty.
Definition: PartonDistributions.h:102
double xfMPI(int idIn, double x, double Q2, xfModPrepData &xfData)
Definition: BeamParticle.h:251
bool isHadron() const
As hadrons here we only count those we know how to handle remnants for.
Definition: BeamParticle.h:208
void pT2gamma2qqbar(double pT2in)
Store the pT2 value of gamma->qqbar splitting.
Definition: BeamParticle.h:406
double xf(int idIn, double x, double Q2)
Standard parton distributions.
Definition: BeamParticle.h:240
void popBack()
Remove the last particle from the beam. Reset companion code if needed.
Definition: BeamParticle.h:315
void calcPDFEnvelope(int idNow, double xNow, double Q2Now, int valSea)
Calculate envelope of PDF predictions.
Definition: BeamParticle.h:276
ResolvedParton & operator[](int i)
Overload index operator to access a resolved parton from the list.
Definition: BeamParticle.h:290
double mQuarkPDF(int idIn)
Return quark masses used in the PDF fit (LHAPDF6 only).
Definition: BeamParticle.h:270
bool insideBounds(double x, double Q2)
Check whether x and Q2 values fall inside the fit bounds (LHAPDF6 only).
Definition: BeamParticle.h:263
double gammaPDFxDependence(int flavour, double x)
Functions to approximate pdfs for ISR.
Definition: BeamParticle.h:422
void initSwitchID(const vector< PDFPtr > &pdfSavePtrsIn)
Initialize array of PDFs for switching between them.
Definition: BeamParticle.h:173
double alphaS(double Q2)
Access the running alpha_s of a PDF set (LHAPDF6 only).
Definition: BeamParticle.h:267
int iPos() const
Get info on initiator or remnant parton.
Definition: BeamParticle.h:71
int id() const
Member functions for output.
Definition: BeamParticle.h:195
Definition: BeamParticle.h:38
double m2(const Vec4 &v1)
The squared invariant mass of one or more four-vectors.
Definition: Basics.cc:598
int nMembers()
Return number of members in PDF family (LHAPDF6 only).
Definition: BeamParticle.h:273
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
double m(const Vec4 &v1)
Invariant mass and its square.
Definition: Basics.cc:588
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:87
void initID(int idIn)
Initialize only the id.
Definition: BeamParticle.h:163
void resetGamma()
Reset variables related to photon beam.
Definition: BeamParticle.h:303
double xfSame(int idIn, double x, double Q2)
Definition: BeamParticle.h:236
double sampleXgamma(double xMinIn)
Sample x and Q2 for emitted photons according to flux.
Definition: BeamParticle.h:457
double xfVal(int idIn, double x, double Q2)
Ditto, split into valence and sea parts (where gluon counts as sea).
Definition: BeamParticle.h:244
double xGammaMin()
Get the kinematic limits for photons emitted by the beam.
Definition: BeamParticle.h:437
void xPom(double xpom=-1.0)
Keep track of pomeron momentum fraction.
Definition: BeamParticle.h:453
Definition: Basics.h:32
int append(int iPos, int idIn, double x, int companion=-1)
Add a resolved parton to list.
Definition: BeamParticle.h:310