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) 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 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 nondiffractive, 2 -> 3, Les Houches.
10 
11 #ifndef Pythia8_PhaseSpace_H
12 #define Pythia8_PhaseSpace_H
13 
14 #include "Pythia8/Basics.h"
15 #include "Pythia8/BeamParticle.h"
16 #include "Pythia8/Info.h"
17 #include "Pythia8/LesHouches.h"
18 #include "Pythia8/MultipartonInteractions.h"
19 #include "Pythia8/ParticleData.h"
20 #include "Pythia8/PartonDistributions.h"
21 #include "Pythia8/PythiaStdlib.h"
22 #include "Pythia8/SigmaProcess.h"
23 #include "Pythia8/SigmaTotal.h"
24 #include "Pythia8/Settings.h"
25 #include "Pythia8/StandardModel.h"
26 #include "Pythia8/UserHooks.h"
27 #include "Pythia8/GammaKinematics.h"
28 
29 namespace Pythia8 {
30 
31 //==========================================================================
32 
33 // Forward reference to the UserHooks class.
34 class UserHooks;
35 
36 //==========================================================================
37 
38 // PhaseSpace is a base class for phase space generators
39 // used in the selection of hard-process kinematics.
40 
41 class PhaseSpace {
42 
43 public:
44 
45  // Destructor.
46  virtual ~PhaseSpace() {}
47 
48  // Perform simple initialization and store pointers.
49  void init(bool isFirst, SigmaProcess* sigmaProcessPtrIn,
50  Info* infoPtrIn, Settings* settingsPtrIn, ParticleData* particleDataPtrIn,
51  Rndm* rndmPtrIn, BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
52  Couplings* couplingsPtrIn, SigmaTotal* sigmaTotPtrIn,
53  UserHooks* userHooksPtrIn);
54 
55  // Update the CM energy of the event.
56  void newECM(double eCMin) {eCM = eCMin; s = eCM * eCM;}
57 
58  // Store or replace Les Houches pointer.
59  void setLHAPtr(LHAup* lhaUpPtrIn) {lhaUpPtr = lhaUpPtrIn;}
60 
61  // A pure virtual method, wherein an optimization procedure
62  // is used to determine how phase space should be sampled.
63  virtual bool setupSampling() = 0;
64 
65  // A pure virtual method, wherein a trial event kinematics
66  // is to be selected in the derived class.
67  virtual bool trialKin(bool inEvent = true, bool repeatSame = false) = 0;
68 
69  // A pure virtual method, wherein the accepted event kinematics
70  // is to be constructed in the derived class.
71  virtual bool finalKin() = 0;
72 
73  // Allow for nonisotropic decays when ME's available.
74  void decayKinematics( Event& process);
75 
76  // Give back current or maximum cross section, or set latter.
77  double sigmaNow() const {return sigmaNw;}
78  double sigmaMax() const {return sigmaMx;}
79  double biasSelectionWeight() const {return biasWt;}
80  bool newSigmaMax() const {return newSigmaMx;}
81  void setSigmaMax(double sigmaMaxIn) {sigmaMx = sigmaMaxIn;}
82 
83  // For Les Houches with negative event weight needs
84  virtual double sigmaSumSigned() const {return sigmaMx;}
85 
86  // Give back constructed four-vectors and known masses.
87  Vec4 p(int i) const {return pH[i];}
88  double m(int i) const {return mH[i];}
89 
90  // Reset the four-momentum.
91  void setP(int i, Vec4 pNew) { pH[i] = pNew; }
92 
93  // Give back other event properties.
94  double ecm() const {return eCM;}
95  double x1() const {return x1H;}
96  double x2() const {return x2H;}
97  double sHat() const {return sH;}
98  double tHat() const {return tH;}
99  double uHat() const {return uH;}
100  double pTHat() const {return pTH;}
101  double thetaHat() const {return theta;}
102  double phiHat() const {return phi;}
103  double runBW3() const {return runBW3H;}
104  double runBW4() const {return runBW4H;}
105  double runBW5() const {return runBW5H;}
106 
107  // Inform whether beam particles are resolved in partons or scatter directly.
108  virtual bool isResolved() const {return true;}
109 
110  // Functions to rescale momenta and cross section for new sHat
111  // Currently implemented only for PhaseSpace2to2tauyz class.
112  virtual void rescaleSigma( double){}
113  virtual void rescaleMomenta( double){}
114 
115  // Calculate the weight for over-estimated cross section.
116  virtual double weightGammaPDFApprox(){ return 1.;}
117 
118  // Set the GammaKinematics pointer. Implemented for non-diffractive gm+gm.
119  virtual void setGammaKinPtr( GammaKinematics*){}
120 
121 protected:
122 
123  // Constructor.
124  PhaseSpace() {}
125 
126  // Constants: could only be changed in the code itself.
127  static const int NMAXTRY, NTRY3BODY;
128  static const double SAFETYMARGIN, TINY, EVENFRAC, SAMESIGMA, MRESMINABS,
129  WIDTHMARGIN, SAMEMASS, MASSMARGIN, EXTRABWWTMAX,
130  THRESHOLDSIZE, THRESHOLDSTEP, YRANGEMARGIN, LEPTONXMIN,
131  LEPTONXMAX, LEPTONXLOGMIN, LEPTONXLOGMAX, LEPTONTAUMIN,
132  SHATMINZ, PT2RATMINZ, WTCORRECTION[11];
133 
134  // Pointer to cross section.
135  SigmaProcess* sigmaProcessPtr;
136 
137  // Pointer to various information on the generation.
138  Info* infoPtr;
139 
140  // Pointer to the settings database.
141  Settings* settingsPtr;
142 
143  // Pointer to the particle data table.
144  ParticleData* particleDataPtr;
145 
146  // Pointer to the random number generator.
147  Rndm* rndmPtr;
148 
149  // Pointers to incoming beams.
150  BeamParticle* beamAPtr;
151  BeamParticle* beamBPtr;
152 
153  // Pointer to Standard Model couplings.
154  Couplings* couplingsPtr;
155 
156  // Pointer to the total/elastic/diffractive cross section object.
157  SigmaTotal* sigmaTotPtr;
158 
159  // Pointer to userHooks object for user interaction with program.
160  UserHooks* userHooksPtr;
161 
162  // Pointer to LHAup for generating external events.
163  LHAup* lhaUpPtr;
164 
165  // Initialization data, normally only set once.
166  bool useBreitWigners, doEnergySpread, showSearch, showViolation,
167  increaseMaximum, hasQ2Min;
168  int gmZmodeGlobal;
169  double mHatGlobalMin, mHatGlobalMax, pTHatGlobalMin, pTHatGlobalMax,
170  Q2GlobalMin, pTHatMinDiverge, minWidthBreitWigners;
171 
172  // Information on incoming beams.
173  int idA, idB;
174  double mA, mB, eCM, s;
175  bool hasLeptonBeamA, hasLeptonBeamB, hasOneLeptonBeam, hasTwoLeptonBeams,
176  hasPointGammaA, hasPointGammaB, hasOnePointParticle,
177  hasTwoPointParticles;
178 
179  // Cross section information.
180  bool newSigmaMx, canModifySigma, canBiasSelection, canBias2Sel;
181  int gmZmode;
182  double bias2SelPow, bias2SelRef, wtBW, sigmaNw, sigmaMx, sigmaPos,
183  sigmaNeg, biasWt;
184 
185  // Process-specific kinematics properties, almost always available.
186  double mHatMin, mHatMax, sHatMin, sHatMax, pTHatMin, pTHatMax,
187  pT2HatMin, pT2HatMax;
188 
189  // Event-specific kinematics properties, almost always available.
190  double x1H, x2H, m3, m4, m5, s3, s4, s5, mHat, sH, tH, uH, pAbs, p2Abs,
191  pTH, theta, phi, betaZ;
192  Vec4 pH[12];
193  double mH[12];
194 
195  // Reselect decay products momenta isotropically in phase space.
196  void decayKinematicsStep( Event& process, int iRes);
197 
198  // Much common code for normal 2 -> 1, 2 -> 2 and 2 -> 3 cases:
199 
200  // Determine how phase space should be sampled.
201  void setup3Body();
202  bool setupSampling123(bool is2, bool is3);
203 
204  // Select a trial kinematics phase space point.
205  bool trialKin123(bool is2, bool is3, bool inEvent = true);
206 
207  // Presence and properties of any s-channel resonances.
208  int idResA, idResB;
209  double mResA, mResB, GammaResA, GammaResB, tauResA, tauResB, widResA,
210  widResB;
211  bool sameResMass;
212 
213  // Kinematics properties specific to 2 -> 1/2/3.
214  bool useMirrorWeight, hasNegZ, hasPosZ;
215  double tau, y, z, tauMin, tauMax, yMax, zMin, zMax, ratio34, unity34,
216  zNeg, zPos, wtTau, wtY, wtZ, wt3Body, runBW3H, runBW4H, runBW5H,
217  intTau0, intTau1, intTau2, intTau3, intTau4, intTau5, intTau6,
218  intY0, intY12, intY34, intY56, mTchan1, sTchan1, mTchan2, sTchan2,
219  frac3Flat, frac3Pow1, frac3Pow2, zNegMin, zNegMax, zPosMin, zPosMax;
220  Vec4 p3cm, p4cm, p5cm;
221 
222  // Coefficients for optimized selection in 2 -> 1/2/3.
223  int nTau, nY, nZ;
224  double tauCoef[8], yCoef[8], zCoef[8], tauCoefSum[8], yCoefSum[8],
225  zCoefSum[8];
226 
227  // Calculate kinematical limits for 2 -> 1/2/3.
228  bool limitTau(bool is2, bool is3);
229  bool limitY();
230  bool limitZ();
231 
232  // Select kinematical variable between defined limits for 2 -> 1/2/3.
233  void selectTau(int iTau, double tauVal, bool is2);
234  void selectY(int iY, double yVal);
235  void selectZ(int iZ, double zVal);
236  bool select3Body();
237 
238  // Solve equation system for better phase space coefficients in 2 -> 1/2/3.
239  void solveSys( int n, int bin[8], double vec[8], double mat[8][8],
240  double coef[8]);
241 
242  // Properties specific to resonance mass selection in 2 -> 2 and 2 -> 3.
243  bool useBW[6];
244  int idMass[6];
245  double mPeak[6], sPeak[6], mWidth[6], mMin[6], mMax[6], mw[6], wmRat[6],
246  mLower[6], mUpper[6], sLower[6], sUpper[6], fracFlatS[6],
247  fracFlatM[6], fracInv[6], fracInv2[6], atanLower[6], atanUpper[6],
248  intBW[6], intFlatS[6], intFlatM[6], intInv[6], intInv2[6];
249 
250  // Setup mass selection for one resonance at a time. Split in two parts.
251  void setupMass1(int iM);
252  void setupMass2(int iM, double distToThresh);
253 
254  // Do mass selection and find the associated weight.
255  void trialMass(int iM);
256  double weightMass(int iM);
257 
258  // Standard methods to find t range of a 2 -> 2 process
259  // and to check whether a given t value is in that range.
260  pair<double,double> tRange( double sIn, double s1In, double s2In,
261  double s3In, double s4In) {
262  double lambda12 = pow2( sIn - s1In - s2In) - 4. * s1In * s2In;
263  double lambda34 = pow2( sIn - s3In - s4In) - 4. * s3In * s4In;
264  if (lambda12 < 0. || lambda34 < 0.) return make_pair( 0., 0.);
265  double tLow = -0.5 * (sIn - (s1In + s2In + s3In + s4In) + (s1In - s2In)
266  * (s3In - s4In) / sIn + sqrtpos(lambda12 * lambda34) / sIn);
267  double tUpp = ( (s3In - s1In) * (s4In - s2In) + (s1In + s4In - s2In - s3In)
268  * (s1In * s4In - s2In * s3In) / sIn ) / tLow;
269  return make_pair( tLow, tUpp); }
270  bool tInRange( double tIn, double sIn, double s1In, double s2In,
271  double s3In, double s4In) {
272  pair<double, double> tRng = tRange( sIn, s1In, s2In, s3In, s4In);
273  return (tIn > tRng.first && tIn < tRng.second); }
274 
275  // The error function erf(x) should normally be in your math library,
276  // but if not uncomment this simple parametrization by Sergei Winitzki.
277  //double erf(double x) { double x2 = x * x; double kx2 = 0.147 * x2;
278  // double tmp = sqrt(1. - exp(-x2 * (4./M_PI + kx2) / (1. + kx2)));
279  // return ((x >= 0.) ? tmp : -tmp); }
280 
281 };
282 
283 //==========================================================================
284 
285 // A derived class with 2 -> 1 kinematics set up in tau, y.
286 
287 class PhaseSpace2to1tauy : public PhaseSpace {
288 
289 public:
290 
291  // Constructor.
292  PhaseSpace2to1tauy() {}
293 
294  // Optimize subsequent kinematics selection.
295  virtual bool setupSampling() {if (!setupMass()) return false;
296  return setupSampling123(false, false);}
297 
298  // Construct the trial kinematics.
299  virtual bool trialKin(bool inEvent = true, bool = false) {wtBW = 1.;
300  return trialKin123(false, false, inEvent);}
301 
302  // Construct the final event kinematics.
303  virtual bool finalKin();
304 
305 private:
306 
307  // Set up allowed mass range.
308  bool setupMass();
309 
310 };
311 
312 //==========================================================================
313 
314 // A derived class with 2 -> 2 kinematics set up in tau, y, z = cos(theta).
315 
316 class PhaseSpace2to2tauyz : public PhaseSpace {
317 
318 public:
319 
320  // Constructor.
321  PhaseSpace2to2tauyz() {}
322 
323  // Optimize subsequent kinematics selection.
324  virtual bool setupSampling() {if (!setupMasses()) return false;
325  return setupSampling123(true, false);}
326 
327  // Construct the trial kinematics.
328  virtual bool trialKin(bool inEvent = true, bool = false) {
329  if (!trialMasses()) return false;
330  return trialKin123(true, false, inEvent);}
331 
332  // Construct the final event kinematics.
333  virtual bool finalKin();
334 
335  // Rescales the momenta of incoming and outgoing partons according to
336  // new sHat.
337  virtual void rescaleMomenta ( double sHatNew);
338 
339  // Recalculates cross section with rescaled sHat.
340  virtual void rescaleSigma ( double sHatNew);
341 
342  // Calculate the weight for over-estimated cross section.
343  virtual double weightGammaPDFApprox();
344 
345 private:
346 
347  // Set up for fixed or Breit-Wigner mass selection.
348  bool setupMasses();
349 
350  // Select fixed or Breit-Wigner-distributed masses.
351  bool trialMasses();
352 
353  // Pick off-shell initialization masses when on-shell not allowed.
354  bool constrainedM3M4();
355  bool constrainedM3();
356  bool constrainedM4();
357 
358 };
359 
360 //==========================================================================
361 
362 // A derived class with 2 -> 2 kinematics set up for elastic scattering.
363 
364 class PhaseSpace2to2elastic : public PhaseSpace {
365 
366 public:
367 
368  // Constructor.
369  PhaseSpace2to2elastic() {}
370 
371  // Construct the trial or final event kinematics.
372  virtual bool setupSampling();
373  virtual bool trialKin(bool inEvent = true, bool = false);
374  virtual bool finalKin();
375 
376  // Are beam particles resolved in partons or scatter directly?
377  virtual bool isResolved() const {return false;}
378 
379 private:
380 
381  // Constants: could only be changed in the code itself.
382  static const int NTRY;
383  static const double HBARC2, BNARROW, BWIDE, WIDEFRAC, TOFFSET;
384 
385  // Kinematics properties specific to 2 -> 2 elastic.
386  bool isOneExp, useCoulomb;
387  double s1, s2, alphaEM0, lambda12S, tLow, tUpp, bSlope1, bSlope2,
388  sigRef1, sigRef2, sigRef, sigNorm1, sigNorm2, sigNorm3,
389  sigNormSum, rel2;
390 
391 };
392 
393 //==========================================================================
394 
395 // A derived class with 2 -> 2 kinematics set up for diffractive scattering.
396 
397 class PhaseSpace2to2diffractive : public PhaseSpace {
398 
399 public:
400 
401  // Constructor.
402  PhaseSpace2to2diffractive(bool isDiffAin = false, bool isDiffBin = false)
403  : isDiffA(isDiffAin), isDiffB(isDiffBin) {isSD = !isDiffA || !isDiffB;}
404 
405  // Construct the trial or final event kinematics.
406  virtual bool setupSampling();
407  virtual bool trialKin(bool inEvent = true, bool = false);
408  virtual bool finalKin();
409 
410  // Are beam particles resolved in partons or scatter directly?
411  virtual bool isResolved() const {return false;}
412 
413 private:
414 
415  // Constants: could only be changed in the code itself.
416  static const int NTRY;
417  static const double BWID1, BWID2, BWID3, BWID4, FWID1SD, FWID2SD, FWID3SD,
418  FWID4SD, FWID1DD, FWID2DD, FWID3DD, FWID4DD, MAXFUDGESD,
419  MAXFUDGEDD, MAXFUDGET, DIFFMASSMARGIN, SPROTON;
420 
421  // Initialization data.
422  bool isDiffA, isDiffB, isSD, splitxit;
423 
424  // Kinematics properties specific to 2 -> 2 diffraction.
425  double m3ElDiff, m4ElDiff, s1, s2, xiMin, xiMax, xiNow, sigNow, sigMax,
426  sigMaxNow, lambda12, lambda34, bNow, tempA, tempB, tempC,
427  tLow, tUpp, tWeight, fWid1, fWid2, fWid3, fWid4, fbWid1, fbWid2,
428  fbWid3, fbWid4, fbWid1234;
429 
430 };
431 
432 //==========================================================================
433 
434 // A derived class with 2 -> 3 kinematics set up for central diffractive
435 // scattering.
436 
437 class PhaseSpace2to3diffractive : public PhaseSpace {
438 
439 public:
440 
441  // Constructor.
442  PhaseSpace2to3diffractive() {}
443 
444  // Construct the trial or final event kinematics.
445  virtual bool setupSampling();
446  virtual bool trialKin(bool inEvent = true, bool = false);
447  virtual bool finalKin();
448 
449  // Are beam particles resolved in partons or scatter directly?
450  virtual bool isResolved() const {return false;}
451 
452  private:
453 
454  // Constants: could only be changed in the code itself.
455  static const int NTRY;
456  static const double BWID1, BWID2, BWID3, BWID4, FWID1, FWID2, FWID3,
457  MAXFUDGECD, MAXFUDGET, DIFFMASSMARGIN;
458 
459  // Initialization data.
460  bool splitxit;
461 
462  // Local variables to calculate DPE kinematics.
463  double s1, s2, m5min, s5min, m5, sigNow, sigMax, sigMaxNow, xiMin, xi1, xi2,
464  fWid1, fWid2, fWid3, fbWid1, fbWid2, fbWid3, fbWid123;
465  Vec4 p1, p2, p3, p4, p5;
466 
467 };
468 
469 //==========================================================================
470 
471 // A derived class for nondiffractive events. Hardly does anything, since
472 // the real action is taken care of by the MultipartonInteractions class.
473 
474 class PhaseSpace2to2nondiffractive : public PhaseSpace {
475 
476 public:
477 
478  // Constructor.
479  PhaseSpace2to2nondiffractive() {}
480 
481  // Construct the trial or final event kinematics.
482  virtual bool setupSampling() {sigmaNw = sigmaProcessPtr->sigmaHat();
483  sigmaMx = sigmaNw; return true;}
484  virtual bool trialKin( bool , bool = false) {return true;}
485  virtual bool finalKin() {return true;}
486 
487 private:
488 
489 };
490 
491 //==========================================================================
492 
493 // A derived class for nondiffractive events in l+l- -> gm+gm.
494 // The process is still generated in MultipartonInteraction but the sampling
495 // of the sub-collision is done here to allow cuts for the phase space.
496 
498 
499 public:
500 
501  // Constructor.
503 
504  // Construct the trial or final event kinematics.
505  virtual bool setupSampling();
506  virtual bool trialKin(bool inEvent = true, bool = false);
507  virtual bool finalKin() {gammaKinPtr->finalize(); return true;}
508 
509  // Set the pointer to GammaKinematics.
510  virtual void setGammaKinPtr( GammaKinematics* gammaKinPtrIn) {
511  gammaKinPtr = gammaKinPtrIn; }
512 
513 private:
514 
515  // Pointer to object that samples photon kinematics from leptons.
516  GammaKinematics* gammaKinPtr;
517 
518  // Parameters.
519  int idAin, idBin;
520  bool gammaA, gammaB, externalFlux, sampleQ2;
521  double Q2maxGamma, Wmin, sigmaNDestimate, sigmaNDmax, sCM, alphaEM,
522  m2BeamA, m2BeamB, m2sA, m2sB, log2xMinA, log2xMaxA, log2xMinB, log2xMaxB,
523  xGamma1, xGamma2, Q2gamma1, Q2gamma2, mGmGm, Q2min1, Q2min2;
524 
525 };
526 
527 //==========================================================================
528 
529 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
530 // tau, y, pT2_4, pT2_5, phi_4, phi_5 and y_3 (partial cylindrical symmetry).
531 
532 class PhaseSpace2to3tauycyl : public PhaseSpace {
533 
534 public:
535 
536  // Constructor.
538 
539  // Optimize subsequent kinematics selection.
540  virtual bool setupSampling() {if (!setupMasses()) return false;
541  setup3Body(); return setupSampling123(false, true);}
542 
543  // Construct the trial kinematics.
544  virtual bool trialKin(bool inEvent = true, bool = false) {
545  if (!trialMasses()) return false;
546  return trialKin123(false, true, inEvent);}
547 
548  // Construct the final event kinematics.
549  virtual bool finalKin();
550 
551 private:
552 
553  // Constants: could only be changed in the code itself.
554  static const int NITERNR;
555 
556  // Set up for fixed or Breit-Wigner mass selection.
557  bool setupMasses();
558 
559  // Select fixed or Breit-Wigner-distributed masses.
560  bool trialMasses();
561 
562 };
563 
564 //==========================================================================
565 
566 // A derived class with 2 -> 3 kinematics 1 + 2 -> 3 + 4 + 5 set up in
567 // y3, y4, y5, pT2_3, pT2_5, phi_3 and phi_5, and with R separation cut.
568 // Intended specifically for (essentially massless) 2 -> 3 QCD processes.
569 
570 class PhaseSpace2to3yyycyl : public PhaseSpace {
571 
572 public:
573 
574  // Constructor.
575  PhaseSpace2to3yyycyl() {}
576 
577  // Optimize subsequent kinematics selection.
578  virtual bool setupSampling();
579 
580  // Construct the trial kinematics.
581  virtual bool trialKin(bool inEvent = true, bool = false);
582 
583  // Construct the final event kinematics.
584  virtual bool finalKin();
585 
586 private:
587 
588  // Phase space cuts specifically for 2 -> 3 QCD processes.
589  double pTHat3Min, pTHat3Max, pTHat5Min, pTHat5Max, RsepMin, R2sepMin;
590  bool hasBaryonBeams;
591 
592  // Event kinematics choices.
593  double pT3Min, pT3Max, pT5Min, pT5Max, y3Max, y4Max, y5Max,
594  pT3, pT4, pT5, phi3, phi4, phi5, y3, y4, y5, dphi;
595  Vec4 pInSum;
596 
597 };
598 
599 //==========================================================================
600 
601 // A derived class for Les Houches events.
602 
603 class PhaseSpaceLHA : public PhaseSpace {
604 
605 public:
606 
607  // Constructor.
608  PhaseSpaceLHA() {idProcSave = 0;}
609 
610  // Find maximal cross section for comparison with internal processes.
611  virtual bool setupSampling();
612 
613  // Construct the next process, by interface to Les Houches class.
614  virtual bool trialKin( bool , bool repeatSame = false);
615 
616  // Set scale, alpha_s and alpha_em if not done.
617  virtual bool finalKin() {sigmaProcessPtr->setScale(); return true;}
618 
619  // For Les Houches with negative event weight needs
620  virtual double sigmaSumSigned() const {return sigmaSgn;}
621 
622 private:
623 
624  // Constants.
625  static const double CONVERTPB2MB;
626 
627  // Local properties.
628  int strategy, stratAbs, nProc, idProcSave;
629  double xMaxAbsSum, xSecSgnSum, sigmaSgn;
630  vector<int> idProc;
631  vector<double> xMaxAbsProc;
632 
633 };
634 
635 //==========================================================================
636 
637 } // end namespace Pythia8
638 
639 #endif // Pythia8_PhaseSpace_H
Definition: AgUStep.h:26