source: HiSusy/trunk/Pythia8/pythia8170/include/ResonanceWidths.h @ 1

Last change on this file since 1 was 1, checked in by zerwas, 12 years ago

first import of structure, PYTHIA8 and DELPHES

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