Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-03-29 08:28:06

0001 // PhaseSpace.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2025 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 
0006 // Header file for phase space generators in kinematics selection.
0007 // PhaseSpace: base class for phase space generators.
0008 // Base class for derived classes> 2 ->1 , 2 -> 2, 2 -> 2 elastic/diffractive,
0009 // 2 -> 2 nondiffractive, 2 -> 3, Les Houches.
0010 
0011 #ifndef Pythia8_PhaseSpace_H
0012 #define Pythia8_PhaseSpace_H
0013 
0014 #include "Pythia8/Basics.h"
0015 #include "Pythia8/BeamParticle.h"
0016 #include "Pythia8/Info.h"
0017 #include "Pythia8/LesHouches.h"
0018 #include "Pythia8/MultipartonInteractions.h"
0019 #include "Pythia8/ParticleData.h"
0020 #include "Pythia8/PartonDistributions.h"
0021 #include "Pythia8/PhysicsBase.h"
0022 #include "Pythia8/PythiaStdlib.h"
0023 #include "Pythia8/SigmaProcess.h"
0024 #include "Pythia8/SigmaTotal.h"
0025 #include "Pythia8/Settings.h"
0026 #include "Pythia8/StandardModel.h"
0027 #include "Pythia8/UserHooks.h"
0028 #include "Pythia8/GammaKinematics.h"
0029 
0030 namespace Pythia8 {
0031 
0032 //==========================================================================
0033 
0034 // Forward reference to the UserHooks class.
0035 class UserHooks;
0036 
0037 //==========================================================================
0038 
0039 // PhaseSpace is a base class for  phase space generators
0040 // used in the selection of hard-process kinematics.
0041 
0042 class PhaseSpace : public PhysicsBase {
0043 
0044 public:
0045 
0046   // Destructor.
0047   virtual ~PhaseSpace() {}
0048 
0049   // Perform simple initialization and store pointers.
0050   void init(bool isFirst, SigmaProcessPtr sigmaProcessPtrIn);
0051 
0052   // Switch to new beam particle identities; for similar hadrons only.
0053   void updateBeamIDs() { idAold = idA; idBold = idB; idA = beamAPtr->id();
0054     idB = beamBPtr->id(); mA = beamAPtr->m(); mB = beamBPtr->m();
0055     sigmaProcessPtr->updateBeamIDs();}
0056 
0057   // Update the CM energy of the event.
0058   void newECM(double eCMin) {eCM = eCMin; s = eCM * eCM;}
0059 
0060   // Store or replace Les Houches pointer.
0061   void setLHAPtr(LHAupPtr lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;}
0062 
0063   // A pure virtual method, wherein an optimization procedure
0064   // is used to determine how phase space should be sampled.
0065   virtual bool setupSampling() = 0;
0066 
0067   // A pure virtual method, wherein a trial event kinematics
0068   // is to be selected in the derived class.
0069   virtual bool trialKin(bool inEvent = true, bool repeatSame = false) = 0;
0070 
0071   // A pure virtual method, wherein the accepted event kinematics
0072   // is to be constructed in the derived class.
0073   virtual bool finalKin() = 0;
0074 
0075   // Allow for nonisotropic decays when ME's available.
0076   void   decayKinematics( Event& process);
0077 
0078   // Give back current or maximum cross section, or set latter.
0079   double sigmaNow() const {return sigmaNw;}
0080   double sigmaMax() const {return sigmaMx;}
0081   double biasSelectionWeight()  const {return biasWt;}
0082   bool   newSigmaMax() const {return newSigmaMx;}
0083   void   setSigmaMax(double sigmaMaxIn) {sigmaMx = sigmaMaxIn;}
0084 
0085   // New value for switched beam identity or energy (for SoftQCD processes).
0086   double sigmaMaxSwitch() {sigmaNw = sigmaProcessPtr->sigmaHatWrap();
0087     sigmaMx = sigmaNw; return sigmaMx;}
0088 
0089   // For Les Houches with negative event weight needs
0090   virtual double sigmaSumSigned() const {return sigmaMx;}
0091 
0092   // Give back constructed four-vectors and known masses.
0093   Vec4   p(int i)   const {return pH[i];}
0094   double m(int i)   const {return mH[i];}
0095 
0096   // Reset the four-momentum.
0097   void   setP(int i, Vec4 pNew) { pH[i] = pNew; }
0098 
0099   // Give back other event properties.
0100   double ecm()      const {return eCM;}
0101   double x1()       const {return x1H;}
0102   double x2()       const {return x2H;}
0103   double sHat()     const {return sH;}
0104   double tHat()     const {return tH;}
0105   double uHat()     const {return uH;}
0106   double pTHat()    const {return pTH;}
0107   double thetaHat() const {return theta;}
0108   double phiHat()   const {return phi;}
0109   double runBW3()   const {return runBW3H;}
0110   double runBW4()   const {return runBW4H;}
0111   double runBW5()   const {return runBW5H;}
0112 
0113   // Inform whether beam particles are resolved in partons or scatter directly.
0114   virtual bool isResolved() const {return true;}
0115 
0116   // Functions to rescale momenta and cross section for new sHat
0117   // Currently implemented only for PhaseSpace2to2tauyz class.
0118   virtual void rescaleSigma( double){}
0119   virtual void rescaleMomenta( double){}
0120 
0121   // Calculate the weight for over-estimated cross section.
0122   virtual double weightGammaPDFApprox(){ return 1.;}
0123 
0124   // Set the GammaKinematics pointer needed for soft photoproduction.
0125   virtual void setGammaKinPtr( GammaKinematics* gammaKinPtrIn) {
0126     gammaKinPtr = gammaKinPtrIn; }
0127 
0128 protected:
0129 
0130   // Constructor.
0131   PhaseSpace() : sigmaProcessPtr(), lhaUpPtr(), gammaKinPtr(),
0132     useBreitWigners(), doEnergySpread(), showSearch(),
0133     showViolation(), increaseMaximum(), hasQ2Min(), gmZmodeGlobal(),
0134     mHatGlobalMin(), mHatGlobalMax(), pTHatGlobalMin(),
0135     pTHatGlobalMax(), Q2GlobalMin(), Q2GlobalMax(), pTHatMinDiverge(),
0136     minWidthBreitWigners(), minWidthNarrowBW(), idA(), idB(),
0137     idAold(), idBold(), idAgm(), idBgm(), mA(), mB(), eCM(), s(),
0138     sigmaMxGm(), hasLeptonBeamA(), hasLeptonBeamB(),
0139     hasOneLeptonBeam(), hasTwoLeptonBeams(), hasPointGammaA(),
0140     hasPointGammaB(), hasOnePointParticle(), hasTwoPointParticles(),
0141     hasGamma(), hasVMD(), newSigmaMx(), canModifySigma(),
0142     canBiasSelection(), canBias2Sel(), gmZmode(), bias2SelPow(),
0143     bias2SelRef(), wtBW(), sigmaNw(), sigmaMx(), sigmaPos(),
0144     sigmaNeg(), biasWt(), mHatMin(), mHatMax(), sHatMin(), sHatMax(),
0145     pTHatMin(), pTHatMax(), pT2HatMin(), pT2HatMax(), x1H(), x2H(),
0146     m3(), m4(), m5(), s3(), s4(), s5(), mHat(), sH(), tH(), uH(),
0147     pAbs(), p2Abs(), pTH(), theta(), phi(), betaZ(), mH(), idResA(),
0148     idResB(), mResA(), mResB(), GammaResA(), GammaResB(), tauResA(),
0149     tauResB(), widResA(), widResB(), sameResMass(), useMirrorWeight(),
0150     hasNegZ(), hasPosZ(), tau(), y(), z(), tauMin(), tauMax(), yMax(),
0151     zMin(), zMax(), ratio34(), unity34(), zNeg(), zPos(), wtTau(),
0152     wtY(), wtZ(), wt3Body(), runBW3H(), runBW4H(), runBW5H(),
0153     intTau0(), intTau1(), intTau2(), intTau3(), intTau4(), intTau5(),
0154     intTau6(), intY0(), intY12(), intY34(), intY56(), mTchan1(),
0155     sTchan1(), mTchan2(), sTchan2(), frac3Flat(), frac3Pow1(),
0156     frac3Pow2(), zNegMin(), zNegMax(), zPosMin(), zPosMax(), nTau(),
0157     nY(), nZ(), tauCoef(), yCoef(), zCoef(), tauCoefSum(), yCoefSum(),
0158     zCoefSum(), useBW(), useNarrowBW(), idMass(), mPeak(), sPeak(),
0159     mWidth(), mMin(), mMax(), mw(), wmRat(), mLower(), mUpper(),
0160     sLower(), sUpper(), fracFlatS(), fracFlatM(), fracInv(),
0161     fracInv2(), atanLower(), atanUpper(), intBW(), intFlatS(),
0162     intFlatM(), intInv(), intInv2(), doTopPair(), topThresholdModel(),
0163     topThresholdWidth(), eThreshold(), m3Threshold(), m4Threshold() {}
0164 
0165   // Constants: could only be changed in the code itself.
0166   static const int    NMAXTRY, NTRY3BODY;
0167   static const double SAFETYMARGIN, TINY, EVENFRAC, SAMESIGMA, MRESMINABS,
0168                       WIDTHMARGIN, SAMEMASS, MASSMARGIN, EXTRABWWTMAX,
0169                       THRESHOLDSIZE, THRESHOLDSTEP, YRANGEMARGIN, LEPTONXMIN,
0170                       LEPTONXMAX, LEPTONXLOGMIN, LEPTONXLOGMAX, LEPTONTAUMIN,
0171                       SHATMINZ, PT2RATMINZ, WTCORRECTION[11];
0172 
0173   // Pointer to cross section.
0174   SigmaProcessPtr sigmaProcessPtr;
0175 
0176   // Pointer to LHAup for generating external events.
0177   LHAupPtr      lhaUpPtr;
0178 
0179   // Pointer to object that samples photon kinematics from leptons.
0180   GammaKinematics* gammaKinPtr;
0181 
0182   // Initialization data, normally only set once.
0183   bool   useBreitWigners, doEnergySpread, showSearch, showViolation,
0184          increaseMaximum, hasQ2Min;
0185   int    gmZmodeGlobal;
0186   double mHatGlobalMin, mHatGlobalMax, pTHatGlobalMin, pTHatGlobalMax,
0187          Q2GlobalMin, Q2GlobalMax, pTHatMinDiverge, minWidthBreitWigners,
0188          minWidthNarrowBW;
0189 
0190   // Information on incoming beams.
0191   int    idA, idB, idAold, idBold, idAgm, idBgm;
0192   double mA, mB, eCM, s, sigmaMxGm;
0193   bool   hasLeptonBeamA, hasLeptonBeamB, hasOneLeptonBeam, hasTwoLeptonBeams,
0194          hasPointGammaA, hasPointGammaB, hasOnePointParticle,
0195          hasTwoPointParticles, hasGamma, hasVMD;
0196 
0197   // Cross section information.
0198   bool   newSigmaMx, canModifySigma, canBiasSelection, canBias2Sel;
0199   int    gmZmode;
0200   double bias2SelPow, bias2SelRef, wtBW, sigmaNw, sigmaMx, sigmaPos,
0201          sigmaNeg, biasWt;
0202 
0203   // Process-specific kinematics properties, almost always available.
0204   double mHatMin, mHatMax, sHatMin, sHatMax, pTHatMin, pTHatMax,
0205          pT2HatMin, pT2HatMax;
0206 
0207   // Event-specific kinematics properties, almost always available.
0208   double x1H, x2H, m3, m4, m5, s3, s4, s5, mHat, sH, tH, uH, pAbs, p2Abs,
0209          pTH, theta, phi, betaZ;
0210   Vec4   pH[12];
0211   double mH[12];
0212 
0213   // Reselect decay products momenta isotropically in phase space.
0214   void decayKinematicsStep( Event& process, int iRes);
0215 
0216   // Much common code for normal 2 -> 1, 2 -> 2 and 2 -> 3 cases:
0217 
0218   // Determine how phase space should be sampled.
0219   void setup3Body();
0220   bool setupSampling123(bool is2, bool is3);
0221 
0222   // Select a trial kinematics phase space point.
0223   bool trialKin123(bool is2, bool is3, bool inEvent = true);
0224 
0225   // Presence and properties of any s-channel resonances.
0226   int    idResA, idResB;
0227   double mResA, mResB, GammaResA, GammaResB, tauResA, tauResB, widResA,
0228          widResB;
0229   bool   sameResMass;
0230 
0231   // Kinematics properties specific to 2 -> 1/2/3.
0232   bool   useMirrorWeight, hasNegZ, hasPosZ;
0233   double tau, y, z, tauMin, tauMax, yMax, zMin, zMax, ratio34, unity34,
0234          zNeg, zPos, wtTau, wtY, wtZ, wt3Body, runBW3H, runBW4H, runBW5H,
0235          intTau0, intTau1, intTau2, intTau3, intTau4, intTau5, intTau6,
0236          intY0, intY12, intY34, intY56, mTchan1, sTchan1, mTchan2, sTchan2,
0237          frac3Flat, frac3Pow1, frac3Pow2, zNegMin, zNegMax, zPosMin, zPosMax;
0238   Vec4   p3cm, p4cm, p5cm;
0239 
0240   // Coefficients for optimized selection in 2 -> 1/2/3.
0241   int    nTau, nY, nZ;
0242   double tauCoef[8], yCoef[8], zCoef[8], tauCoefSum[8], yCoefSum[8],
0243          zCoefSum[8];
0244 
0245   // Calculate kinematical limits for 2 -> 1/2/3.
0246   bool limitTau(bool is2, bool is3);
0247   bool limitY();
0248   bool limitZ();
0249 
0250   // Select kinematical variable between defined limits for 2 -> 1/2/3.
0251   void selectTau(int iTau, double tauVal, bool is2);
0252   void selectY(int iY, double yVal);
0253   void selectZ(int iZ, double zVal);
0254   bool select3Body();
0255 
0256   // Solve equation system for better phase space coefficients in 2 -> 1/2/3.
0257   void solveSys( int n, int bin[8], double vec[8], double mat[8][8],
0258     double coef[8]);
0259 
0260   // Properties specific to resonance mass selection in 2 -> 2 and 2 -> 3.
0261   bool   useBW[6], useNarrowBW[6];
0262   int    idMass[6];
0263   double mPeak[6], sPeak[6], mWidth[6], mMin[6], mMax[6], mw[6], wmRat[6],
0264          mLower[6], mUpper[6], sLower[6], sUpper[6], fracFlatS[6],
0265          fracFlatM[6], fracInv[6], fracInv2[6], atanLower[6], atanUpper[6],
0266          intBW[6], intFlatS[6], intFlatM[6], intInv[6], intInv2[6];
0267 
0268   // Properties specific to top threshold enhancement.
0269   bool   doTopPair;
0270   int    topThresholdModel;
0271   double topThresholdWidth, eThreshold, m3Threshold, m4Threshold;
0272 
0273   // Setup mass selection for one resonance at a time. Split in two parts.
0274   void   setupMass1(int iM);
0275   void   setupMass2(int iM, double distToThresh);
0276 
0277   // Do mass selection and find the associated weight.
0278   void   trialMass(int iM);
0279   double weightMass(int iM);
0280 
0281   // Standard methods to find t range of a 2 -> 2 process
0282   // and to check whether a given t value is in that range.
0283   pair<double,double> tRange( double sIn, double s1In, double s2In,
0284     double s3In,  double s4In) {
0285     double lambda12 = pow2( sIn - s1In - s2In) - 4. * s1In * s2In;
0286     double lambda34 = pow2( sIn - s3In - s4In) - 4. * s3In * s4In;
0287     if (lambda12 < 0. || lambda34 < 0.) return make_pair( 0., 0.);
0288     double tLow = -0.5 * (sIn - (s1In + s2In + s3In + s4In) + (s1In - s2In)
0289       * (s3In - s4In) / sIn + sqrtpos(lambda12 *  lambda34) / sIn);
0290     double tUpp = ( (s3In - s1In) * (s4In - s2In) + (s1In + s4In - s2In - s3In)
0291       * (s1In * s4In - s2In * s3In) / sIn ) / tLow;
0292     return make_pair( tLow, tUpp); }
0293   bool tInRange( double tIn, double sIn, double s1In, double s2In,
0294     double s3In,  double s4In) {
0295     pair<double, double> tRng = tRange( sIn, s1In, s2In, s3In, s4In);
0296     return (tIn > tRng.first && tIn < tRng.second); }
0297 
0298   // The error function erf(x) should normally be in your math library,
0299   // but if not uncomment this simple parametrization by Sergei Winitzki.
0300   //double erf(double x) { double x2 = x * x; double kx2 = 0.147 * x2;
0301   //  double tmp = sqrt(1. - exp(-x2 * (4./M_PI + kx2) / (1. + kx2)));
0302   //  return ((x >= 0.) ? tmp : -tmp); }
0303 
0304 };
0305 
0306 //==========================================================================
0307 
0308 // A derived class with 2 -> 1 kinematics set up in tau, y.
0309 
0310 class PhaseSpace2to1tauy : public PhaseSpace {
0311 
0312 public:
0313 
0314   // Constructor.
0315   PhaseSpace2to1tauy() {}
0316 
0317   // Optimize subsequent kinematics selection.
0318   virtual bool setupSampling() {if (!setupMass()) return false;
0319     return setupSampling123(false, false);}
0320 
0321   // Construct the trial kinematics.
0322   virtual bool trialKin(bool inEvent = true, bool = false) {wtBW = 1.;
0323     return trialKin123(false, false, inEvent);}
0324 
0325   // Construct the final event kinematics.
0326   virtual bool finalKin();
0327 
0328 private:
0329 
0330   // Set up allowed mass range.
0331   bool setupMass();
0332 
0333 };
0334 
0335 //==========================================================================
0336 
0337 // A derived class with 2 -> 2 kinematics set up in tau, y, z = cos(theta).
0338 
0339 class PhaseSpace2to2tauyz : public PhaseSpace {
0340 
0341 public:
0342 
0343   // Constructor.
0344   PhaseSpace2to2tauyz() {}
0345 
0346   // Optimize subsequent kinematics selection.
0347   virtual bool setupSampling() {if (!setupMasses()) return false;
0348     return setupSampling123(true, false);}
0349 
0350   // Construct the trial kinematics.
0351   virtual bool trialKin(bool inEvent = true, bool = false) {
0352     if (!trialMasses()) return false;
0353     return trialKin123(true, false, inEvent);}
0354 
0355   // Construct the final event kinematics.
0356   virtual bool finalKin();
0357 
0358   // Rescales the momenta of incoming and outgoing partons according to
0359   // new sHat.
0360   virtual void rescaleMomenta ( double sHatNew);
0361 
0362   // Recalculates cross section with rescaled sHat.
0363   virtual void rescaleSigma ( double sHatNew);
0364 
0365   // Calculate the weight for over-estimated cross section.
0366   virtual double weightGammaPDFApprox();
0367 
0368 private:
0369 
0370   // Set up for fixed or Breit-Wigner mass selection.
0371   bool setupMasses();
0372 
0373   // Select fixed or Breit-Wigner-distributed masses.
0374   bool trialMasses();
0375 
0376   // Pick off-shell initialization masses when on-shell not allowed.
0377   bool constrainedM3M4();
0378   bool constrainedM3();
0379   bool constrainedM4();
0380 
0381 };
0382 
0383 //==========================================================================
0384 
0385 // A derived class with 2 -> 2 kinematics set up for elastic scattering.
0386 
0387 class PhaseSpace2to2elastic : public PhaseSpace {
0388 
0389 public:
0390 
0391   // Constructor.
0392   PhaseSpace2to2elastic() : isOneExp(), useCoulomb(), s1(), s2(), alphaEM0(),
0393     lambda12S(), lambda12(), lambda34(), tempA(), tempB(), tempC(),
0394     tLow(), tUpp(), bSlope1(), bSlope2(), sigRef1(), sigRef2(),
0395     sigRef(), sigNorm1(), sigNorm2(), sigNorm3(), sigNormSum(), rel2() {}
0396 
0397   // Construct the trial or final event kinematics.
0398   virtual bool setupSampling();
0399   virtual bool trialKin(bool inEvent = true, bool = false);
0400   virtual bool finalKin();
0401 
0402   // Are beam particles resolved in partons or scatter directly?
0403   virtual bool isResolved() const {return false;}
0404 
0405 private:
0406 
0407   // Constants: could only be changed in the code itself.
0408   static const int    NTRY;
0409   static const double BNARROW, BWIDE, WIDEFRAC, TOFFSET;
0410 
0411   // Kinematics properties specific to 2 -> 2 elastic.
0412   bool   isOneExp, useCoulomb;
0413   double s1, s2, alphaEM0, lambda12S, lambda12, lambda34, tempA, tempB, tempC,
0414          tLow, tUpp, bSlope1, bSlope2, sigRef1, sigRef2, sigRef,
0415          sigNorm1, sigNorm2, sigNorm3, sigNormSum, rel2;
0416 
0417 };
0418 
0419 //==========================================================================
0420 
0421 // A derived class with 2 -> 2 kinematics set up for diffractive scattering.
0422 
0423 class PhaseSpace2to2diffractive : public PhaseSpace {
0424 
0425 public:
0426 
0427   // Constructor.
0428   PhaseSpace2to2diffractive(bool isDiffAin = false, bool isDiffBin = false)
0429     : isDiffA(isDiffAin), isDiffB(isDiffBin), splitxit(), m3ElDiff(),
0430     m4ElDiff(), s1(), s2(), xiMin(), xiMax(), xiNow(), sigNow(), sigMax(),
0431     sigMaxNow(), lambda12(), lambda34(), bNow(), tempA(), tempB(), tempC(),
0432     tLow(), tUpp(), tWeight(), fWid1(), fWid2(), fWid3(), fWid4(), fbWid1(),
0433     fbWid2(), fbWid3(), fbWid4(), fbWid1234() {isSD = !isDiffA || !isDiffB;}
0434 
0435   // Construct the trial or final event kinematics.
0436   virtual bool setupSampling();
0437   virtual bool trialKin(bool inEvent = true, bool = false);
0438   virtual bool finalKin();
0439 
0440   // Are beam particles resolved in partons or scatter directly?
0441   virtual bool isResolved() const {return false;}
0442 
0443 private:
0444 
0445   // Constants: could only be changed in the code itself.
0446   static const int    NTRY;
0447   static const double BWID1, BWID2, BWID3, BWID4, FWID1SD, FWID2SD, FWID3SD,
0448                       FWID4SD, FWID1DD, FWID2DD, FWID3DD, FWID4DD, MAXFUDGESD,
0449                       MAXFUDGEDD, MAXFUDGET, DIFFMASSMARGIN, SPROTON;
0450 
0451   // Initialization data.
0452   bool   isDiffA, isDiffB, isSD, splitxit;
0453 
0454   // Pion mass, cached for efficiency.
0455   double mPi;
0456 
0457   // Kinematics properties specific to 2 -> 2 diffraction.
0458   double m3ElDiff, m4ElDiff, s1, s2, xiMin, xiMax, xiNow, sigNow, sigMax,
0459          sigMaxNow, lambda12, lambda34, bNow, tempA, tempB, tempC,
0460          tLow, tUpp, tWeight, fWid1, fWid2, fWid3, fWid4, fbWid1, fbWid2,
0461          fbWid3, fbWid4, fbWid1234;
0462 
0463 };
0464 
0465 //==========================================================================
0466 
0467 // A derived class with 2 -> 3 kinematics set up for central diffractive
0468 // scattering.
0469 
0470 class PhaseSpace2to3diffractive : public PhaseSpace {
0471 
0472 public:
0473 
0474   // Constructor.
0475  PhaseSpace2to3diffractive() : PhaseSpace(), splitxit(), s1(), s2(), m5min(),
0476     s5min(), sigNow(), sigMax(), sigMaxNow(), xiMin(), xi1(), xi2(), fWid1(),
0477     fWid2(), fWid3(), fbWid1(), fbWid2(), fbWid3(), fbWid123() {}
0478 
0479   // Construct the trial or final event kinematics.
0480   virtual bool setupSampling();
0481   virtual bool trialKin(bool inEvent = true, bool = false);
0482   virtual bool finalKin();
0483 
0484   // Are beam particles resolved in partons or scatter directly?
0485   virtual bool isResolved() const {return false;}
0486 
0487  private:
0488 
0489   // Constants: could only be changed in the code itself.
0490   static const int    NTRY;
0491   static const double BWID1, BWID2, BWID3, BWID4, FWID1, FWID2, FWID3,
0492                       MAXFUDGECD, MAXFUDGET, DIFFMASSMARGIN;
0493 
0494   // Initialization data.
0495   bool   splitxit;
0496 
0497   // Local variables to calculate DPE kinematics.
0498   double s1, s2, m5min, s5min, sigNow, sigMax, sigMaxNow, xiMin, xi1, xi2,
0499          fWid1, fWid2, fWid3, fbWid1, fbWid2, fbWid3, fbWid123;
0500   Vec4   p1, p2, p3, p4, p5;
0501 
0502 };
0503 
0504 //==========================================================================
0505 
0506 class PhaseSpace2to2nondiffractive : public PhaseSpace {
0507 
0508 // A derived class for nondiffractive events. Hardly does anything, since
0509 // the real action is taken care of by the MultipartonInteractions class.
0510 
0511 public:
0512 
0513   // Constructor.
0514   PhaseSpace2to2nondiffractive() {}
0515 
0516   // Construct the trial or final event kinematics.
0517   virtual bool setupSampling();
0518   virtual bool trialKin( bool , bool = false);
0519   virtual bool finalKin() { if (hasGamma) gammaKinPtr->finalize();
0520     return true;}
0521 
0522 private:
0523 
0524 };
0525 
0526 //==========================================================================
0527 
0528 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
0529 // tau, y, pT2_4, pT2_5, phi_4, phi_5 and y_3 (partial cylindrical symmetry).
0530 
0531 class PhaseSpace2to3tauycyl : public PhaseSpace {
0532 
0533 public:
0534 
0535   // Constructor.
0536   PhaseSpace2to3tauycyl() {}
0537 
0538   // Optimize subsequent kinematics selection.
0539   virtual bool setupSampling() {if (!setupMasses()) return false;
0540     setup3Body(); return setupSampling123(false, true);}
0541 
0542   // Construct the trial kinematics.
0543   virtual bool trialKin(bool inEvent = true, bool = false) {
0544     if (!trialMasses()) return false;
0545     return trialKin123(false, true, inEvent);}
0546 
0547   // Construct the final event kinematics.
0548   virtual bool finalKin();
0549 
0550 private:
0551 
0552   // Constants: could only be changed in the code itself.
0553   static const int    NITERNR;
0554 
0555   // Set up for fixed or Breit-Wigner mass selection.
0556   bool setupMasses();
0557 
0558   // Select fixed or Breit-Wigner-distributed masses.
0559   bool trialMasses();
0560 
0561 };
0562 
0563 //==========================================================================
0564 
0565 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
0566 // y3, y4, y5, pT2_3, pT2_5, phi_3 and phi_5, and with R separation cut.
0567 // Intended specifically for (essentially massless) 2 -> 3 QCD processes.
0568 
0569 class PhaseSpace2to3yyycyl : public PhaseSpace {
0570 
0571 public:
0572 
0573   // Constructor.
0574   PhaseSpace2to3yyycyl() : pTHat3Min(), pTHat3Max(), pTHat5Min(), pTHat5Max(),
0575     RsepMin(), R2sepMin(), hasBaryonBeams(), pT3Min(), pT3Max(), pT5Min(),
0576     pT5Max(), y3Max(), y4Max(), y5Max(), pT3(), pT4(), pT5(), phi3(), phi4(),
0577     phi5(), y3(), y4(), y5(), dphi() {}
0578 
0579   // Optimize subsequent kinematics selection.
0580   virtual bool setupSampling();
0581 
0582   // Construct the trial kinematics.
0583   virtual bool trialKin(bool inEvent = true, bool = false);
0584 
0585   // Construct the final event kinematics.
0586   virtual bool finalKin();
0587 
0588 private:
0589 
0590   // Phase space cuts specifically for 2 -> 3 QCD processes.
0591   double pTHat3Min, pTHat3Max, pTHat5Min, pTHat5Max, RsepMin, R2sepMin;
0592   bool   hasBaryonBeams;
0593 
0594   // Event kinematics choices.
0595   double pT3Min, pT3Max, pT5Min, pT5Max, y3Max, y4Max, y5Max,
0596          pT3, pT4, pT5, phi3, phi4, phi5, y3, y4, y5, dphi;
0597   Vec4   pInSum;
0598 
0599 };
0600 
0601 //==========================================================================
0602 
0603 // A derived class for Les Houches events.
0604 
0605 class PhaseSpaceLHA : public PhaseSpace {
0606 
0607 public:
0608 
0609   // Constructor.
0610   PhaseSpaceLHA() : strategy(), stratAbs(), nProc(), idProcSave(0),
0611     xMaxAbsSum(), xSecSgnSum(), sigmaSgn() {}
0612 
0613   // Find maximal cross section for comparison with internal processes.
0614   virtual bool setupSampling();
0615 
0616   // Construct the next process, by interface to Les Houches class.
0617   virtual bool trialKin( bool , bool repeatSame = false);
0618 
0619   // Set scale, alpha_s and alpha_em if not done.
0620   virtual bool finalKin() {sigmaProcessPtr->setScale(); return true;}
0621 
0622   // For Les Houches with negative event weight needs
0623   virtual double sigmaSumSigned() const {return sigmaSgn;}
0624 
0625 private:
0626 
0627   // Constants.
0628   static const double CONVERTPB2MB;
0629 
0630   // Local properties.
0631   int    strategy, stratAbs, nProc, idProcSave;
0632   double xMaxAbsSum, xSecSgnSum, sigmaSgn;
0633   vector<int>    idProc;
0634   vector<double> xMaxAbsProc;
0635 
0636 };
0637 
0638 //==========================================================================
0639 
0640 // Rambo flat phase-space generator.
0641 
0642 // This is an implementation of the Rambo phase-space generator as
0643 // presented in A New Monte Carlo Treatment Of Multiparticle Phase
0644 // Space At High-Energies, R. Kleiss, W.J. Stirling, S.D. Ellis, CPC40
0645 // (1986) 359.
0646 
0647 class Rambo {
0648 
0649  public:
0650 
0651   // Deafult constructor.
0652   Rambo() { rndmPtr=nullptr; isInitPtr=false;}
0653 
0654   // Initializing constructor.
0655   Rambo(Rndm* rndmPtrIn) { initPtr(rndmPtrIn); }
0656 
0657   // Destructor.
0658   virtual ~Rambo() {}
0659 
0660   // Initialize pointers.
0661   void initPtr(Rndm* rndmPtrIn) {rndmPtr = rndmPtrIn; isInitPtr = true;}
0662 
0663   // Rambo phase space generator. Generates nOut uniformly distributed
0664   // massless 4-vectors with sqrt(s) = eCM. Output in pOut.
0665   double genPoint(double eCM, int nOut, vector<Vec4>& pOut);
0666 
0667   // Massive generalisation, weights NOT 1 anymore - literal implementation
0668   // of original RAMBO paper by Ellis, Kleiss and Stirling. Number of particles
0669   // determined from size of mIn vector.
0670   double genPoint(double eCM, vector<double> mIn, vector<Vec4>& pOut);
0671 
0672  private:
0673 
0674   // Is initialized.
0675   bool isInitPtr;
0676 
0677   // Pointer to the random number generator.
0678   Rndm*  rndmPtr;
0679 
0680 };
0681 
0682 //==========================================================================
0683 
0684 } // end namespace Pythia8
0685 
0686 #endif // Pythia8_PhaseSpace_H