StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
PhaseSpace.h
1 // PhaseSpace.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 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 phase space generators in kinematics selection.
7 // PhaseSpace: base class for phase space generators.
8 // Base class for derived classes> 2 ->1 , 2 -> 2, 2 -> 2 elastic/diffractive,
9 // 2 -> 2 minbias, 2 -> 3, Les Houches.
10 
11 #ifndef Pythia8_PhaseSpace_H
12 #define Pythia8_PhaseSpace_H
13 
14 #include "Basics.h"
15 #include "BeamParticle.h"
16 #include "Info.h"
17 #include "LesHouches.h"
18 #include "MultipartonInteractions.h"
19 #include "ParticleData.h"
20 #include "PartonDistributions.h"
21 #include "PythiaStdlib.h"
22 #include "SigmaProcess.h"
23 #include "SigmaTotal.h"
24 #include "Settings.h"
25 #include "StandardModel.h"
26 #include "UserHooks.h"
27 
28 namespace Pythia8 {
29 
30 //==========================================================================
31 
32 // Forward reference to the UserHooks class.
33 class UserHooks;
34 
35 //==========================================================================
36 
37 // PhaseSpace is a base class for phase space generators
38 // used in the selection of hard-process kinematics.
39 
40 class PhaseSpace {
41 
42 public:
43 
44  // Destructor.
45  virtual ~PhaseSpace() {}
46 
47  // Perform simple initialization and store pointers.
48  void init(bool isFirst, SigmaProcess* sigmaProcessPtrIn,
49  Info* infoPtrIn, Settings* settingsPtrIn, ParticleData* particleDataPtrIn,
50  Rndm* rndmPtrIn, BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
51  Couplings* couplingsPtrIn, SigmaTotal* sigmaTotPtrIn,
52  UserHooks* userHooksPtrIn);
53 
54  // Update the CM energy of the event.
55  void newECM(double eCMin) {eCM = eCMin; s = eCM * eCM;}
56 
57  // Store or replace Les Houches pointer.
58  void setLHAPtr(LHAup* lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;}
59 
60  // A pure virtual method, wherein an optimization procedure
61  // is used to determine how phase space should be sampled.
62  virtual bool setupSampling() = 0;
63 
64  // A pure virtual method, wherein a trial event kinematics
65  // is to be selected in the derived class.
66  virtual bool trialKin(bool inEvent = true, bool repeatSame = false) = 0;
67 
68  // A pure virtual method, wherein the accepted event kinematics
69  // is to be constructed in the derived class.
70  virtual bool finalKin() = 0;
71 
72  // Allow for nonisotropic decays when ME's available.
73  void decayKinematics( Event& process);
74 
75  // Give back current or maximum cross section, or set latter.
76  double sigmaNow() const {return sigmaNw;}
77  double sigmaMax() const {return sigmaMx;}
78  double biasSelectionWeight() const {return biasWt;}
79  bool newSigmaMax() const {return newSigmaMx;}
80  void setSigmaMax(double sigmaMaxIn) {sigmaMx = sigmaMaxIn;}
81 
82  // For Les Houches with negative event weight needs
83  virtual double sigmaSumSigned() const {return sigmaMx;}
84 
85  // Give back constructed four-vectors and known masses.
86  Vec4 p(int i) const {return pH[i];}
87  double m(int i) const {return mH[i];}
88 
89  // Give back other event properties.
90  double ecm() const {return eCM;}
91  double x1() const {return x1H;}
92  double x2() const {return x2H;}
93  double sHat() const {return sH;}
94  double tHat() const {return tH;}
95  double uHat() const {return uH;}
96  double pTHat() const {return pTH;}
97  double thetaHat() const {return theta;}
98  double phiHat() const {return phi;}
99  double runBW3() const {return runBW3H;}
100  double runBW4() const {return runBW4H;}
101  double runBW5() const {return runBW5H;}
102 
103  // Inform whether beam particles are resolved in partons or scatter directly.
104  virtual bool isResolved() const {return true;}
105 
106 protected:
107 
108  // Constructor.
109  PhaseSpace() {}
110 
111  // Constants: could only be changed in the code itself.
112  static const int NMAXTRY, NTRY3BODY;
113  static const double SAFETYMARGIN, TINY, EVENFRAC, SAMESIGMA, WIDTHMARGIN,
114  SAMEMASS, MASSMARGIN, EXTRABWWTMAX, THRESHOLDSIZE,
115  THRESHOLDSTEP, YRANGEMARGIN, LEPTONXMIN, LEPTONXMAX,
116  LEPTONXLOGMIN, LEPTONXLOGMAX, LEPTONTAUMIN,
117  SHATMINZ, PT2RATMINZ, WTCORRECTION[11];
118 
119  // Pointer to cross section.
120  SigmaProcess* sigmaProcessPtr;
121 
122  // Pointer to various information on the generation.
123  Info* infoPtr;
124 
125  // Pointer to the settings database.
126  Settings* settingsPtr;
127 
128  // Pointer to the particle data table.
129  ParticleData* particleDataPtr;
130 
131  // Pointer to the random number generator.
132  Rndm* rndmPtr;
133 
134  // Pointers to incoming beams.
135  BeamParticle* beamAPtr;
136  BeamParticle* beamBPtr;
137 
138  // Pointer to Standard Model couplings.
139  Couplings* couplingsPtr;
140 
141  // Pointer to the total/elastic/diffractive cross section object.
142  SigmaTotal* sigmaTotPtr;
143 
144  // Pointer to userHooks object for user interaction with program.
145  UserHooks* userHooksPtr;
146 
147  // Pointer to LHAup for generating external events.
148  LHAup* lhaUpPtr;
149 
150  // Initialization data, normally only set once.
151  bool useBreitWigners, doEnergySpread, showSearch, showViolation,
152  increaseMaximum;
153  int gmZmodeGlobal;
154  double mHatGlobalMin, mHatGlobalMax, pTHatGlobalMin, pTHatGlobalMax,
155  pTHatMinDiverge, minWidthBreitWigners;
156 
157  // Information on incoming beams.
158  int idA, idB;
159  double mA, mB, eCM, s;
160  bool hasLeptonBeams, hasPointLeptons;
161 
162  // Cross section information.
163  bool newSigmaMx, canModifySigma, canBiasSelection;
164  int gmZmode;
165  double wtBW, sigmaNw, sigmaMx, sigmaPos, sigmaNeg, biasWt;
166 
167  // Process-specific kinematics properties, almost always available.
168  double mHatMin, mHatMax, sHatMin, sHatMax, pTHatMin, pTHatMax,
169  pT2HatMin, pT2HatMax;
170 
171  // Event-specific kinematics properties, almost always available.
172  double x1H, x2H, m3, m4, m5, s3, s4, s5, mHat, sH, tH, uH, pAbs, p2Abs,
173  pTH, theta, phi, betaZ;
174  Vec4 pH[6];
175  double mH[6];
176 
177  // Reselect decay products momenta isotropically in phase space.
178  void decayKinematicsStep( Event& process, int iRes);
179 
180  // Much common code for normal 2 -> 1, 2 -> 2 and 2 -> 3 cases:
181 
182  // Determine how phase space should be sampled.
183  void setup3Body();
184  bool setupSampling123(bool is2, bool is3, ostream& os = cout);
185 
186  // Select a trial kinematics phase space point.
187  bool trialKin123(bool is2, bool is3, bool inEvent = true,
188  ostream& os = cout);
189 
190  // Presence and properties of any s-channel resonances.
191  int idResA, idResB;
192  double mResA, mResB, GammaResA, GammaResB, tauResA, tauResB, widResA,
193  widResB;
194  bool sameResMass;
195 
196  // Kinematics properties specific to 2 -> 1/2/3.
197  bool useMirrorWeight;
198  double tau, y, z, tauMin, tauMax, yMax, zMin, zMax, ratio34, unity34,
199  zNeg, zPos, wtTau, wtY, wtZ, wt3Body, runBW3H, runBW4H, runBW5H,
200  intTau0, intTau1, intTau2, intTau3, intTau4, intTau5, intTau6,
201  intY0, intY12, intY34, intY56, mTchan1, sTchan1, mTchan2, sTchan2,
202  frac3Flat, frac3Pow1, frac3Pow2;
203  Vec4 p3cm, p4cm, p5cm;
204 
205  // Coefficients for optimized selection in 2 -> 1/2/3.
206  int nTau, nY, nZ;
207  double tauCoef[8], yCoef[8], zCoef[8], tauCoefSum[8], yCoefSum[8],
208  zCoefSum[8];
209 
210  // Calculate kinematical limits for 2 -> 1/2/3.
211  bool limitTau(bool is2, bool is3);
212  bool limitY();
213  bool limitZ();
214 
215  // Select kinematical variable between defined limits for 2 -> 1/2/3.
216  void selectTau(int iTau, double tauVal, bool is2);
217  void selectY(int iY, double yVal);
218  void selectZ(int iZ, double zVal);
219  bool select3Body();
220 
221  // Solve equation system for better phase space coefficients in 2 -> 1/2/3.
222  void solveSys( int n, int bin[8], double vec[8], double mat[8][8],
223  double coef[8], ostream& os = cout);
224 
225  // Properties specific to resonance mass selection in 2 -> 2 and 2 -> 3.
226  bool useBW[6];
227  int idMass[6];
228  double mPeak[6], sPeak[6], mWidth[6], mMin[6], mMax[6], mw[6], wmRat[6],
229  mLower[6], mUpper[6], sLower[6], sUpper[6], fracFlat[6], fracInv[6],
230  fracInv2[6], atanLower[6], atanUpper[6], intBW[6], intFlat[6],
231  intInv[6], intInv2[6];
232 
233  // Setup mass selection for one resonance at a time. Split in two parts.
234  void setupMass1(int iM);
235  void setupMass2(int iM, double distToThresh);
236 
237  // Do mass selection and find the associated weight.
238  void trialMass(int iM);
239  double weightMass(int iM);
240 
241 };
242 
243 //==========================================================================
244 
245 // A derived class with 2 -> 1 kinematics set up in tau, y.
246 
248 
249 public:
250 
251  // Constructor.
252  PhaseSpace2to1tauy() {}
253 
254  // Optimize subsequent kinematics selection.
255  virtual bool setupSampling() {if (!setupMass()) return false;
256  return setupSampling123(false, false);}
257 
258  // Construct the trial kinematics.
259  virtual bool trialKin(bool inEvent = true, bool = false) {wtBW = 1.;
260  return trialKin123(false, false, inEvent);}
261 
262  // Construct the final event kinematics.
263  virtual bool finalKin();
264 
265 private:
266 
267  // Set up allowed mass range.
268  bool setupMass();
269 
270 };
271 
272 //==========================================================================
273 
274 // A derived class with 2 -> 2 kinematics set up in tau, y, z = cos(theta).
275 
277 
278 public:
279 
280  // Constructor.
282 
283  // Optimize subsequent kinematics selection.
284  virtual bool setupSampling() {if (!setupMasses()) return false;
285  return setupSampling123(true, false);}
286 
287  // Construct the trial kinematics.
288  virtual bool trialKin(bool inEvent = true, bool = false) {
289  if (!trialMasses()) return false;
290  return trialKin123(true, false, inEvent);}
291 
292  // Construct the final event kinematics.
293  virtual bool finalKin();
294 
295 private:
296 
297  // Set up for fixed or Breit-Wigner mass selection.
298  bool setupMasses();
299 
300  // Select fixed or Breit-Wigner-distributed masses.
301  bool trialMasses();
302 
303  // Pick off-shell initialization masses when on-shell not allowed.
304  bool constrainedM3M4();
305  bool constrainedM3();
306  bool constrainedM4();
307 
308 };
309 
310 //==========================================================================
311 
312 // A derived class with 2 -> 2 kinematics set up for elastic scattering.
313 
315 
316 public:
317 
318  // Constructor.
320 
321  // Construct the trial or final event kinematics.
322  virtual bool setupSampling();
323  virtual bool trialKin(bool inEvent = true, bool = false);
324  virtual bool finalKin();
325 
326  // Are beam particles resolved in partons or scatter directly?
327  virtual bool isResolved() const {return false;}
328 
329 private:
330 
331  // Constants: could only be changed in the code itself.
332  static const double EXPMAX, CONVERTEL;
333 
334  // Kinematics properties specific to 2 -> 2 elastic.
335  bool useCoulomb;
336  double s1, s2, bSlope, lambda12S, tLow, tUpp, tAux, sigmaTot, rho,
337  lambda, tAbsMin, phaseCst, alphaEM0, sigmaNuc, sigmaCou, signCou;
338 
339  // Calculation of alphaElectromagnetic.
340  AlphaEM alphaEM;
341 
342 };
343 
344 //==========================================================================
345 
346 // A derived class with 2 -> 2 kinematics set up for diffractive scattering.
347 
349 
350 public:
351 
352  // Constructor.
353  PhaseSpace2to2diffractive(bool isDiffAin = false, bool isDiffBin = false)
354  : isDiffA(isDiffAin), isDiffB(isDiffBin) {}
355 
356  // Construct the trial or final event kinematics.
357  virtual bool setupSampling();
358  virtual bool trialKin(bool inEvent = true, bool = false);
359  virtual bool finalKin();
360 
361  // Are beam particles resolved in partons or scatter directly?
362  virtual bool isResolved() const {return false;}
363 
364 private:
365 
366  // Constants: could only be changed in the code itself.
367  static const int NTRY;
368  static const double EXPMAX, DIFFMASSMAX;
369 
370  // Initialization data, in constructor or read from Settings.
371  bool isDiffA, isDiffB;
372  int PomFlux;
373  double epsilonPF, alphaPrimePF;
374 
375  // Initialization: kinematics properties specific to 2 -> 2 diffractive.
376  double m3ElDiff, m4ElDiff, s1, s2, lambda12, lambda34, tLow, tUpp,
377  cRes, sResXB, sResAX, sProton, bMin, bSlope, bSlope1, bSlope2,
378  probSlope1, xIntPF, xtCorPF, mp24DL, coefDL, tAux, tAux1, tAux2;
379 
380 };
381 
382 //==========================================================================
383 
384 // A derived class for minumum bias events. Hardly does anything, since
385 // the real action is taken care of by the MultipartonInteractions class.
386 
388 
389 public:
390 
391  // Constructor.
393 
394  // Construct the trial or final event kinematics.
395  virtual bool setupSampling() {sigmaNw = sigmaProcessPtr->sigmaHat();
396  sigmaMx = sigmaNw; return true;}
397  virtual bool trialKin( bool , bool = false) {return true;}
398  virtual bool finalKin() {return true;}
399 
400 private:
401 
402 };
403 
404 //==========================================================================
405 
406 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
407 // tau, y, pT2_4, pT2_5, phi_4, phi_5 and y_3 (partial cylindrical symmetry).
408 
410 
411 public:
412 
413  // Constructor.
415 
416  // Optimize subsequent kinematics selection.
417  virtual bool setupSampling() {if (!setupMasses()) return false;
418  setup3Body(); return setupSampling123(false, true);}
419 
420  // Construct the trial kinematics.
421  virtual bool trialKin(bool inEvent = true, bool = false) {
422  if (!trialMasses()) return false;
423  return trialKin123(false, true, inEvent);}
424 
425  // Construct the final event kinematics.
426  virtual bool finalKin();
427 
428 private:
429 
430  // Constants: could only be changed in the code itself.
431  static const int NITERNR;
432 
433  // Set up for fixed or Breit-Wigner mass selection.
434  bool setupMasses();
435 
436  // Select fixed or Breit-Wigner-distributed masses.
437  bool trialMasses();
438 
439 };
440 
441 //==========================================================================
442 
443 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
444 // y3, y4, y5, pT2_3, pT2_5, phi_3 and phi_5, and with R separation cut.
445 // Intended specifically for (essentially massless) 2 -> 3 QCD processes.
446 
448 
449 public:
450 
451  // Constructor.
453 
454  // Optimize subsequent kinematics selection.
455  virtual bool setupSampling();
456 
457  // Construct the trial kinematics.
458  virtual bool trialKin(bool inEvent = true, bool = false);
459 
460  // Construct the final event kinematics.
461  virtual bool finalKin();
462 
463 private:
464 
465  // Phase space cuts specifically for 2 -> 3 QCD processes.
466  double pTHat3Min, pTHat3Max, pTHat5Min, pTHat5Max, RsepMin, R2sepMin;
467  bool hasBaryonBeams;
468 
469  // Event kinematics choices.
470  double pT3Min, pT3Max, pT5Min, pT5Max, y3Max, y4Max, y5Max,
471  pT3, pT4, pT5, phi3, phi4, phi5, y3, y4, y5, dphi;
472  Vec4 pInSum;
473 
474 };
475 
476 //==========================================================================
477 
478 // A derived class for Les Houches events.
479 
480 class PhaseSpaceLHA : public PhaseSpace {
481 
482 public:
483 
484  // Constructor.
485  PhaseSpaceLHA() {idProcSave = 0;}
486 
487  // Find maximal cross section for comparison with internal processes.
488  virtual bool setupSampling();
489 
490  // Construct the next process, by interface to Les Houches class.
491  virtual bool trialKin( bool , bool repeatSame = false);
492 
493  // Set scale, alpha_s and alpha_em if not done.
494  virtual bool finalKin() {sigmaProcessPtr->setScale(); return true;}
495 
496  // For Les Houches with negative event weight needs
497  virtual double sigmaSumSigned() const {return sigmaSgn;}
498 
499 private:
500 
501  // Constants.
502  static const double CONVERTPB2MB;
503 
504  // Local properties.
505  int strategy, stratAbs, nProc, idProcSave;
506  double xMaxAbsSum, xSecSgnSum, sigmaSgn;
507  vector<int> idProc;
508  vector<double> xMaxAbsProc;
509 
510 };
511 
512 //==========================================================================
513 
514 } // end namespace Pythia8
515 
516 #endif // Pythia8_PhaseSpace_H
517