Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:29

0001 // ResonanceWidths.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2024 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 resonance properties: dynamical widths etc.
0007 // ResonanceWidths: base class for all resonances.
0008 // ResonanceGmZ, ...: derived classes for individual resonances.
0009 
0010 #ifndef Pythia8_ResonanceWidths_H
0011 #define Pythia8_ResonanceWidths_H
0012 
0013 #include "Pythia8/Basics.h"
0014 #include "Pythia8/Info.h"
0015 #include "Pythia8/PythiaStdlib.h"
0016 #include "Pythia8/Settings.h"
0017 #include "Pythia8/StandardModel.h"
0018 
0019 namespace Pythia8 {
0020 
0021 //==========================================================================
0022 
0023 // Forward references to ParticleData and StandardModel classes.
0024 class DecayChannel;
0025 class ParticleData;
0026 class ParticleDataEntry;
0027 class CoupSM;
0028 class CoupSUSY;
0029 
0030 //==========================================================================
0031 
0032 // The ResonanceWidths is the base class. Also used for generic resonaces.
0033 
0034 class ResonanceWidths {
0035 
0036 public:
0037 
0038   // Destructor.
0039   virtual ~ResonanceWidths() {}
0040 
0041   // Set up standard properties.
0042   void initBasic(int idResIn, bool isGenericIn = false) {
0043     idRes = idResIn; isGeneric = isGenericIn;}
0044 
0045   // Calculate and store partial and total widths at the nominal mass.
0046   virtual bool init(Info* infoPtrIn);
0047 
0048   // Return identity of particle species.
0049   int id() const {return idRes;}
0050 
0051   // Calculate the total/open width for given mass, charge and instate.
0052   double width(int idSgn, double mHatIn, int idInFlavIn = 0,
0053     bool openOnly = false, bool setBR = false, int idOutFlav1 = 0,
0054     int idOutFlav2 = 0);
0055 
0056   // Special case to calculate open final-state width.
0057   double widthOpen(int idSgn, double mHatIn, int idIn = 0) {
0058     return width( idSgn, mHatIn, idIn, true, false);}
0059 
0060   // Special case to store open final-state widths for channel selection.
0061   double widthStore(int idSgn, double mHatIn, int idIn = 0) {
0062     return width( idSgn, mHatIn, idIn, true, true);}
0063 
0064   // Return fraction of width open for particle and antiparticle.
0065   double openFrac(int idSgn) {return (idSgn > 0) ? openPos : openNeg;}
0066 
0067   // Return forced rescaling factor of resonance width.
0068   double widthRescaleFactor() {return forceFactor;}
0069 
0070   // Special case to calculate one final-state width.
0071   // Currently only used for Higgs -> qqbar, g g or gamma gamma.
0072   double widthChan(double mHatIn, int idOutFlav1, int idOutFlav2) {
0073     return width( 1, mHatIn, 0, false, false, idOutFlav1, idOutFlav2);}
0074 
0075 protected:
0076 
0077   // Constructor.
0078   ResonanceWidths() : idRes(), hasAntiRes(), doForceWidth(), isGeneric(),
0079     allowCalcWidth(), minWidth(), minThreshold(), mRes(), GammaRes(), m2Res(),
0080     GamMRat(), openPos(), openNeg(), forceFactor(), iChannel(), onMode(),
0081     meMode(), mult(), id1(), id2(), id3(), id1Abs(), id2Abs(), id3Abs(),
0082     idInFlav(), widNow(), mHat(), mf1(), mf2(), mf3(), mr1(), mr2(), mr3(),
0083     ps(), kinFac(), alpEM(), alpS(), colQ(), preFac(), particlePtr(),
0084     infoPtr(), settingsPtr(), particleDataPtr(), coupSMPtr(),
0085     coupSUSYPtr(){}
0086 
0087   // Constants: could only be changed in the code itself.
0088   static const int    NPOINT;
0089   static const double MASSMIN, MASSMARGIN;
0090 
0091   // Particle properties always present.
0092   int    idRes, hasAntiRes;
0093   bool   doForceWidth, isGeneric, allowCalcWidth;
0094   double minWidth, minThreshold, mRes, GammaRes, m2Res, GamMRat,
0095          openPos, openNeg, forceFactor;
0096 
0097   // Properties for currently studied decay channel(s).
0098   int    iChannel, onMode, meMode, mult, id1, id2, id3, id1Abs,
0099          id2Abs, id3Abs, idInFlav;
0100   double widNow, mHat, mf1, mf2, mf3, mr1, mr2, mr3, ps, kinFac,
0101          alpEM, alpS, colQ, preFac;
0102 
0103   // Pointer to properties of the particle species.
0104   weak_ptr<ParticleDataEntry> particlePtr;
0105 
0106   // Pointer to various information on the generation.
0107   Info*         infoPtr;
0108 
0109   // Pointer to the settings database.
0110   Settings*     settingsPtr;
0111 
0112   // Pointer to the particle data table.
0113   ParticleData* particleDataPtr;
0114 
0115   // Pointer to the logger.
0116   Logger*       loggerPtr;
0117 
0118   // Pointers to Standard Model and SUSY couplings.
0119   CoupSM*       coupSMPtr;
0120   CoupSUSY*     coupSUSYPtr;
0121 
0122   // Initialize constants.
0123   virtual void initConstants() {}
0124 
0125   // Virtual methods to handle model-specific (non-SM) part of initialization
0126   // for use by derived classes that implement additional models (eg SUSY).
0127   virtual bool initBSM() {return true;}
0128   virtual bool allowCalc() {return true;}
0129 
0130   // Calculate various common prefactors for the current mass.
0131   // Optional argument calledFromInit only used for Z0.
0132   virtual void calcPreFac(bool = false) {}
0133 
0134   // Calculate width for currently considered channel.
0135   // Optional argument calledFromInit only used for Z0.
0136   virtual void calcWidth(bool = false) {}
0137 
0138   // Simple routines for matrix-element integration over Breit-Wigners.
0139   double numInt1BW(double mHatIn, double m1, double Gamma1, double mMin1,
0140     double m2, int psMode = 1);
0141   double numInt2BW(double mHatIn, double m1, double Gamma1, double mMin1,
0142     double m2, double Gamma2, double mMin2, int psMode = 1);
0143 
0144 };
0145 
0146 //==========================================================================
0147 
0148 // The ResonanceGeneric class handles a generic resonance.
0149 // Only needs a constructor and allowCalc = false; for the rest uses
0150 // defaults in base class.
0151 
0152 class ResonanceGeneric : public ResonanceWidths {
0153 
0154 public:
0155 
0156   // Constructor.
0157   ResonanceGeneric(int idResIn) {initBasic(idResIn, true);}
0158 
0159   // By default, assume no dedicated code exists to compute width.
0160   virtual bool allowCalc() override {return false;}
0161 
0162 };
0163 
0164 //==========================================================================
0165 
0166 // The ResonanceGmZ class handles the gamma*/Z0 resonance.
0167 
0168 class ResonanceGmZ : public ResonanceWidths {
0169 
0170 public:
0171 
0172   // Constructor.
0173   ResonanceGmZ(int idResIn) : gmZmode(), thetaWRat(), ei2(), eivi(), vi2ai2(),
0174     gamNorm(), intNorm(), resNorm() {initBasic(idResIn);}
0175 
0176 private:
0177 
0178   // Locally stored properties and couplings.
0179   int    gmZmode;
0180   double thetaWRat, ei2, eivi, vi2ai2, gamNorm, intNorm, resNorm;
0181 
0182   // Initialize constants.
0183   virtual void initConstants() override;
0184 
0185   // Calculate various common prefactors for the current mass.
0186   virtual void calcPreFac(bool = false) override;
0187 
0188   // Caclulate width for currently considered channel.
0189   virtual void calcWidth(bool calledFromInit = false) override;
0190 
0191 };
0192 
0193 //==========================================================================
0194 
0195 // The ResonanceW class handles the W+- resonance.
0196 
0197 class ResonanceW : public ResonanceWidths {
0198 
0199 public:
0200 
0201   // Constructor.
0202  ResonanceW(int idResIn) : ResonanceWidths(), thetaWRat() {initBasic(idResIn);}
0203 
0204 private:
0205 
0206   // Locally stored properties and couplings.
0207   double thetaWRat;
0208 
0209   // Initialize constants.
0210   virtual void initConstants() override;
0211 
0212   // Calculate various common prefactors for the current mass.
0213   virtual void calcPreFac(bool = false) override;
0214 
0215   // Caclulate width for currently considered channel.
0216   virtual void calcWidth(bool = false) override;
0217 
0218 };
0219 
0220 //==========================================================================
0221 
0222 // The ResonanceTop class handles the top/antitop resonance.
0223 
0224 class ResonanceTop : public ResonanceWidths {
0225 
0226 public:
0227 
0228   // Constructor.
0229   ResonanceTop(int idResIn) : thetaWRat(), m2W(), tanBeta(), tan2Beta(),
0230     mbRun() {initBasic(idResIn);}
0231 
0232 private:
0233 
0234   // Locally stored properties and couplings.
0235   double thetaWRat, m2W, tanBeta, tan2Beta, mbRun;
0236 
0237   // Initialize constants.
0238   virtual void initConstants() override;
0239 
0240   // Calculate various common prefactors for the current mass.
0241   virtual void calcPreFac(bool = false) override;
0242 
0243   // Caclulate width for currently considered channel.
0244   virtual void calcWidth(bool = false) override;
0245 
0246 };
0247 
0248 //==========================================================================
0249 
0250 // The ResonanceFour class handles fourth-generation resonances.
0251 
0252 class ResonanceFour : public ResonanceWidths {
0253 
0254 public:
0255 
0256   // Constructor.
0257   ResonanceFour(int idResIn) : thetaWRat(), m2W() {initBasic(idResIn);}
0258 
0259 private:
0260 
0261   // Locally stored properties and couplings.
0262   double thetaWRat, m2W;
0263 
0264   // Initialize constants.
0265   virtual void initConstants() override;
0266 
0267   // Calculate various common prefactors for the current mass.
0268   virtual void calcPreFac(bool = false) override;
0269 
0270   // Caclulate width for currently considered channel.
0271   virtual void calcWidth(bool = false) override;
0272 
0273 };
0274 
0275 //==========================================================================
0276 
0277 // The ResonanceH class handles the SM and BSM Higgs resonance.
0278 // higgsType = 0 : SM H; = 1: h^0/H_1; = 2 : H^0/H_2; = 3 : A^0/A_3.
0279 
0280 class ResonanceH : public ResonanceWidths {
0281 
0282 public:
0283 
0284   // Constructor.
0285   ResonanceH(int higgsTypeIn, int idResIn) : higgsType(higgsTypeIn),
0286     useCubicWidth(), useRunLoopMass(), useNLOWidths(), sin2tW(), cos2tW(),
0287     mT(), mZ(), mW(), mHchg(), GammaT(), GammaZ(), GammaW(), rescAlpS(),
0288     rescColQ(), coup2d(), coup2u(), coup2l(), coup2Z(), coup2W(), coup2Hchg(),
0289     coup2H1H1(), coup2A3A3(), coup2H1Z(), coup2A3Z(), coup2A3H1(),
0290     coup2HchgW(), mLowT(), mStepT(), mLowZ(), mStepZ(), mLowW(), mStepW(),
0291     kinFacT(), kinFacZ(), kinFacW() {initBasic(idResIn);}
0292 
0293 private:
0294 
0295   // Constants: could only be changed in the code itself.
0296   static const double MASSMINWZ, MASSMINT, GAMMAMARGIN;
0297 
0298   // Higgs type in current instance.
0299   int    higgsType;
0300 
0301   // Locally stored properties and couplings.
0302   bool   useCubicWidth, useRunLoopMass, useNLOWidths;
0303   double sin2tW, cos2tW, mT, mZ, mW, mHchg, GammaT, GammaZ, GammaW,
0304          rescAlpS, rescColQ, coup2d, coup2u, coup2l, coup2Z, coup2W,
0305          coup2Hchg, coup2H1H1, coup2A3A3, coup2H1Z, coup2A3Z, coup2A3H1,
0306          coup2HchgW, mLowT, mStepT, mLowZ, mStepZ, mLowW, mStepW,
0307          kinFacT[101], kinFacZ[101], kinFacW[101];
0308 
0309   // Initialize constants.
0310   virtual void initConstants() override;
0311 
0312   // Calculate various common prefactors for the current mass.
0313   virtual void calcPreFac(bool = false) override;
0314 
0315   // Caclulate width for currently considered channel.
0316   virtual void calcWidth(bool = false) override;
0317 
0318   // Sum up loop contributions in Higgs -> g + g.
0319   double eta2gg();
0320 
0321   // Sum up loop contributions in Higgs -> gamma + gamma.
0322   double eta2gaga();
0323 
0324   // Sum up loop contributions in Higgs -> gamma + Z0.
0325   double eta2gaZ();
0326 
0327 };
0328 
0329 //==========================================================================
0330 
0331 // The ResonanceHchg class handles the H+- resonance.
0332 
0333 class ResonanceHchg : public ResonanceWidths {
0334 
0335 public:
0336 
0337   // Constructor.
0338   ResonanceHchg(int idResIn) : useCubicWidth(), thetaWRat(), mW(), tanBeta(),
0339     tan2Beta(), coup2H1W() {initBasic(idResIn);}
0340 
0341 private:
0342 
0343   // Locally stored properties and couplings.
0344   bool   useCubicWidth;
0345   double thetaWRat, mW, tanBeta, tan2Beta, coup2H1W;
0346 
0347   // Initialize constants.
0348   virtual void initConstants() override;
0349 
0350   // Calculate various common prefactors for the current mass.
0351   virtual void calcPreFac(bool = false) override;
0352 
0353   // Caclulate width for currently considered channel.
0354   virtual void calcWidth(bool = false) override;
0355 
0356 };
0357 
0358 //==========================================================================
0359 
0360 // The ResonanceZprime class handles the gamma*/Z0 /Z'^0 resonance.
0361 
0362 class ResonanceZprime : public ResonanceWidths {
0363 
0364 public:
0365 
0366   // Constructor.
0367   ResonanceZprime(int idResIn) : gmZmode(), maxZpGen(), sin2tW(), cos2tW(),
0368     thetaWRat(), mZ(), GammaZ(), m2Z(), GamMRatZ(), afZp(), vfZp(), coupZpWW(),
0369     ei2(), eivi(), vai2(), eivpi(), vaivapi(), vapi2(), gamNorm(), gamZNorm(),
0370     ZNorm(), gamZpNorm(), ZZpNorm(), ZpNorm() {initBasic(idResIn);}
0371 
0372 private:
0373 
0374   // Locally stored properties and couplings.
0375   int    gmZmode, maxZpGen;
0376   double sin2tW, cos2tW, thetaWRat, mZ, GammaZ, m2Z, GamMRatZ, afZp[20],
0377          vfZp[20], coupZpWW, ei2, eivi, vai2, eivpi, vaivapi, vapi2,
0378          gamNorm, gamZNorm, ZNorm, gamZpNorm, ZZpNorm, ZpNorm;
0379 
0380   // Initialize constants.
0381   virtual void initConstants() override;
0382 
0383   // Calculate various common prefactors for the current mass.
0384   virtual void calcPreFac(bool = false) override;
0385 
0386   // Caclulate width for currently considered channel.
0387   virtual void calcWidth(bool calledFromInit = false) override;
0388 
0389 };
0390 
0391 //==========================================================================
0392 
0393 // The ResonanceWprime class handles the W'+- resonance.
0394 
0395 class ResonanceWprime : public ResonanceWidths {
0396 
0397 public:
0398 
0399   // Constructor.
0400  ResonanceWprime(int idResIn) : ResonanceWidths(), thetaWRat(), cos2tW(),
0401     aqWp(), vqWp(), alWp(), vlWp(), coupWpWZ() {initBasic(idResIn);}
0402 
0403 private:
0404 
0405   // Locally stored properties and couplings.
0406   double thetaWRat, cos2tW, aqWp, vqWp, alWp, vlWp, coupWpWZ;
0407 
0408   // Initialize constants.
0409   virtual void initConstants() override;
0410 
0411   // Calculate various common prefactors for the current mass.
0412   virtual void calcPreFac(bool = false) override;
0413 
0414   // Caclulate width for currently considered channel.
0415   virtual void calcWidth(bool = false) override;
0416 
0417 };
0418 
0419 //==========================================================================
0420 
0421 // The ResonanceRhorizontal class handles the R^0 resonance.
0422 
0423 class ResonanceRhorizontal : public ResonanceWidths {
0424 
0425 public:
0426 
0427   // Constructor.
0428   ResonanceRhorizontal(int idResIn) : thetaWRat() {initBasic(idResIn);}
0429 
0430 private:
0431 
0432   // Locally stored properties and couplings.
0433   double thetaWRat;
0434 
0435   // Initialize constants.
0436   virtual void initConstants() override;
0437 
0438   // Calculate various common prefactors for the current mass.
0439   virtual void calcPreFac(bool = false) override;
0440 
0441   // Caclulate width for currently considered channel.
0442   virtual void calcWidth(bool = false) override;
0443 
0444 };
0445 
0446 //==========================================================================
0447 
0448 // The ResonanceExcited class handles excited-fermion resonances.
0449 
0450 class ResonanceExcited : public ResonanceWidths {
0451 
0452 public:
0453 
0454   // Constructor.
0455   ResonanceExcited(int idResIn) : Lambda(), coupF(), coupFprime(), coupFcol(),
0456     contactDec(), sin2tW(), cos2tW() {initBasic(idResIn);}
0457 
0458 private:
0459 
0460   // Locally stored properties and couplings.
0461   double Lambda, coupF, coupFprime, coupFcol, contactDec, sin2tW, cos2tW;
0462 
0463   // Initialize constants.
0464   virtual void initConstants() override;
0465 
0466   // Calculate various common prefactors for the current mass.
0467   virtual void calcPreFac(bool = false) override;
0468 
0469   // Caclulate width for currently considered channel.
0470   virtual void calcWidth(bool = false) override;
0471 
0472 };
0473 
0474 //==========================================================================
0475 
0476 // The ResonanceGraviton class handles the excited Graviton resonance.
0477 
0478 class ResonanceGraviton : public ResonanceWidths {
0479 
0480 public:
0481 
0482   // Constructor.
0483   ResonanceGraviton(int idResIn) : eDsmbulk(), eDvlvl(), kappaMG(),
0484     eDcoupling() {initBasic(idResIn);}
0485 
0486 private:
0487 
0488   // Locally stored properties and couplings.
0489   bool   eDsmbulk, eDvlvl;
0490   double kappaMG;
0491 
0492   // Couplings between graviton and SM (map from particle id to coupling).
0493   double eDcoupling[27];
0494 
0495   // Initialize constants.
0496   virtual void initConstants() override;
0497 
0498   // Calculate various common prefactors for the current mass.
0499   virtual void calcPreFac(bool = false) override;
0500 
0501   // Caclulate width for currently considered channel.
0502   virtual void calcWidth(bool = false) override;
0503 
0504 };
0505 
0506 //==========================================================================
0507 
0508 // The ResonanceKKgluon class handles the g^*/KK-gluon^* resonance.
0509 
0510 class ResonanceKKgluon : public ResonanceWidths {
0511 
0512 public:
0513 
0514   // Constructor.
0515   ResonanceKKgluon(int idResIn) : normSM(), normInt(), normKK(), eDgv(),
0516     eDga(), interfMode() {initBasic(idResIn);}
0517 
0518 private:
0519 
0520   // Locally stored properties.
0521   double normSM, normInt, normKK;
0522 
0523   // Couplings between kk gluon and SM (indexed by particle id).
0524   // Helicity dependent couplings. Use vector/axial-vector
0525   // couplings internally, gv/ga = 0.5 * (gL +/- gR).
0526   double eDgv[10], eDga[10];
0527 
0528   // Interference parameter.
0529   int interfMode;
0530 
0531   // Initialize constants.
0532   virtual void initConstants() override;
0533 
0534   // Calculate various common prefactors for the current mass.
0535   virtual void calcPreFac(bool calledFromInit = false) override;
0536 
0537   // Caclulate width for currently considered channel.
0538   virtual void calcWidth(bool calledFromInit = false) override;
0539 
0540 };
0541 
0542 //==========================================================================
0543 
0544 // The ResonanceLeptoquark class handles the LQ/LQbar resonance.
0545 
0546 class ResonanceLeptoquark : public ResonanceWidths {
0547 
0548 public:
0549 
0550   // Constructor.
0551   ResonanceLeptoquark(int idResIn) : kCoup() {initBasic(idResIn);}
0552 
0553 private:
0554 
0555   // Locally stored properties and couplings.
0556   double kCoup;
0557 
0558   // Initialize constants.
0559   virtual void initConstants() override;
0560 
0561   // Calculate various common prefactors for the current mass.
0562   virtual void calcPreFac(bool = false) override;
0563 
0564   // Caclulate width for currently considered channel.
0565   virtual void calcWidth(bool = false) override;
0566 
0567 };
0568 
0569 //==========================================================================
0570 
0571 // The ResonanceNuRight class handles righthanded Majorana neutrinos.
0572 
0573 class ResonanceNuRight : public ResonanceWidths {
0574 
0575 public:
0576 
0577   // Constructor.
0578   ResonanceNuRight(int idResIn) : thetaWRat(), mWR() {initBasic(idResIn);}
0579 
0580 private:
0581 
0582   // Locally stored properties and couplings.
0583   double thetaWRat, mWR;
0584 
0585   // Initialize constants.
0586   virtual void initConstants() override;
0587 
0588   // Calculate various common prefactors for the current mass.
0589   virtual void calcPreFac(bool = false) override;
0590 
0591   // Caclulate width for currently considered channel.
0592   virtual void calcWidth(bool = false) override;
0593 
0594 };
0595 
0596 //==========================================================================
0597 
0598 // The ResonanceZRight class handles the Z_R^0 resonance.
0599 
0600 class ResonanceZRight : public ResonanceWidths {
0601 
0602 public:
0603 
0604   // Constructor.
0605   ResonanceZRight(int idResIn) : sin2tW(), thetaWRat() {initBasic(idResIn);}
0606 
0607 private:
0608 
0609   // Locally stored properties and couplings.
0610   double sin2tW, thetaWRat;
0611 
0612   // Initialize constants.
0613   virtual void initConstants() override;
0614 
0615   // Calculate various common prefactors for the current mass.
0616   virtual void calcPreFac(bool = false) override;
0617 
0618   // Caclulate width for currently considered channel.
0619   virtual void calcWidth(bool = false) override;
0620 
0621 };
0622 
0623 //==========================================================================
0624 
0625 // The ResonanceWRight class handles the W_R+- resonance.
0626 
0627 class ResonanceWRight : public ResonanceWidths {
0628 
0629 public:
0630 
0631   // Constructor.
0632   ResonanceWRight(int idResIn) : thetaWRat() {initBasic(idResIn);}
0633 
0634 private:
0635 
0636   // Locally stored properties and couplings.
0637   double thetaWRat;
0638 
0639   // Initialize constants.
0640   virtual void initConstants() override;
0641 
0642   // Calculate various common prefactors for the current mass.
0643   virtual void calcPreFac(bool = false) override;
0644 
0645   // Caclulate width for currently considered channel.
0646   virtual void calcWidth(bool = false) override;
0647 
0648 };
0649 
0650 //==========================================================================
0651 
0652 // The ResonanceHchgchgLeft class handles the H++/H-- (left) resonance.
0653 
0654 class ResonanceHchgchgLeft : public ResonanceWidths {
0655 
0656 public:
0657 
0658   // Constructor.
0659   ResonanceHchgchgLeft(int idResIn) : yukawa(), gL(), vL(),
0660     mW() {initBasic(idResIn);}
0661 
0662 private:
0663 
0664   // Locally stored properties and couplings.
0665   double yukawa[4][4], gL, vL, mW;
0666 
0667   // Initialize constants.
0668   virtual void initConstants() override;
0669 
0670   // Calculate various common prefactors for the current mass.
0671   virtual void calcPreFac(bool = false) override;
0672 
0673   // Caclulate width for currently considered channel.
0674   virtual void calcWidth(bool = false) override;
0675 
0676 };
0677 
0678 //==========================================================================
0679 
0680 // The ResonanceHchgchgRight class handles the H++/H-- (right) resonance.
0681 
0682 class ResonanceHchgchgRight : public ResonanceWidths {
0683 
0684 public:
0685 
0686   // Constructor.
0687   ResonanceHchgchgRight(int idResIn) : idWR(), yukawa(),
0688     gR() {initBasic(idResIn);}
0689 
0690 private:
0691 
0692   // Locally stored properties and couplings.
0693   int    idWR;
0694   double yukawa[4][4], gR;
0695 
0696   // Initialize constants.
0697   virtual void initConstants() override;
0698 
0699   // Calculate various common prefactors for the current mass.
0700   virtual void calcPreFac(bool = false) override;
0701 
0702   // Caclulate width for currently considered channel.
0703   virtual void calcWidth(bool = false) override;
0704 
0705 };
0706 
0707 //==========================================================================
0708 
0709 } // end namespace Pythia8
0710 
0711 #endif // Pythia8_ResonanceWidths_H