StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ResonanceWidths.h
1 // ResonanceWidths.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 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 resonance properties: dynamical widths etc.
7 // ResonanceWidths: base class for all resonances.
8 // ResonanceGmZ, ...: derived classes for individual resonances.
9 
10 #ifndef Pythia8_ResonanceWidths_H
11 #define Pythia8_ResonanceWidths_H
12 
13 #include "Pythia8/Basics.h"
14 #include "Pythia8/Info.h"
15 #include "Pythia8/PythiaStdlib.h"
16 #include "Pythia8/Settings.h"
17 #include "Pythia8/StandardModel.h"
18 
19 namespace Pythia8 {
20 
21 //==========================================================================
22 
23 // Forward references to ParticleData and StandardModel classes.
24 class DecayChannel;
25 class ParticleData;
26 class ParticleDataEntry;
27 class CoupSM;
28 class CoupSUSY;
29 
30 //==========================================================================
31 
32 // The ResonanceWidths is the base class. Also used for generic resonaces.
33 
34 class ResonanceWidths {
35 
36 public:
37 
38  // Destructor.
39  virtual ~ResonanceWidths() {}
40 
41  // Set up standard properties.
42  void initBasic(int idResIn, bool isGenericIn = false) {
43  idRes = idResIn; isGeneric = isGenericIn;}
44 
45  // Calculate and store partial and total widths at the nominal mass.
46  virtual bool init(Info* infoPtrIn);
47 
48  // Return identity of particle species.
49  int id() const {return idRes;}
50 
51  // Calculate the total/open width for given mass, charge and instate.
52  double width(int idSgn, double mHatIn, int idInFlavIn = 0,
53  bool openOnly = false, bool setBR = false, int idOutFlav1 = 0,
54  int idOutFlav2 = 0);
55 
56  // Special case to calculate open final-state width.
57  double widthOpen(int idSgn, double mHatIn, int idIn = 0) {
58  return width( idSgn, mHatIn, idIn, true, false);}
59 
60  // Special case to store open final-state widths for channel selection.
61  double widthStore(int idSgn, double mHatIn, int idIn = 0) {
62  return width( idSgn, mHatIn, idIn, true, true);}
63 
64  // Return fraction of width open for particle and antiparticle.
65  double openFrac(int idSgn) {return (idSgn > 0) ? openPos : openNeg;}
66 
67  // Return forced rescaling factor of resonance width.
68  double widthRescaleFactor() {return forceFactor;}
69 
70  // Special case to calculate one final-state width.
71  // Currently only used for Higgs -> qqbar, g g or gamma gamma.
72  double widthChan(double mHatIn, int idOutFlav1, int idOutFlav2) {
73  return width( 1, mHatIn, 0, false, false, idOutFlav1, idOutFlav2);}
74 
75 protected:
76 
77  // Constructor.
78  ResonanceWidths() : idRes(), hasAntiRes(), doForceWidth(), isGeneric(),
79  allowCalcWidth(), minWidth(), minThreshold(), mRes(), GammaRes(), m2Res(),
80  GamMRat(), openPos(), openNeg(), forceFactor(), iChannel(), onMode(),
81  meMode(), mult(), id1(), id2(), id3(), id1Abs(), id2Abs(), id3Abs(),
82  idInFlav(), widNow(), mHat(), mf1(), mf2(), mf3(), mr1(), mr2(), mr3(),
83  ps(), kinFac(), alpEM(), alpS(), colQ(), preFac(), particlePtr(),
84  infoPtr(), settingsPtr(), particleDataPtr(), coupSMPtr(),
85  coupSUSYPtr(){}
86 
87  // Constants: could only be changed in the code itself.
88  static const int NPOINT;
89  static const double MASSMIN, MASSMARGIN;
90 
91  // Particle properties always present.
92  int idRes, hasAntiRes;
93  bool doForceWidth, isGeneric, allowCalcWidth;
94  double minWidth, minThreshold, mRes, GammaRes, m2Res, GamMRat,
95  openPos, openNeg, forceFactor;
96 
97  // Properties for currently studied decay channel(s).
98  int iChannel, onMode, meMode, mult, id1, id2, id3, id1Abs,
99  id2Abs, id3Abs, idInFlav;
100  double widNow, mHat, mf1, mf2, mf3, mr1, mr2, mr3, ps, kinFac,
101  alpEM, alpS, colQ, preFac;
102 
103  // Pointer to properties of the particle species.
104  ParticleDataEntry* particlePtr;
105 
106  // Pointer to various information on the generation.
107  Info* infoPtr;
108 
109  // Pointer to the settings database.
110  Settings* settingsPtr;
111 
112  // Pointer to the particle data table.
113  ParticleData* particleDataPtr;
114 
115  // Pointers to Standard Model and SUSY couplings.
116  CoupSM* coupSMPtr;
117  CoupSUSY* coupSUSYPtr;
118 
119  // Initialize constants.
120  virtual void initConstants() {}
121 
122  // Virtual methods to handle model-specific (non-SM) part of initialization
123  // for use by derived classes that implement additional models (eg SUSY).
124  virtual bool initBSM() {return true;}
125  virtual bool allowCalc() {return true;}
126 
127  // Calculate various common prefactors for the current mass.
128  // Optional argument calledFromInit only used for Z0.
129  virtual void calcPreFac(bool = false) {}
130 
131  // Calculate width for currently considered channel.
132  // Optional argument calledFromInit only used for Z0.
133  virtual void calcWidth(bool = false) {}
134 
135  // Simple routines for matrix-element integration over Breit-Wigners.
136  double numInt1BW(double mHatIn, double m1, double Gamma1, double mMin1,
137  double m2, int psMode = 1);
138  double numInt2BW(double mHatIn, double m1, double Gamma1, double mMin1,
139  double m2, double Gamma2, double mMin2, int psMode = 1);
140 
141 };
142 
143 //==========================================================================
144 
145 // The ResonanceGeneric class handles a generic resonance.
146 // Only needs a constructor and allowCalc = false; for the rest uses
147 // defaults in base class.
148 
149 class ResonanceGeneric : public ResonanceWidths {
150 
151 public:
152 
153  // Constructor.
154  ResonanceGeneric(int idResIn) {initBasic(idResIn, true);}
155 
156  // By default, assume no dedicated code exists to compute width.
157  virtual bool allowCalc() override {return false;}
158 
159 };
160 
161 //==========================================================================
162 
163 // The ResonanceGmZ class handles the gamma*/Z0 resonance.
164 
165 class ResonanceGmZ : public ResonanceWidths {
166 
167 public:
168 
169  // Constructor.
170  ResonanceGmZ(int idResIn) : gmZmode(), thetaWRat(), ei2(), eivi(), vi2ai2(),
171  gamNorm(), intNorm(), resNorm() {initBasic(idResIn);}
172 
173 private:
174 
175  // Locally stored properties and couplings.
176  int gmZmode;
177  double thetaWRat, ei2, eivi, vi2ai2, gamNorm, intNorm, resNorm;
178 
179  // Initialize constants.
180  virtual void initConstants() override;
181 
182  // Calculate various common prefactors for the current mass.
183  virtual void calcPreFac(bool = false) override;
184 
185  // Caclulate width for currently considered channel.
186  virtual void calcWidth(bool calledFromInit = false) override;
187 
188 };
189 
190 //==========================================================================
191 
192 // The ResonanceW class handles the W+- resonance.
193 
194 class ResonanceW : public ResonanceWidths {
195 
196 public:
197 
198  // Constructor.
199  ResonanceW(int idResIn) : ResonanceWidths(), thetaWRat() {initBasic(idResIn);}
200 
201 private:
202 
203  // Locally stored properties and couplings.
204  double thetaWRat;
205 
206  // Initialize constants.
207  virtual void initConstants() override;
208 
209  // Calculate various common prefactors for the current mass.
210  virtual void calcPreFac(bool = false) override;
211 
212  // Caclulate width for currently considered channel.
213  virtual void calcWidth(bool = false) override;
214 
215 };
216 
217 //==========================================================================
218 
219 // The ResonanceTop class handles the top/antitop resonance.
220 
221 class ResonanceTop : public ResonanceWidths {
222 
223 public:
224 
225  // Constructor.
226  ResonanceTop(int idResIn) : thetaWRat(), m2W(), tanBeta(), tan2Beta(),
227  mbRun() {initBasic(idResIn);}
228 
229 private:
230 
231  // Locally stored properties and couplings.
232  double thetaWRat, m2W, tanBeta, tan2Beta, mbRun;
233 
234  // Initialize constants.
235  virtual void initConstants() override;
236 
237  // Calculate various common prefactors for the current mass.
238  virtual void calcPreFac(bool = false) override;
239 
240  // Caclulate width for currently considered channel.
241  virtual void calcWidth(bool = false) override;
242 
243 };
244 
245 //==========================================================================
246 
247 // The ResonanceFour class handles fourth-generation resonances.
248 
249 class ResonanceFour : public ResonanceWidths {
250 
251 public:
252 
253  // Constructor.
254  ResonanceFour(int idResIn) : thetaWRat(), m2W() {initBasic(idResIn);}
255 
256 private:
257 
258  // Locally stored properties and couplings.
259  double thetaWRat, m2W;
260 
261  // Initialize constants.
262  virtual void initConstants() override;
263 
264  // Calculate various common prefactors for the current mass.
265  virtual void calcPreFac(bool = false) override;
266 
267  // Caclulate width for currently considered channel.
268  virtual void calcWidth(bool = false) override;
269 
270 };
271 
272 //==========================================================================
273 
274 // The ResonanceH class handles the SM and BSM Higgs resonance.
275 // higgsType = 0 : SM H; = 1: h^0/H_1; = 2 : H^0/H_2; = 3 : A^0/A_3.
276 
277 class ResonanceH : public ResonanceWidths {
278 
279 public:
280 
281  // Constructor.
282  ResonanceH(int higgsTypeIn, int idResIn) : higgsType(higgsTypeIn),
283  useCubicWidth(), useRunLoopMass(), useNLOWidths(), sin2tW(), cos2tW(),
284  mT(), mZ(), mW(), mHchg(), GammaT(), GammaZ(), GammaW(), rescAlpS(),
285  rescColQ(), coup2d(), coup2u(), coup2l(), coup2Z(), coup2W(), coup2Hchg(),
286  coup2H1H1(), coup2A3A3(), coup2H1Z(), coup2A3Z(), coup2A3H1(),
287  coup2HchgW(), mLowT(), mStepT(), mLowZ(), mStepZ(), mLowW(), mStepW(),
288  kinFacT(), kinFacZ(), kinFacW() {initBasic(idResIn);}
289 
290 private:
291 
292  // Constants: could only be changed in the code itself.
293  static const double MASSMINWZ, MASSMINT, GAMMAMARGIN;
294 
295  // Higgs type in current instance.
296  int higgsType;
297 
298  // Locally stored properties and couplings.
299  bool useCubicWidth, useRunLoopMass, useNLOWidths;
300  double sin2tW, cos2tW, mT, mZ, mW, mHchg, GammaT, GammaZ, GammaW,
301  rescAlpS, rescColQ, coup2d, coup2u, coup2l, coup2Z, coup2W,
302  coup2Hchg, coup2H1H1, coup2A3A3, coup2H1Z, coup2A3Z, coup2A3H1,
303  coup2HchgW, mLowT, mStepT, mLowZ, mStepZ, mLowW, mStepW,
304  kinFacT[101], kinFacZ[101], kinFacW[101];
305 
306  // Initialize constants.
307  virtual void initConstants() override;
308 
309  // Calculate various common prefactors for the current mass.
310  virtual void calcPreFac(bool = false) override;
311 
312  // Caclulate width for currently considered channel.
313  virtual void calcWidth(bool = false) override;
314 
315  // Sum up loop contributions in Higgs -> g + g.
316  double eta2gg();
317 
318  // Sum up loop contributions in Higgs -> gamma + gamma.
319  double eta2gaga();
320 
321  // Sum up loop contributions in Higgs -> gamma + Z0.
322  double eta2gaZ();
323 
324 };
325 
326 //==========================================================================
327 
328 // The ResonanceHchg class handles the H+- resonance.
329 
330 class ResonanceHchg : public ResonanceWidths {
331 
332 public:
333 
334  // Constructor.
335  ResonanceHchg(int idResIn) : useCubicWidth(), thetaWRat(), mW(), tanBeta(),
336  tan2Beta(), coup2H1W() {initBasic(idResIn);}
337 
338 private:
339 
340  // Locally stored properties and couplings.
341  bool useCubicWidth;
342  double thetaWRat, mW, tanBeta, tan2Beta, coup2H1W;
343 
344  // Initialize constants.
345  virtual void initConstants() override;
346 
347  // Calculate various common prefactors for the current mass.
348  virtual void calcPreFac(bool = false) override;
349 
350  // Caclulate width for currently considered channel.
351  virtual void calcWidth(bool = false) override;
352 
353 };
354 
355 //==========================================================================
356 
357 // The ResonanceZprime class handles the gamma*/Z0 /Z'^0 resonance.
358 
359 class ResonanceZprime : public ResonanceWidths {
360 
361 public:
362 
363  // Constructor.
364  ResonanceZprime(int idResIn) : gmZmode(), maxZpGen(), sin2tW(), cos2tW(),
365  thetaWRat(), mZ(), GammaZ(), m2Z(), GamMRatZ(), afZp(), vfZp(), coupZpWW(),
366  ei2(), eivi(), vai2(), eivpi(), vaivapi(), vapi2(), gamNorm(), gamZNorm(),
367  ZNorm(), gamZpNorm(), ZZpNorm(), ZpNorm() {initBasic(idResIn);}
368 
369 private:
370 
371  // Locally stored properties and couplings.
372  int gmZmode, maxZpGen;
373  double sin2tW, cos2tW, thetaWRat, mZ, GammaZ, m2Z, GamMRatZ, afZp[20],
374  vfZp[20], coupZpWW, ei2, eivi, vai2, eivpi, vaivapi, vapi2,
375  gamNorm, gamZNorm, ZNorm, gamZpNorm, ZZpNorm, ZpNorm;
376 
377  // Initialize constants.
378  virtual void initConstants() override;
379 
380  // Calculate various common prefactors for the current mass.
381  virtual void calcPreFac(bool = false) override;
382 
383  // Caclulate width for currently considered channel.
384  virtual void calcWidth(bool calledFromInit = false) override;
385 
386 };
387 
388 //==========================================================================
389 
390 // The ResonanceWprime class handles the W'+- resonance.
391 
392 class ResonanceWprime : public ResonanceWidths {
393 
394 public:
395 
396  // Constructor.
397  ResonanceWprime(int idResIn) : ResonanceWidths(), thetaWRat(), cos2tW(),
398  aqWp(), vqWp(), alWp(), vlWp(), coupWpWZ() {initBasic(idResIn);}
399 
400 private:
401 
402  // Locally stored properties and couplings.
403  double thetaWRat, cos2tW, aqWp, vqWp, alWp, vlWp, coupWpWZ;
404 
405  // Initialize constants.
406  virtual void initConstants() override;
407 
408  // Calculate various common prefactors for the current mass.
409  virtual void calcPreFac(bool = false) override;
410 
411  // Caclulate width for currently considered channel.
412  virtual void calcWidth(bool = false) override;
413 
414 };
415 
416 //==========================================================================
417 
418 // The ResonanceRhorizontal class handles the R^0 resonance.
419 
420 class ResonanceRhorizontal : public ResonanceWidths {
421 
422 public:
423 
424  // Constructor.
425  ResonanceRhorizontal(int idResIn) : thetaWRat() {initBasic(idResIn);}
426 
427 private:
428 
429  // Locally stored properties and couplings.
430  double thetaWRat;
431 
432  // Initialize constants.
433  virtual void initConstants() override;
434 
435  // Calculate various common prefactors for the current mass.
436  virtual void calcPreFac(bool = false) override;
437 
438  // Caclulate width for currently considered channel.
439  virtual void calcWidth(bool = false) override;
440 
441 };
442 
443 //==========================================================================
444 
445 // The ResonanceExcited class handles excited-fermion resonances.
446 
447 class ResonanceExcited : public ResonanceWidths {
448 
449 public:
450 
451  // Constructor.
452  ResonanceExcited(int idResIn) : Lambda(), coupF(), coupFprime(), coupFcol(),
453  contactDec(), sin2tW(), cos2tW() {initBasic(idResIn);}
454 
455 private:
456 
457  // Locally stored properties and couplings.
458  double Lambda, coupF, coupFprime, coupFcol, contactDec, sin2tW, cos2tW;
459 
460  // Initialize constants.
461  virtual void initConstants() override;
462 
463  // Calculate various common prefactors for the current mass.
464  virtual void calcPreFac(bool = false) override;
465 
466  // Caclulate width for currently considered channel.
467  virtual void calcWidth(bool = false) override;
468 
469 };
470 
471 //==========================================================================
472 
473 // The ResonanceGraviton class handles the excited Graviton resonance.
474 
475 class ResonanceGraviton : public ResonanceWidths {
476 
477 public:
478 
479  // Constructor.
480  ResonanceGraviton(int idResIn) : eDsmbulk(), eDvlvl(), kappaMG(),
481  eDcoupling() {initBasic(idResIn);}
482 
483 private:
484 
485  // Locally stored properties and couplings.
486  bool eDsmbulk, eDvlvl;
487  double kappaMG;
488 
489  // Couplings between graviton and SM (map from particle id to coupling).
490  double eDcoupling[27];
491 
492  // Initialize constants.
493  virtual void initConstants() override;
494 
495  // Calculate various common prefactors for the current mass.
496  virtual void calcPreFac(bool = false) override;
497 
498  // Caclulate width for currently considered channel.
499  virtual void calcWidth(bool = false) override;
500 
501 };
502 
503 //==========================================================================
504 
505 // The ResonanceKKgluon class handles the g^*/KK-gluon^* resonance.
506 
507 class ResonanceKKgluon : public ResonanceWidths {
508 
509 public:
510 
511  // Constructor.
512  ResonanceKKgluon(int idResIn) : normSM(), normInt(), normKK(), eDgv(),
513  eDga(), interfMode() {initBasic(idResIn);}
514 
515 private:
516 
517  // Locally stored properties.
518  double normSM, normInt, normKK;
519 
520  // Couplings between kk gluon and SM (indexed by particle id).
521  // Helicity dependent couplings. Use vector/axial-vector
522  // couplings internally, gv/ga = 0.5 * (gL +/- gR).
523  double eDgv[10], eDga[10];
524 
525  // Interference parameter.
526  int interfMode;
527 
528  // Initialize constants.
529  virtual void initConstants() override;
530 
531  // Calculate various common prefactors for the current mass.
532  virtual void calcPreFac(bool calledFromInit = false) override;
533 
534  // Caclulate width for currently considered channel.
535  virtual void calcWidth(bool calledFromInit = false) override;
536 
537 };
538 
539 //==========================================================================
540 
541 // The ResonanceLeptoquark class handles the LQ/LQbar resonance.
542 
543 class ResonanceLeptoquark : public ResonanceWidths {
544 
545 public:
546 
547  // Constructor.
548  ResonanceLeptoquark(int idResIn) : kCoup() {initBasic(idResIn);}
549 
550 private:
551 
552  // Locally stored properties and couplings.
553  double kCoup;
554 
555  // Initialize constants.
556  virtual void initConstants() override;
557 
558  // Calculate various common prefactors for the current mass.
559  virtual void calcPreFac(bool = false) override;
560 
561  // Caclulate width for currently considered channel.
562  virtual void calcWidth(bool = false) override;
563 
564 };
565 
566 //==========================================================================
567 
568 // The ResonanceNuRight class handles righthanded Majorana neutrinos.
569 
570 class ResonanceNuRight : public ResonanceWidths {
571 
572 public:
573 
574  // Constructor.
575  ResonanceNuRight(int idResIn) : thetaWRat(), mWR() {initBasic(idResIn);}
576 
577 private:
578 
579  // Locally stored properties and couplings.
580  double thetaWRat, mWR;
581 
582  // Initialize constants.
583  virtual void initConstants() override;
584 
585  // Calculate various common prefactors for the current mass.
586  virtual void calcPreFac(bool = false) override;
587 
588  // Caclulate width for currently considered channel.
589  virtual void calcWidth(bool = false) override;
590 
591 };
592 
593 //==========================================================================
594 
595 // The ResonanceZRight class handles the Z_R^0 resonance.
596 
597 class ResonanceZRight : public ResonanceWidths {
598 
599 public:
600 
601  // Constructor.
602  ResonanceZRight(int idResIn) : sin2tW(), thetaWRat() {initBasic(idResIn);}
603 
604 private:
605 
606  // Locally stored properties and couplings.
607  double sin2tW, thetaWRat;
608 
609  // Initialize constants.
610  virtual void initConstants() override;
611 
612  // Calculate various common prefactors for the current mass.
613  virtual void calcPreFac(bool = false) override;
614 
615  // Caclulate width for currently considered channel.
616  virtual void calcWidth(bool = false) override;
617 
618 };
619 
620 //==========================================================================
621 
622 // The ResonanceWRight class handles the W_R+- resonance.
623 
624 class ResonanceWRight : public ResonanceWidths {
625 
626 public:
627 
628  // Constructor.
629  ResonanceWRight(int idResIn) : thetaWRat() {initBasic(idResIn);}
630 
631 private:
632 
633  // Locally stored properties and couplings.
634  double thetaWRat;
635 
636  // Initialize constants.
637  virtual void initConstants() override;
638 
639  // Calculate various common prefactors for the current mass.
640  virtual void calcPreFac(bool = false) override;
641 
642  // Caclulate width for currently considered channel.
643  virtual void calcWidth(bool = false) override;
644 
645 };
646 
647 //==========================================================================
648 
649 // The ResonanceHchgchgLeft class handles the H++/H-- (left) resonance.
650 
651 class ResonanceHchgchgLeft : public ResonanceWidths {
652 
653 public:
654 
655  // Constructor.
656  ResonanceHchgchgLeft(int idResIn) : yukawa(), gL(), vL(),
657  mW() {initBasic(idResIn);}
658 
659 private:
660 
661  // Locally stored properties and couplings.
662  double yukawa[4][4], gL, vL, mW;
663 
664  // Initialize constants.
665  virtual void initConstants() override;
666 
667  // Calculate various common prefactors for the current mass.
668  virtual void calcPreFac(bool = false) override;
669 
670  // Caclulate width for currently considered channel.
671  virtual void calcWidth(bool = false) override;
672 
673 };
674 
675 //==========================================================================
676 
677 // The ResonanceHchgchgRight class handles the H++/H-- (right) resonance.
678 
679 class ResonanceHchgchgRight : public ResonanceWidths {
680 
681 public:
682 
683  // Constructor.
684  ResonanceHchgchgRight(int idResIn) : idWR(), yukawa(),
685  gR() {initBasic(idResIn);}
686 
687 private:
688 
689  // Locally stored properties and couplings.
690  int idWR;
691  double yukawa[4][4], gR;
692 
693  // Initialize constants.
694  virtual void initConstants() override;
695 
696  // Calculate various common prefactors for the current mass.
697  virtual void calcPreFac(bool = false) override;
698 
699  // Caclulate width for currently considered channel.
700  virtual void calcWidth(bool = false) override;
701 
702 };
703 
704 //==========================================================================
705 
706 } // end namespace Pythia8
707 
708 #endif // Pythia8_ResonanceWidths_H