StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
SpaceShower.h
1 // SpaceShower.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2018 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 spacelike initial-state showers.
7 // SpaceDipoleEnd: radiating dipole end in ISR.
8 // SpaceShower: handles the showering description.
9 
10 #ifndef Pythia8_SpaceShower_H
11 #define Pythia8_SpaceShower_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/BeamParticle.h"
15 #include "Pythia8/Event.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/ParticleData.h"
18 #include "Pythia8/PartonSystems.h"
19 #include "Pythia8/PartonVertex.h"
20 #include "Pythia8/PythiaStdlib.h"
21 #include "Pythia8/Settings.h"
22 #include "Pythia8/StandardModel.h"
23 #include "Pythia8/UserHooks.h"
24 #include "Pythia8/MergingHooks.h"
25 #include "Pythia8/WeakShowerMEs.h"
26 
27 
28 namespace Pythia8 {
29 
30 //==========================================================================
31 
32 // Data on radiating dipole ends, only used inside SpaceShower.
33 
34 class SpaceDipoleEnd {
35 
36 public:
37 
38  // Constructor.
39  SpaceDipoleEnd( int systemIn = 0, int sideIn = 0, int iRadiatorIn = 0,
40  int iRecoilerIn = 0, double pTmaxIn = 0., int colTypeIn = 0,
41  int chgTypeIn = 0, int weakTypeIn = 0, int MEtypeIn = 0,
42  bool normalRecoilIn = true, int weakPolIn = 0,
43  int iColPartnerIn = 0, int idColPartnerIn = 0) :
44  system(systemIn), side(sideIn), iRadiator(iRadiatorIn),
45  iRecoiler(iRecoilerIn), pTmax(pTmaxIn), colType(colTypeIn),
46  chgType(chgTypeIn), weakType(weakTypeIn), MEtype(MEtypeIn),
47  normalRecoil(normalRecoilIn), weakPol(weakPolIn),
48  iColPartner(iColPartnerIn), idColPartner(idColPartnerIn),
49  nBranch(0), pT2Old(0.), zOld(0.5) { }
50 
51  // Store values for trial emission.
52  void store( int idDaughterIn, int idMotherIn, int idSisterIn,
53  double x1In, double x2In, double m2DipIn, double pT2In, double zIn,
54  double xMoIn, double Q2In, double mSisterIn, double m2SisterIn,
55  double pT2corrIn, int iColPartnerIn, double m2IFIn, double mColPartnerIn)
56  {idDaughter = idDaughterIn; idMother = idMotherIn;
57  idSister = idSisterIn; x1 = x1In; x2 = x2In; m2Dip = m2DipIn;
58  pT2 = pT2In; z = zIn; xMo = xMoIn; Q2 = Q2In; mSister = mSisterIn;
59  m2Sister = m2SisterIn; pT2corr = pT2corrIn; iColPartner = iColPartnerIn;
60  m2IF = m2IFIn; mColPartner = mColPartnerIn;}
61 
62  // Basic properties related to evolution and matrix element corrections.
63  int system, side, iRadiator, iRecoiler;
64  double pTmax;
65  int colType, chgType, weakType, MEtype;
66  bool normalRecoil;
67  int weakPol, iColPartner, idColPartner;
68 
69  // Properties specific to current trial emission.
70  int nBranch, idDaughter, idMother, idSister, iFinPol;
71  double x1, x2, m2Dip, pT2, z, xMo, Q2, mSister, m2Sister, pT2corr,
72  pT2Old, zOld, asymPol, m2IF, mColPartner;
73 
74  // Properties needed for the evaluation of parameter variations
75  double pAccept;
76 
77 } ;
78 
79 //==========================================================================
80 
81 // The SpaceShower class does spacelike showers.
82 
83 class SpaceShower {
84 
85 public:
86 
87  // Constructor.
88  SpaceShower() {beamOffset = 0;}
89 
90  // Destructor.
91  virtual ~SpaceShower() {}
92 
93  // Initialize various pointers.
94  // (Separated from rest of init since not virtual.)
95  void initPtr(Info* infoPtrIn, Settings* settingsPtrIn,
96  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn, CoupSM* coupSMPtrIn,
97  PartonSystems* partonSystemsPtrIn, UserHooks* userHooksPtrIn,
98  MergingHooks* mergingHooksPtrIn, PartonVertex* partonVertexPtrIn) {
99  infoPtr = infoPtrIn; settingsPtr = settingsPtrIn;
100  particleDataPtr = particleDataPtrIn; rndmPtr = rndmPtrIn;
101  coupSMPtr = coupSMPtrIn; partonSystemsPtr = partonSystemsPtrIn;
102  userHooksPtr = userHooksPtrIn; mergingHooksPtr = mergingHooksPtrIn;
103  partonVertexPtr = partonVertexPtrIn; }
104 
105  // Initialize generation. Possibility to force re-initialization by hand.
106  virtual void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn);
107 
108  // New beams possible for handling of hard diffraction. (Not virtual.)
109  void reassignBeamPtrs( BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
110  int beamOffsetIn = 0) {beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;
111  beamOffset = beamOffsetIn;}
112 
113  // Find whether to limit maximum scale of emissions, and whether to dampen.
114  virtual bool limitPTmax( Event& event, double Q2Fac = 0.,
115  double Q2Ren = 0.);
116 
117  // Potential enhancement factor of pTmax scale for hardest emission.
118  virtual double enhancePTmax() const {return pTmaxFudge;}
119 
120  // Prepare system for evolution; identify ME.
121  virtual void prepare( int iSys, Event& event, bool limitPTmaxIn = true);
122 
123  // Update dipole list after each FSR emission.
124  // Usage: update( iSys, event).
125  virtual void update( int , Event&, bool hasWeakRad = false);
126 
127  // Select next pT in downwards evolution.
128  virtual double pTnext( Event& event, double pTbegAll, double pTendAll,
129  int nRadIn = -1, bool doTrialIn = false);
130 
131  // ME corrections and kinematics that may give failure.
132  virtual bool branch( Event& event);
133 
134  // Initialize data members for calculation of uncertainty bands.
135  bool initUncertainties();
136 
137  // Calculate uncertainty-band weights for accepted/rejected trial branching.
138  void calcUncertainties(bool accept, double pAcceptIn, double pT20in,
139  double enhance, double vp, SpaceDipoleEnd* dip, Particle* motherPtr,
140  Particle* sisterPtr);
141 
142  // Tell if latest scattering was a gamma->qqbar.
143  bool wasGamma2qqbar() { return gamma2qqbar; }
144 
145  // Tell which system was the last processed one.
146  virtual int system() const {return iSysSel;}
147 
148  // Flag for failure in branch(...) that will force a retry of parton level.
149  bool doRestart() const {return rescatterFail;}
150 
151  // Tell whether ISR has done a weak emission.
152  bool getHasWeaklyRadiated() {return hasWeaklyRadiated;}
153 
154  // Print dipole list; for debug mainly.
155  virtual void list() const;
156 
157  // Functions to allow usage of shower kinematics, evolution variables,
158  // and splitting probabilities outside of shower.
159  // Virtual so that shower plugins can overwrite these functions.
160  // This makes it possible for another piece of the code to request
161  // these - which is very convenient for merging.
162  // Function variable names are not included to avoid compiler warnings.
163  // Please see the documentation under "Implement New Showers" for details.
164 
165  // Return clustering kinematics - as needed form merging.
166  virtual Event clustered( const Event& , int , int , int , string )
167  { return Event();}
168 
169  // Return the evolution variable(s).
170  // Important note: this map must contain the following entries
171  // - a key "t" for the value of the shower evolution variable;
172  // - a key "tRS" for the value of the shower evolution variable
173  // from which the shower would be restarted after a branching;
174  // - a key "scaleAS" for the argument of alpha_s used for the branching;
175  // - a key "scalePDF" for the argument of the PDFs used for the branching.
176  // Usage: getStateVariables( event, iRad, iEmt, iRec, name)
177  virtual map<string, double> getStateVariables (const Event& , int , int ,
178  int , string ) { return map<string,double>();}
179 
180  // Check if attempted clustering is handled by spacelike shower.
181  // Usage: isSpacelike( event, iRad, iEmt, iRec, name)
182  virtual bool isSpacelike(const Event&, int, int, int, string)
183  { return false; }
184 
185  // Return a string identifier of a splitting.
186  // Usage: getSplittingName( event, iRad, iEmt, iRec)
187  virtual vector<string> getSplittingName( const Event& , int , int , int )
188  { return vector<string>();}
189 
190  // Return the splitting probability.
191  // Usage: getSplittingProb( event, iRad, iEmt, iRec)
192  virtual double getSplittingProb( const Event& , int , int , int , string )
193  { return 0.;}
194 
195  virtual bool allowedSplitting( const Event& , int , int)
196  { return true;}
197  virtual vector<int> getRecoilers( const Event&, int, int, string)
198  { return vector<int>(); }
199 
200  // Pointer to MergingHooks object for NLO merging.
201  MergingHooks* mergingHooksPtr;
202 
203 protected:
204 
205  // Pointer to various information on the generation.
206  Info* infoPtr;
207 
208  // Pointer to the settings database.
209  Settings* settingsPtr;
210 
211  // Pointer to the particle data table.
212  ParticleData* particleDataPtr;
213 
214  // Pointer to the random number generator.
215  Rndm* rndmPtr;
216 
217  // Pointer to Standard Model couplings.
218  CoupSM* coupSMPtr;
219 
220  // Pointers to the two incoming beams. Offset their location in event.
221  BeamParticle* beamAPtr;
222  BeamParticle* beamBPtr;
223  int beamOffset;
224 
225  // Pointer to information on subcollision parton locations.
226  PartonSystems* partonSystemsPtr;
227 
228  // Pointer to userHooks object for user interaction with program.
229  UserHooks* userHooksPtr;
230 
231  // Pointer to assign space-time vertices during parton evolution.
232  PartonVertex* partonVertexPtr;
233 
234  // Weak matrix elements used for corrections both of ISR and FSR.
235  WeakShowerMEs weakShowerMEs;
236 
237  // Store properties to be returned by methods.
238  bool rescatterFail;
239  int iSysSel;
240  double pTmaxFudge;
241 
242 private:
243 
244  // Constants: could only be changed in the code itself.
245  static const int MAXLOOPTINYPDF;
246  static const double MCMIN, MBMIN, CTHRESHOLD, BTHRESHOLD, EVALPDFSTEP,
247  TINYPDF, TINYKERNELPDF, TINYPT2, HEAVYPT2EVOL, HEAVYXEVOL,
248  EXTRASPACEQ, LAMBDA3MARGIN, PT2MINWARN, LEPTONXMIN, LEPTONXMAX,
249  LEPTONPT2MIN, LEPTONFUDGE, WEAKPSWEIGHT, HEADROOMQ2Q, HEADROOMQ2G,
250  HEADROOMG2G, HEADROOMG2Q, HEADROOMHQG, REJECTFACTOR, PROBLIMIT;
251 
252  // Initialization data, normally only set once.
253  bool doQCDshower, doQEDshowerByQ, doQEDshowerByL, useSamePTasMPI,
254  doWeakShower, doMEcorrections, doMEafterFirst, doPhiPolAsym,
255  doPhiPolAsymHard, doPhiIntAsym, doRapidityOrder, useFixedFacScale,
256  doSecondHard, canVetoEmission, hasUserHooks, alphaSuseCMW,
257  singleWeakEmission, vetoWeakJets, weakExternal, doRapidityOrderMPI,
258  doUncertainties, uVarMuSoftCorr, uVarMPIshowers, doMPI, gamma2qqbar,
259  doDipoleRecoil, doPartonVertex;
260  int pTmaxMatch, pTdampMatch, alphaSorder, alphaSnfmax, alphaEMorder,
261  nQuarkIn, enhanceScreening, weakMode, pT0paramMode;
262  double pTdampFudge, mc, mb, m2c, m2b, renormMultFac, factorMultFac,
263  fixedFacScale2, alphaSvalue, alphaS2pi, Lambda3flav, Lambda4flav,
264  Lambda5flav, Lambda3flav2, Lambda4flav2, Lambda5flav2, pT0Ref,
265  ecmRef, ecmPow, pTmin, sCM, eCM, pT0, pTminChgQ, pTminChgL, pT20,
266  pT2min, pT2minChgQ, pT2minChgL, pTweakCut, pT2weakCut, pTmaxFudgeMPI,
267  strengthIntAsym, weakEnhancement, mZ, gammaZ, thetaWRat, mW, gammaW,
268  weakMaxWt, vetoWeakDeltaR2, dASmax, cNSpTmin, uVarpTmin2, overFactor;
269 
270  // alphaStrong and alphaEM calculations.
271  AlphaStrong alphaS;
272  AlphaEM alphaEM;
273 
274  // Some current values.
275  bool sideA, dopTlimit1, dopTlimit2, dopTdamp, hasWeaklyRadiated, tChannel,
276  doUncertaintiesNow;
277  int iNow, iRec, idDaughter, nRad, idResFirst, idResSecond;
278  double xDaughter, x1Now, x2Now, m2ColPair, mColPartner, m2ColPartner,
279  m2Dip, m2Rec, pT2damp, pTbegRef, pdfScale2;
280 
281  // Bookkeeping of enhanced actual or trial emissions (see EPJC (2013) 73).
282  bool doTrialNow, canEnhanceEmission, canEnhanceTrial, canEnhanceET;
283  string splittingNameNow, splittingNameSel;
284  map< double, pair<string,double> > enhanceFactors;
285  void storeEnhanceFactor(double pT2, string name, double enhanceFactorIn)
286  { enhanceFactors.insert(make_pair(pT2,make_pair(name,enhanceFactorIn)));}
287 
288  // List of emissions in different sides in different systems:
289  vector<int> nRadA,nRadB;
290 
291  // All dipole ends
292  vector<SpaceDipoleEnd> dipEnd;
293 
294  // List of 2 -> 2 momenta for external weak setup.
295  vector<Vec4> weakMomenta;
296 
297  // Pointers to the current and hardest (so far) dipole ends.
298  int iDipNow, iSysNow;
299  SpaceDipoleEnd* dipEndNow;
300  int iDipSel;
301  SpaceDipoleEnd* dipEndSel;
302 
303  // Evolve a QCD dipole end.
304  void pT2nextQCD( double pT2begDip, double pT2endDip);
305 
306  // Evolve a QCD and QED dipole end near heavy quark threshold region.
307  void pT2nearThreshold( BeamParticle& beam, double m2Massive,
308  double m2Threshold, double xMaxAbs, double zMinAbs,
309  double zMaxMassive, int iColPartner);
310 
311  // Evolve a QED dipole end.
312  void pT2nextQED( double pT2begDip, double pT2endDip);
313 
314  // Evolve a Weak dipole end.
315  void pT2nextWeak( double pT2begDip, double pT2endDip);
316 
317  // Find class of ME correction.
318  int findMEtype( int iSys, Event& event, bool weakRadiation = false);
319 
320  // Provide maximum of expected ME weight; for preweighting of evolution.
321  double calcMEmax( int MEtype, int idMother, int idDaughterIn);
322 
323  // Provide actual ME weight for current branching.
324  double calcMEcorr(int MEtype, int idMother, int idDaughterIn, double M2,
325  double z, double Q2,double m2Sister);
326 
327  // Provide actual ME weight for t-channel weak emissions.
328  double calcMEcorrWeak(int MEtype, double m2, double z,
329  double pT2, Vec4 pMother, Vec4 pB, Vec4 pDaughter,
330  Vec4 pB0, Vec4 p1, Vec4 p2, Vec4 pSister);
331 
332  // Find coefficient of azimuthal asymmetry from gluon polarization.
333  void findAsymPol( Event& event, SpaceDipoleEnd* dip);
334 
335  // Store uncertainty variations relevant to TimeShower.
336  int nUncertaintyVariations, nVarQCD, uVarNflavQ;
337  map<int,double> varG2GGmuRfac, varQ2QGmuRfac, varQ2GQmuRfac, varG2QQmuRfac,
338  varX2XGmuRfac;
339  map<int,double> varG2GGcNS, varQ2QGcNS, varQ2GQcNS, varG2QQcNS, varX2XGcNS;
340  map<int,double>* varPDFplus;
341  map<int,double>* varPDFminus;
342  map<int,double>* varPDFmember;
343 
344  // Find a possible colour partner in the case of dipole recoil.
345  int findColPartner(Event& event, int iSideA, int iSideB, int iSystem);
346 
347 };
348 
349 //==========================================================================
350 
351 } // end namespace Pythia8
352 
353 #endif // Pythia8_SpaceShower_H
Definition: beam.h:43
Definition: AgUStep.h:26