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) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, 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/PythiaStdlib.h"
20 #include "Pythia8/Settings.h"
21 #include "Pythia8/StandardModel.h"
22 #include "Pythia8/UserHooks.h"
23 #include "Pythia8/MergingHooks.h"
24 #include "Pythia8/WeakShowerMEs.h"
25 
26 namespace Pythia8 {
27 
28 //==========================================================================
29 
30 // Data on radiating dipole ends, only used inside SpaceShower.
31 
32 class SpaceDipoleEnd {
33 
34 public:
35 
36  // Constructor.
37  SpaceDipoleEnd( int systemIn = 0, int sideIn = 0, int iRadiatorIn = 0,
38  int iRecoilerIn = 0, double pTmaxIn = 0., int colTypeIn = 0,
39  int chgTypeIn = 0, int weakTypeIn = 0, int MEtypeIn = 0,
40  bool normalRecoilIn = true, int weakPolIn = 0) :
41  system(systemIn), side(sideIn), iRadiator(iRadiatorIn),
42  iRecoiler(iRecoilerIn), pTmax(pTmaxIn), colType(colTypeIn),
43  chgType(chgTypeIn), weakType(weakTypeIn), MEtype(MEtypeIn),
44  normalRecoil(normalRecoilIn), weakPol(weakPolIn), nBranch(0),
45  pT2Old(0.), zOld(0.5) { }
46 
47  // Store values for trial emission.
48  void store( int idDaughterIn, int idMotherIn, int idSisterIn,
49  double x1In, double x2In, double m2DipIn, double pT2In, double zIn,
50  double xMoIn, double Q2In, double mSisterIn, double m2SisterIn,
51  double pT2corrIn) {idDaughter = idDaughterIn; idMother = idMotherIn;
52  idSister = idSisterIn; x1 = x1In; x2 = x2In; m2Dip = m2DipIn;
53  pT2 = pT2In; z = zIn; xMo = xMoIn; Q2 = Q2In; mSister = mSisterIn;
54  m2Sister = m2SisterIn; pT2corr = pT2corrIn;}
55 
56  // Basic properties related to evolution and matrix element corrections.
57  int system, side, iRadiator, iRecoiler;
58  double pTmax;
59  int colType, chgType, weakType, MEtype;
60  bool normalRecoil;
61  int weakPol;
62 
63  // Properties specific to current trial emission.
64  int nBranch, idDaughter, idMother, idSister, iFinPol;
65  double x1, x2, m2Dip, pT2, z, xMo, Q2, mSister, m2Sister, pT2corr,
66  pT2Old, zOld, asymPol;
67 
68 } ;
69 
70 //==========================================================================
71 
72 // The SpaceShower class does spacelike showers.
73 
74 class SpaceShower {
75 
76 public:
77 
78  // Constructor.
79  SpaceShower() {beamOffset = 0;}
80 
81  // Destructor.
82  virtual ~SpaceShower() {}
83 
84  // Initialize various pointers.
85  // (Separated from rest of init since not virtual.)
86  void initPtr(Info* infoPtrIn, Settings* settingsPtrIn,
87  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn, CoupSM* coupSMPtrIn,
88  PartonSystems* partonSystemsPtrIn, UserHooks* userHooksPtrIn,
89  MergingHooks* mergingHooksPtrIn = 0) {
90  infoPtr = infoPtrIn; settingsPtr = settingsPtrIn;
91  particleDataPtr = particleDataPtrIn; rndmPtr = rndmPtrIn;
92  coupSMPtr = coupSMPtrIn; partonSystemsPtr = partonSystemsPtrIn;
93  userHooksPtr = userHooksPtrIn; mergingHooksPtr = mergingHooksPtrIn;}
94 
95  // Initialize generation. Possibility to force re-initialization by hand.
96  virtual void init(BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn);
97 
98  // New beams possible for handling of hard diffraction. (Not virtual.)
99  void reassignBeamPtrs( BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
100  int beamOffsetIn = 0) {beamAPtr = beamAPtrIn; beamBPtr = beamBPtrIn;
101  beamOffset = beamOffsetIn;}
102 
103  // Find whether to limit maximum scale of emissions, and whether to dampen.
104  virtual bool limitPTmax( Event& event, double Q2Fac = 0.,
105  double Q2Ren = 0.);
106 
107  // Potential enhancement factor of pTmax scale for hardest emission.
108  virtual double enhancePTmax() const {return pTmaxFudge;}
109 
110  // Prepare system for evolution; identify ME.
111  virtual void prepare( int iSys, Event& event, bool limitPTmaxIn = true);
112 
113  // Update dipole list after each FSR emission.
114  // Usage: update( iSys, event).
115  virtual void update( int , Event&, bool hasWeakRad = false);
116 
117  // Select next pT in downwards evolution.
118  virtual double pTnext( Event& event, double pTbegAll, double pTendAll,
119  int nRadIn = -1);
120 
121  // ME corrections and kinematics that may give failure.
122  virtual bool branch( Event& event);
123 
124  // Tell which system was the last processed one.
125  int system() const {return iSysSel;}
126 
127  // Flag for failure in branch(...) that will force a retry of parton level.
128  bool doRestart() const {return rescatterFail;}
129 
130  // Tell whether ISR has done a weak emission.
131  bool getHasWeaklyRadiated() {return hasWeaklyRadiated;}
132 
133  // Print dipole list; for debug mainly.
134  virtual void list(ostream& os = cout) const;
135 
136 protected:
137 
138  // Pointer to various information on the generation.
139  Info* infoPtr;
140 
141  // Pointer to the settings database.
142  Settings* settingsPtr;
143 
144  // Pointer to the particle data table.
145  ParticleData* particleDataPtr;
146 
147  // Pointer to the random number generator.
148  Rndm* rndmPtr;
149 
150  // Pointer to Standard Model couplings.
151  CoupSM* coupSMPtr;
152 
153  // Pointers to the two incoming beams. Offset their location in event.
154  BeamParticle* beamAPtr;
155  BeamParticle* beamBPtr;
156  int beamOffset;
157 
158  // Pointer to information on subcollision parton locations.
159  PartonSystems* partonSystemsPtr;
160 
161  // Pointer to userHooks object for user interaction with program.
162  UserHooks* userHooksPtr;
163 
164  // Weak matrix elements used for corrections both of ISR and FSR.
165  WeakShowerMEs weakShowerMEs;
166 
167  // Store properties to be returned by methods.
168  bool rescatterFail;
169  int iSysSel;
170  double pTmaxFudge;
171 
172 private:
173 
174  // Constants: could only be changed in the code itself.
175  static const int MAXLOOPTINYPDF;
176  static const double MCMIN, MBMIN, CTHRESHOLD, BTHRESHOLD, EVALPDFSTEP,
177  TINYPDF, TINYKERNELPDF, TINYPT2, HEAVYPT2EVOL, HEAVYXEVOL,
178  EXTRASPACEQ, LAMBDA3MARGIN, PT2MINWARN, LEPTONXMIN, LEPTONXMAX,
179  LEPTONPT2MIN, LEPTONFUDGE, WEAKPSWEIGHT, HEADROOMQ2Q, HEADROOMQ2G,
180  HEADROOMG2G, HEADROOMG2Q;
181 
182  // Initialization data, normally only set once.
183  bool doQCDshower, doQEDshowerByQ, doQEDshowerByL, useSamePTasMPI,
184  doWeakShower, doMEcorrections, doMEafterFirst, doPhiPolAsym,
185  doPhiIntAsym, doRapidityOrder, useFixedFacScale, doSecondHard,
186  canVetoEmission, alphaSuseCMW, singleWeakEmission, vetoWeakJets;
187  int pTmaxMatch, pTdampMatch, alphaSorder, alphaSnfmax, alphaEMorder,
188  nQuarkIn, enhanceScreening, weakMode;
189  double pTdampFudge, mc, mb, m2c, m2b, renormMultFac, factorMultFac,
190  fixedFacScale2, alphaSvalue, alphaS2pi, Lambda3flav, Lambda4flav,
191  Lambda5flav, Lambda3flav2, Lambda4flav2, Lambda5flav2, pT0Ref,
192  ecmRef, ecmPow, pTmin, sCM, eCM, pT0, pTminChgQ, pTminChgL, pT20,
193  pT2min, pT2minChgQ, pT2minChgL, pTweakCut, pT2weakCut, pTmaxFudgeMPI,
194  strengthIntAsym, weakEnhancement, mZ, gammaZ, thetaWRat, mW, gammaW,
195  weakMaxWt, vetoWeakDeltaR2;
196 
197  // alphaStrong and alphaEM calculations.
198  AlphaStrong alphaS;
199  AlphaEM alphaEM;
200 
201  // Some current values.
202  bool sideA, dopTlimit1, dopTlimit2, dopTdamp, hasWeaklyRadiated, tChannel;
203  int iNow, iRec, idDaughter, nRad, idResFirst, idResSecond;
204  double xDaughter, x1Now, x2Now, m2Dip, m2Rec, pT2damp, pTbegRef, pdfScale2;
205 
206  // List of emissions in different sides in different systems:
207  vector<int> nRadA,nRadB;
208 
209  // All dipole ends
210  vector<SpaceDipoleEnd> dipEnd;
211 
212  // Pointers to the current and hardest (so far) dipole ends.
213  int iDipNow, iSysNow;
214  SpaceDipoleEnd* dipEndNow;
215  int iDipSel;
216  SpaceDipoleEnd* dipEndSel;
217 
218  // Evolve a QCD dipole end.
219  void pT2nextQCD( double pT2begDip, double pT2endDip);
220 
221  // Evolve a QCD dipole end near heavy quark threshold region.
222  void pT2nearQCDthreshold( BeamParticle& beam, double m2Massive,
223  double m2Threshold, double xMaxAbs, double zMinAbs,
224  double zMaxMassive);
225 
226  // Evolve a QED dipole end.
227  void pT2nextQED( double pT2begDip, double pT2endDip);
228 
229  // Evolve a Weak dipole end.
230  void pT2nextWeak( double pT2begDip, double pT2endDip);
231 
232  // Find class of ME correction.
233  int findMEtype( int iSys, Event& event, bool weakRadiation = false);
234 
235  // Provide maximum of expected ME weight; for preweighting of evolution.
236  double calcMEmax( int MEtype, int idMother, int idDaughterIn);
237 
238  // Provide actual ME weight for current branching.
239  double calcMEcorr(int MEtype, int idMother, int idDaughterIn, double M2,
240  double z, double Q2,double m2Sister);
241 
242  // Provide actual ME weight for t-channel weak emissions.
243  double calcMEcorrWeak(int MEtype, double m2, double z,
244  double pT2, Vec4 pMother, Vec4 pB, Vec4 pDaughter,
245  Vec4 pB0, Vec4 p1, Vec4 p2, Vec4 pSister);
246 
247  // Find coefficient of azimuthal asymmetry from gluon polarization.
248  void findAsymPol( Event& event, SpaceDipoleEnd* dip);
249 
250  // Pointer to MergingHooks object for NLO merging.
251  MergingHooks* mergingHooksPtr;
252 
253 };
254 
255 //==========================================================================
256 
257 } // end namespace Pythia8
258 
259 #endif // Pythia8_SpaceShower_H
Definition: beam.h:43
Definition: AgUStep.h:26