StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SimpleSpaceShower.h
1 // SimpleSpaceShower.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 the original simple spacelike initial-state showers.
7 // SpaceDipoleEnd: data on a radiating dipole end in ISR.
8 // SimpleSpaceShower: handles the showering description.
9 
10 #ifndef Pythia8_SimpleSpaceShower_H
11 #define Pythia8_SimpleSpaceShower_H
12 
13 #include "Pythia8/SpaceShower.h"
14 #include "Pythia8/SimpleWeakShowerMEs.h"
15 
16 namespace Pythia8 {
17 
18 //==========================================================================
19 
20 // Data on radiating dipole ends, only used inside SimpleSpaceShower.
21 
22 class SpaceDipoleEnd {
23 
24 public:
25 
26  // Constructor.
27  SpaceDipoleEnd( int systemIn = 0, int sideIn = 0, int iRadiatorIn = 0,
28  int iRecoilerIn = 0, double pTmaxIn = 0., int colTypeIn = 0,
29  int chgTypeIn = 0, int weakTypeIn = 0, int MEtypeIn = 0,
30  bool normalRecoilIn = true, int weakPolIn = 0,
31  int iColPartnerIn = 0, int idColPartnerIn = 0) :
32  system(systemIn), side(sideIn), iRadiator(iRadiatorIn),
33  iRecoiler(iRecoilerIn), pTmax(pTmaxIn), colType(colTypeIn),
34  chgType(chgTypeIn), weakType(weakTypeIn), MEtype(MEtypeIn),
35  normalRecoil(normalRecoilIn), weakPol(weakPolIn),
36  iColPartner(iColPartnerIn), idColPartner(idColPartnerIn),
37  nBranch(0), idDaughter(), idMother(), idSister(), iFinPol(), x1(), x2(),
38  m2Dip(), pT2(), z(), xMo(), Q2(), mSister(), m2Sister(), pT2corr(),
39  pT2Old(0.), zOld(0.5), asymPol(), m2IF(), mColPartner(),
40  pAccept() { }
41 
42  // Store values for trial emission.
43  void store( int idDaughterIn, int idMotherIn, int idSisterIn,
44  double x1In, double x2In, double m2DipIn, double pT2In, double zIn,
45  double xMoIn, double Q2In, double mSisterIn, double m2SisterIn,
46  double pT2corrIn, int iColPartnerIn, double m2IFIn, double mColPartnerIn)
47  {idDaughter = idDaughterIn; idMother = idMotherIn;
48  idSister = idSisterIn; x1 = x1In; x2 = x2In; m2Dip = m2DipIn;
49  pT2 = pT2In; z = zIn; xMo = xMoIn; Q2 = Q2In; mSister = mSisterIn;
50  m2Sister = m2SisterIn; pT2corr = pT2corrIn; iColPartner = iColPartnerIn;
51  m2IF = m2IFIn; mColPartner = mColPartnerIn;}
52 
53  // Basic properties related to evolution and matrix element corrections.
54  int system, side, iRadiator, iRecoiler;
55  double pTmax;
56  int colType, chgType, weakType, MEtype;
57  bool normalRecoil;
58  int weakPol, iColPartner, idColPartner;
59 
60  // Properties specific to current trial emission.
61  int nBranch, idDaughter, idMother, idSister, iFinPol;
62  double x1, x2, m2Dip, pT2, z, xMo, Q2, mSister, m2Sister, pT2corr,
63  pT2Old, zOld, asymPol, m2IF, mColPartner;
64 
65  // Properties needed for the evaluation of parameter variations
66  double pAccept;
67 
68 } ;
69 
70 //==========================================================================
71 
72 // The SimpleSpaceShower class does spacelike showers.
73 
75 
76 public:
77 
78  // Constructor.
79  SimpleSpaceShower() : rescatterFail(), gamma2qqbar(), hasWeaklyRadiated(),
80  iSysSel(), pTmaxFudge(), doQCDshower(), doQEDshowerByQ(), doQEDshowerByL(),
81  useSamePTasMPI(), doWeakShower(), doMEcorrections(), doMEafterFirst(),
82  doPhiPolAsym(), doPhiPolAsymHard(), doPhiIntAsym(), doRapidityOrder(),
83  useFixedFacScale(), doSecondHard(), canVetoEmission(), hasUserHooks(),
84  alphaSuseCMW(), singleWeakEmission(), vetoWeakJets(), weakExternal(),
85  doRapidityOrderMPI(), doMPI(), doDipoleRecoil(), doPartonVertex(),
86  pTmaxMatch(), pTdampMatch(), alphaSorder(), alphaSnfmax(), alphaEMorder(),
87  nQuarkIn(), enhanceScreening(), weakMode(), pT0paramMode(), pTdampFudge(),
88  mc(), mb(), m2c(), m2b(), renormMultFac(), factorMultFac(),
89  fixedFacScale2(), alphaSvalue(), alphaS2pi(), Lambda3flav(), Lambda4flav(),
90  Lambda5flav(), Lambda3flav2(), Lambda4flav2(), Lambda5flav2(), pT0Ref(),
91  ecmRef(), ecmPow(), pTmin(), sCM(), eCM(), pT0(), pTminChgQ(), pTminChgL(),
92  pT20(), pT2min(), pT2minChgQ(), pT2minChgL(), pTweakCut(), pT2weakCut(),
93  pTmaxFudgeMPI(), strengthIntAsym(), weakEnhancement(), mZ(), gammaZ(),
94  thetaWRat(), mW(), gammaW(), weakMaxWt(), vetoWeakDeltaR2(), sideA(),
95  twoHard(), dopTlimit1(), dopTlimit2(), dopTdamp(), tChannel(),
96  doUncertaintiesNow(), iNow(), iRec(), idDaughter(), nRad(), idResFirst(),
97  idResSecond(), xDaughter(), x1Now(), x2Now(), m2ColPair(), mColPartner(),
98  m2ColPartner(), m2Dip(), m2Rec(), pT2damp(), pTbegRef(), pdfScale2(),
99  doTrialNow(), canEnhanceEmission(), canEnhanceTrial(), canEnhanceET(),
100  iDipNow(), iSysNow(), dipEndNow(), iDipSel(), dipEndSel() {
101  beamOffset = 0;}
102 
103  // Destructor.
104  virtual ~SimpleSpaceShower() {}
105 
106  // Initialize generation. Possibility to force re-initialization by hand.
107  virtual void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn);
108 
109  // Find whether to limit maximum scale of emissions, and whether to dampen.
110  virtual bool limitPTmax( Event& event, double Q2Fac = 0.,
111  double Q2Ren = 0.);
112 
113  // Prepare system for evolution; identify ME.
114  virtual void prepare( int iSys, Event& event, bool limitPTmaxIn = true);
115 
116  // Update dipole list after each FSR emission.
117  virtual void update( int iSys, Event& event, bool hasWeakRad = false);
118 
119  // Select next pT in downwards evolution.
120  virtual double pTnext( Event& event, double pTbegAll, double pTendAll,
121  int nRadIn = -1, bool doTrialIn = false);
122 
123  // ME corrections and kinematics that may give failure.
124  virtual bool branch( Event& event);
125 
126  // Print dipole list; for debug mainly.
127  virtual void list() const;
128 
129  // Initialize data members for calculation of uncertainty bands.
130  virtual bool initUncertainties();
131 
132  // Flag for failure in branch(...) that will force a retry of parton level.
133  virtual bool doRestart() const {return rescatterFail;}
134 
135  // Tell if latest scattering was a gamma->qqbar.
136  virtual bool wasGamma2qqbar() { return gamma2qqbar; }
137 
138  // Tell whether ISR has done a weak emission.
139  virtual bool getHasWeaklyRadiated() {return hasWeaklyRadiated;}
140 
141  // Tell which system was the last processed one.
142  virtual int system() const {return iSysSel;}
143 
144  // Potential enhancement factor of pTmax scale for hardest emission.
145  virtual double enhancePTmax() const {return pTmaxFudge;}
146 
147 private:
148 
149  // Constants: could only be changed in the code itself.
150  static const int MAXLOOPTINYPDF;
151  static const double MCMIN, MBMIN, CTHRESHOLD, BTHRESHOLD, EVALPDFSTEP,
152  TINYPDF, TINYKERNELPDF, TINYPT2, HEAVYPT2EVOL, HEAVYXEVOL,
153  EXTRASPACEQ, LAMBDA3MARGIN, PT2MINWARN, LEPTONXMIN, LEPTONXMAX,
154  LEPTONPT2MIN, LEPTONFUDGE, WEAKPSWEIGHT, HEADROOMQ2Q, HEADROOMQ2G,
155  HEADROOMG2G, HEADROOMG2Q, HEADROOMHQG, REJECTFACTOR, PROBLIMIT;
156 
157  // Store properties to be returned by methods.
158  bool rescatterFail, gamma2qqbar, hasWeaklyRadiated;
159  int iSysSel;
160  double pTmaxFudge;
161 
162  // Initialization data, normally only set once.
163  bool doQCDshower, doQEDshowerByQ, doQEDshowerByL, useSamePTasMPI,
164  doWeakShower, doMEcorrections, doMEafterFirst, doPhiPolAsym,
165  doPhiPolAsymHard, doPhiIntAsym, doRapidityOrder, useFixedFacScale,
166  doSecondHard, canVetoEmission, hasUserHooks, alphaSuseCMW,
167  singleWeakEmission, vetoWeakJets, weakExternal, doRapidityOrderMPI,
168  doMPI, doDipoleRecoil, doPartonVertex;
169  int pTmaxMatch, pTdampMatch, alphaSorder, alphaSnfmax, alphaEMorder,
170  nQuarkIn, enhanceScreening, weakMode, pT0paramMode;
171  double pTdampFudge, mc, mb, m2c, m2b, renormMultFac, factorMultFac,
172  fixedFacScale2, alphaSvalue, alphaS2pi, Lambda3flav, Lambda4flav,
173  Lambda5flav, Lambda3flav2, Lambda4flav2, Lambda5flav2, pT0Ref,
174  ecmRef, ecmPow, pTmin, sCM, eCM, pT0, pTminChgQ, pTminChgL, pT20,
175  pT2min, pT2minChgQ, pT2minChgL, pTweakCut, pT2weakCut, pTmaxFudgeMPI,
176  strengthIntAsym, weakEnhancement, mZ, gammaZ, thetaWRat, mW, gammaW,
177  weakMaxWt, vetoWeakDeltaR2;
178 
179  // alphaStrong and alphaEM calculations.
180  AlphaStrong alphaS;
181  AlphaEM alphaEM;
182 
183  // Weak matrix elements used for corrections both of ISR and FSR.
184  SimpleWeakShowerMEs simpleWeakShowerMEs;
185 
186  // Some current values.
187  bool sideA, twoHard, dopTlimit1, dopTlimit2, dopTdamp, tChannel,
188  doUncertaintiesNow;
189  int iNow, iRec, idDaughter, nRad, idResFirst, idResSecond;
190  double xDaughter, x1Now, x2Now, m2ColPair, mColPartner, m2ColPartner,
191  m2Dip, m2Rec, pT2damp, pTbegRef, pdfScale2;
192 
193  // Bookkeeping of enhanced actual or trial emissions (see EPJC (2013) 73).
194  bool doTrialNow, canEnhanceEmission, canEnhanceTrial, canEnhanceET;
195  string splittingNameNow, splittingNameSel;
196  map< double, pair<string,double> > enhanceFactors;
197  void storeEnhanceFactor(double pT2, string name, double enhanceFactorIn)
198  { enhanceFactors.insert(make_pair(pT2,make_pair(name,enhanceFactorIn)));}
199 
200  // List of emissions in different sides in different systems:
201  vector<int> nRadA,nRadB;
202 
203  // All dipole ends
204  vector<SpaceDipoleEnd> dipEnd;
205 
206  // List of 2 -> 2 momenta for external weak setup.
207  vector<Vec4> weakMomenta;
208 
209  // Pointers to the current and hardest (so far) dipole ends.
210  int iDipNow, iSysNow;
211  SpaceDipoleEnd* dipEndNow;
212  int iDipSel;
213  SpaceDipoleEnd* dipEndSel;
214 
215  // Evolve a QCD dipole end.
216  void pT2nextQCD( double pT2begDip, double pT2endDip);
217 
218  // Evolve a QCD and QED dipole end near heavy quark threshold region.
219  void pT2nearThreshold( BeamParticle& beam, double m2Massive,
220  double m2Threshold, double xMaxAbs, double zMinAbs,
221  double zMaxMassive, int iColPartner);
222 
223  // Evolve a QED dipole end.
224  void pT2nextQED( double pT2begDip, double pT2endDip);
225 
226  // Evolve a Weak dipole end.
227  void pT2nextWeak( double pT2begDip, double pT2endDip);
228 
229  // Find class of ME correction.
230  int findMEtype( int iSys, Event& event, bool weakRadiation = false);
231 
232  // Provide maximum of expected ME weight; for preweighting of evolution.
233  double calcMEmax( int MEtype, int idMother, int idDaughterIn);
234 
235  // Provide actual ME weight for current branching.
236  double calcMEcorr(int MEtype, int idMother, int idDaughterIn, double M2,
237  double z, double Q2,double m2Sister);
238 
239  // Provide actual ME weight for t-channel weak emissions.
240  double calcMEcorrWeak(int MEtype, double m2, double z,
241  double pT2, Vec4 pMother, Vec4 pB, Vec4 pDaughter,
242  Vec4 pB0, Vec4 p1, Vec4 p2, Vec4 pSister);
243 
244  // Find coefficient of azimuthal asymmetry from gluon polarization.
245  void findAsymPol( Event& event, SpaceDipoleEnd* dip);
246 
247  // Find a possible colour partner in the case of dipole recoil.
248  int findColPartner(Event& event, int iSideA, int iSideB, int iSystem);
249 
250  // Calculate uncertainty-band weights for accepted/rejected trial branching.
251  void calcUncertainties(bool accept, double pAcceptIn, double pT20in,
252  double enhance, double vp, SpaceDipoleEnd* dip, Particle* motherPtr,
253  Particle* sisterPtr);
254 
255 };
256 
257 //==========================================================================
258 
259 } // end namespace Pythia8
260 
261 #endif // Pythia8_SimpleSpaceShower_H
Definition: beam.h:43