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