StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
HelicityMatrixElements.h
1 // HelicityMatrixElements.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Philip Ilten, 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 a number of physics classes used in tau decays.
7 
8 #ifndef Pythia8_HelicityMatrixElements_H
9 #define Pythia8_HelicityMatrixElements_H
10 
11 #include "Pythia8/Basics.h"
12 #include "Pythia8/Event.h"
13 #include "Pythia8/HelicityBasics.h"
14 #include "Pythia8/PythiaComplex.h"
15 #include "Pythia8/PythiaStdlib.h"
16 #include "Pythia8/StandardModel.h"
17 
18 namespace Pythia8 {
19 
20 //==========================================================================
21 
22 // The helicity matrix element class.
23 
24 class HelicityMatrixElement {
25 
26 public:
27 
28  // Constructor and destructor.
29  HelicityMatrixElement() : DECAYWEIGHTMAX(), particleDataPtr(),
30  coupSMPtr(), settingsPtr() {};
31  virtual ~HelicityMatrixElement() {};
32 
33  // Initialize the physics matrices and pointers.
34  virtual void initPointers(ParticleData*, CoupSM* , Settings* = 0);
35 
36  // Initialize the channel.
37  virtual HelicityMatrixElement* initChannel(vector<HelicityParticle>&);
38 
39  // Calculate the matrix element weight for a decay.
40  virtual double decayWeight(vector<HelicityParticle>&);
41 
42  // Calculate the maximum matrix element decay weight.
43  virtual double decayWeightMax(vector<HelicityParticle>&)
44  {return DECAYWEIGHTMAX;}
45 
46  // Calculate the helicity matrix element.
47  virtual complex calculateME(vector<int>){return complex(0,0);}
48 
49  // Calculate the decay matrix for a particle.
50  virtual void calculateD(vector<HelicityParticle>&);
51 
52  // Calculate the density matrix for a particle.
53  virtual void calculateRho(unsigned int, vector<HelicityParticle>&);
54 
55  // Set a fermion line.
56  void setFermionLine(int, HelicityParticle&, HelicityParticle&);
57 
58  // Calculate Breit-Wigner's with running widths and fixed.
59  virtual complex breitWigner(double s, double M, double G);
60  virtual complex sBreitWigner(double m0, double m1, double s,
61  double M, double G);
62  virtual complex pBreitWigner(double m0, double m1, double s,
63  double M, double G);
64  virtual complex dBreitWigner(double m0, double m1, double s,
65  double M, double G);
66 
67 protected:
68 
69  // Maximum decay weight. WARNING: hardcoded constant.
70  double DECAYWEIGHTMAX;
71 
72  // Physics matrices.
73  vector< GammaMatrix > gamma;
74 
75  // Particle map vector.
76  vector< int > pMap;
77 
78  // Particle ID vector.
79  vector< int > pID;
80 
81  // Particle mass vector.
82  vector< double > pM;
83 
84  // Wave functions.
85  vector< vector< Wave4 > > u;
86 
87  // Initialize the constants for the matrix element (called by initChannel).
88  virtual void initConstants() {};
89 
90  // Initialize the wave functions (called by decayWeight and calculateRho/D).
91  virtual void initWaves(vector<HelicityParticle>&) {};
92 
93  // Pointer to particle data.
94  ParticleData* particleDataPtr;
95 
96  // Pointer to Standard Model constants.
97  CoupSM* coupSMPtr;
98 
99  // Pointer to Settings.
100  Settings* settingsPtr;
101 
102 private:
103 
104  // Recursive sub-method to calculate the density matrix for a particle.
105  void calculateRho(unsigned int, vector<HelicityParticle>&,
106  vector<int>&, vector<int>&, unsigned int);
107 
108  // Recursive sub-method to calculate the decay matrix for a particle.
109  void calculateD(vector<HelicityParticle>&, vector<int>&, vector<int>&,
110  unsigned int);
111 
112  // Recursive sub-method to calculate the matrix element weight for a decay.
113  void decayWeight(vector<HelicityParticle>&, vector<int>&, vector<int>&,
114  complex&, unsigned int);
115 
116  // Calculate the product of the decay matrices for a hard process.
117  complex calculateProductD(unsigned int, unsigned int,
118  vector<HelicityParticle>&, vector<int>&, vector<int>&);
119 
120  // Calculate the product of the decay matrices for a decay process.
121  complex calculateProductD(vector<HelicityParticle>&,
122  vector<int>&, vector<int>&);
123 
124 };
125 
126 //==========================================================================
127 
128 // Helicity matrix element for the hard process of two fermions -> W/W' ->
129 // two fermions.
130 
131 class HMETwoFermions2W2TwoFermions : public HelicityMatrixElement {
132 
133 public:
134 
135  HMETwoFermions2W2TwoFermions() : p0CA(), p2CA(), p0CV(), p2CV() {}
136 
137  void initConstants();
138 
139  void initWaves(vector<HelicityParticle>&);
140 
141  complex calculateME(vector<int>);
142 
143 private:
144 
145  // Vector and axial couplings.
146  double p0CA, p2CA, p0CV, p2CV;
147 
148 };
149 
150 //==========================================================================
151 
152 // Helicity matrix element for the hard process of two fermions ->
153 // photon/Z/Z' -> two fermions.
154 
155 class HMETwoFermions2GammaZ2TwoFermions : public HelicityMatrixElement {
156 
157 public:
158 
159  HMETwoFermions2GammaZ2TwoFermions() : p0CAZ(), p2CAZ(), p0CVZ(), p2CVZ(),
160  p0CAZp(), p2CAZp(), p0CVZp(), p2CVZp(), cos2W(), sin2W(), zG(), zM(),
161  zpG(), zpM(), s(), p0Q(), p2Q(), zaxis(), includeGamma(), includeZ(),
162  includeZp() {}
163 
164  void initConstants();
165 
166  void initWaves(vector<HelicityParticle>&);
167 
168  complex calculateME(vector<int>);
169 
170 private:
171 
172  // Return gamma element for the helicity matrix element.
173  complex calculateGammaME(vector<int>);
174 
175  // Return Z/Z' element for helicity matrix element.
176  complex calculateZME(vector<int>, double, double, double, double,
177  double, double);
178 
179  // Return the Z' vector or axial coupling for a fermion.
180  double zpCoupling(int id, string type);
181 
182  // Vector and axial couplings.
183  double p0CAZ, p2CAZ, p0CVZ, p2CVZ, p0CAZp, p2CAZp, p0CVZp, p2CVZp;
184 
185  // Weinberg angle, Z/Z' width (on-shell), Z/Z' mass (on-shell), and s.
186  double cos2W, sin2W, zG, zM, zpG, zpM, s;
187 
188  // Fermion line charge.
189  double p0Q, p2Q;
190 
191  // Bool whether the incoming fermions are oriented with the z-axis.
192  bool zaxis, includeGamma, includeZ, includeZp;
193 
194 };
195 
196 //==========================================================================
197 
198 // Helicity matrix element for the hard process of X -> two fermions.
199 
200 class HMEX2TwoFermions : public HelicityMatrixElement {
201 
202 public:
203 
204  void initWaves(vector<HelicityParticle>&);
205 
206 };
207 
208 //==========================================================================
209 
210 // Helicity matrix element for the hard process of W/W' -> two fermions.
211 
212 class HMEW2TwoFermions : public HMEX2TwoFermions {
213 
214 public:
215 
216  HMEW2TwoFermions() : p2CA(), p2CV() {}
217 
218  void initConstants();
219 
220  complex calculateME(vector<int>);
221 
222 private:
223 
224  // Vector and axial couplings.
225  double p2CA, p2CV;
226 
227 };
228 
229 //==========================================================================
230 
231 // Helicity matrix element for the hard process of photon -> two fermions.
232 
233 class HMEGamma2TwoFermions : public HMEX2TwoFermions {
234 
235 public:
236 
237  complex calculateME(vector<int>);
238 
239 };
240 
241 //==========================================================================
242 
243 // Helicity matrix element for the hard process of Z/Z' -> two fermions.
244 
245 class HMEZ2TwoFermions : public HMEX2TwoFermions {
246 
247 public:
248 
249  HMEZ2TwoFermions() : p2CA(), p2CV() {}
250 
251  void initConstants();
252 
253  complex calculateME(vector<int>);
254 
255 private:
256 
257  // Return the Z' vector or axial coupling for a fermion.
258  double zpCoupling(int id, string type);
259 
260  // Vector and axial couplings.
261  double p2CA, p2CV;
262 
263 };
264 
265 //==========================================================================
266 
267 // Helicity matrix element for the decay of a Higgs -> two fermions.
268 
269 // Because the Higgs is spin zero the Higgs production mechanism is not
270 // needed for calculating helicity density matrices. However, the CP mixing
271 // is needed.
272 
273 class HMEHiggs2TwoFermions : public HelicityMatrixElement {
274 
275 public:
276 
277  void initConstants();
278 
279  void initWaves(vector<HelicityParticle>&);
280 
281  complex calculateME(vector<int>);
282 
283 private:
284 
285  // Coupling constants of the fermions with the Higgs.
286  complex p2CA, p2CV;
287 
288 };
289 
290 //==========================================================================
291 
292 // Base class for all tau decay helicity matrix elements.
293 
294 class HMETauDecay : public HelicityMatrixElement {
295 
296 public:
297 
298  virtual void initWaves(vector<HelicityParticle>&);
299 
300  virtual complex calculateME(vector<int>);
301 
302  virtual double decayWeightMax(vector<HelicityParticle>&);
303 
304 protected:
305 
306  virtual void initHadronicCurrent(vector<HelicityParticle>&) {};
307 
308  virtual void calculateResonanceWeights(vector<double>&, vector<double>&,
309  vector<complex>&);
310 
311 };
312 
313 //==========================================================================
314 
315 // Helicity matrix element for a tau decaying into a single scalar meson.
316 
317 class HMETau2Meson : public HMETauDecay {
318 
319 public:
320 
321  void initConstants();
322 
323  void initHadronicCurrent(vector<HelicityParticle>&);
324 
325 };
326 
327 //==========================================================================
328 
329 // Helicity matrix element for a tau decaying into two leptons.
330 
331 class HMETau2TwoLeptons : public HMETauDecay {
332 
333 public:
334 
335  void initConstants();
336 
337  void initWaves(vector<HelicityParticle>&);
338 
339  complex calculateME(vector<int>);
340 
341 };
342 
343 //==========================================================================
344 
345 // Helicity matrix element for a tau decaying into two mesons through a
346 // vector meson resonance.
347 
348 class HMETau2TwoMesonsViaVector : public HMETauDecay {
349 
350 public:
351 
352  void initConstants();
353 
354  void initHadronicCurrent(vector<HelicityParticle>&);
355 
356 private:
357 
358  // Resonance masses, widths, and weights.
359  vector<double> vecM, vecG, vecP, vecA;
360  vector<complex> vecW;
361 
362 };
363 
364 //==========================================================================
365 
366 // Helicity matrix element for a tau decay into two mesons through a vector
367 // or scalar meson resonance.
368 
369 class HMETau2TwoMesonsViaVectorScalar : public HMETauDecay {
370 
371 public:
372 
373  HMETau2TwoMesonsViaVectorScalar() : scaC(), vecC() {}
374 
375  void initConstants();
376 
377  void initHadronicCurrent(vector<HelicityParticle>&);
378 
379 private:
380 
381  // Coupling to vector and scalar resonances.
382  double scaC, vecC;
383 
384  // Resonance masses, widths, and weights.
385  vector<double> scaM, scaG, scaP, scaA, vecM, vecG, vecP, vecA;
386  vector<complex> scaW, vecW;
387 
388 };
389 
390 //==========================================================================
391 
392 // Helicity matrix element for a tau decay into three mesons (base class).
393 
394 class HMETau2ThreeMesons : public HMETauDecay {
395 
396 public:
397 
398  HMETau2ThreeMesons() = default;
399 
400  void initConstants();
401 
402  void initHadronicCurrent(vector<HelicityParticle>&);
403 
404 protected:
405 
406  // Decay mode of the tau.
407  enum Mode{Pi0Pi0Pim, PimPimPip, Pi0PimK0b, PimPipKm, Pi0PimEta, PimKmKp,
408  Pi0K0Km, KlPimKs, Pi0Pi0Km, KlKlPim, PimKsKs, PimK0bK0, Uknown};
409  Mode mode{};
410 
411  // Initialize decay mode and resonance constants (called by initConstants).
412  virtual void initMode();
413  virtual void initResonances() {;}
414 
415  // Initialize the momenta.
416  virtual void initMomenta(vector<HelicityParticle>&);
417 
418  // Center of mass energies and momenta.
419  double s1{}, s2{}, s3{}, s4{};
420  Wave4 q{}, q2{}, q3{}, q4{};
421 
422  // Stored a1 Breit-Wigner (for speed).
423  complex a1BW{};
424 
425  // Form factors.
426  virtual complex F1() {return complex(0, 0);}
427  virtual complex F2() {return complex(0, 0);}
428  virtual complex F3() {return complex(0, 0);}
429  virtual complex F4() {return complex(0, 0);}
430 
431  // Phase space and Breit-Wigner for the a1.
432  virtual double a1PhaseSpace(double);
433  virtual complex a1BreitWigner(double);
434 
435  // Sum running p and fixed width Breit-Wigner resonances.
436  complex T(double m0, double m1, double s,
437  vector<double>& M, vector<double>& G, vector<double>& W);
438  complex T(double s, vector<double>& M, vector<double>& G, vector<double>& W);
439 
440 };
441 
442 //==========================================================================
443 
444 // Helicity matrix element for a tau decay into three pions.
445 
446 class HMETau2ThreePions : public HMETau2ThreeMesons {
447 
448 public:
449 
450  HMETau2ThreePions() : f0M(), f0G(), f0P(), f0A(), f2M(), f2G(), f2P(),
451  f2A(), sigM(), sigG(), sigP(), sigA() {}
452 
453 private:
454 
455  void initResonances();
456 
457  // Resonance masses, widths, and weights.
458  vector<double> rhoM, rhoG, rhoPp, rhoAp, rhoPd, rhoAd;
459  double f0M, f0G, f0P, f0A, f2M, f2G, f2P, f2A;
460  double sigM, sigG, sigP, sigA;
461  vector<complex> rhoWp, rhoWd;
462  complex f0W, f2W, sigW;
463 
464  // Form factors.
465  complex F1();
466  complex F2();
467  complex F3();
468 
469  // Running width and Breit-Wigner for the a1.
470  double a1PhaseSpace(double);
471  complex a1BreitWigner(double);
472 
473 };
474 
475 //==========================================================================
476 
477 // Helicity matrix element for a tau decay into three mesons with kaons.
478 
479 class HMETau2ThreeMesonsWithKaons : public HMETau2ThreeMesons {
480 
481 public:
482 
483  HMETau2ThreeMesonsWithKaons() : kM(), piM(), piW() {}
484 
485 private:
486 
487  void initResonances();
488 
489  // Resonance masses, widths, and weights.
490  vector<double> rhoMa, rhoGa, rhoWa, rhoMv, rhoGv, rhoWv;
491  vector<double> kstarMa, kstarGa, kstarWa, kstarMv, kstarGv, kstarWv;
492  vector<double> k1Ma, k1Ga, k1Wa, k1Mb, k1Gb, k1Wb;
493  vector<double> omegaM, omegaG, omegaW;
494  double kM, piM, piW;
495 
496  // Form factors.
497  complex F1();
498  complex F2();
499  complex F4();
500 
501 };
502 
503 //==========================================================================
504 
505 // Helicity matrix element for a tau decay into generic three mesons.
506 
507 class HMETau2ThreeMesonsGeneric : public HMETau2ThreeMesons {
508 
509 public:
510 
511  HMETau2ThreeMesonsGeneric() : kM(), piM(), piW() {}
512 
513 private:
514 
515  void initResonances();
516 
517  // Resonance masses, widths, and weights.
518  vector<double> rhoMa, rhoGa, rhoWa, rhoMv, rhoGv, rhoWv;
519  vector<double> kstarM, kstarG, kstarW, k1M, k1G, k1W;
520  double kM, piM, piW;
521 
522  // Form factors.
523  complex F1();
524  complex F2();
525  complex F4();
526 
527 };
528 
529 //==========================================================================
530 
531 // Helicity matrix element for a tau decay into two pions and a photon.
532 
533 class HMETau2TwoPionsGamma : public HMETauDecay {
534 
535 public:
536 
537  HMETau2TwoPionsGamma() : piM() {}
538 
539  void initConstants();
540 
541  void initWaves(vector<HelicityParticle>&);
542 
543  complex calculateME(vector<int>);
544 
545 protected:
546 
547  // Resonance masses, widths, and weights.
548  vector<double> rhoM, rhoG, rhoW, omegaM, omegaG, omegaW;
549  double piM;
550 
551  // Form factor.
552  complex F(double s, vector<double> M, vector<double> G, vector<double> W);
553 
554 };
555 
556 //==========================================================================
557 
558 // Helicity matrix element for a tau decay into four pions.
559 
560 class HMETau2FourPions : public HMETauDecay {
561 
562 public:
563 
564  HMETau2FourPions() : a1M(), a1G(), rhoM(), rhoG(), sigM(), sigG(), omeM(),
565  omeG(), picM(), pinM(), sigA(), sigP(), omeA(), omeP(), lambda2() {}
566 
567  void initConstants();
568 
569  void initHadronicCurrent(vector<HelicityParticle>& p);
570 
571 private:
572 
573  // G-function form factors (fits).
574  double G(int i, double s);
575 
576  // T-vector functions.
577  Wave4 t1(Wave4&, Wave4&, Wave4&, Wave4&, Wave4&);
578  Wave4 t2(Wave4&, Wave4&, Wave4&, Wave4&, Wave4&);
579  Wave4 t3(Wave4&, Wave4&, Wave4&, Wave4&, Wave4&);
580 
581  // Breit-Wigner denominators for the intermediate mesons.
582  complex a1D(double s);
583  complex rhoD(double s);
584  complex sigD(double s);
585  complex omeD(double s);
586 
587  // Form factors needed for the a1, rho, and omega.
588  double a1FormFactor(double s);
589  double rhoFormFactor1(double s);
590  double rhoFormFactor2(double s);
591  double omeFormFactor(double s);
592 
593  // Masses and widths of the intermediate mesons.
594  double a1M, a1G, rhoM, rhoG, sigM, sigG, omeM, omeG;
595 
596  // Masses for the pions (charged and neutral).
597  double picM, pinM;
598 
599  // Amplitude, phases, and weights for mixing.
600  double sigA, sigP, omeA, omeP;
601  complex sigW, omeW;
602 
603  // Cut-off for a1 form factor.
604  double lambda2;
605 
606 };
607 
608 //==========================================================================
609 
610 // Helicity matrix element for a tau decaying into five pions.
611 
612 class HMETau2FivePions : public HMETauDecay {
613 
614 public:
615 
616  HMETau2FivePions() : a1M(), a1G(), rhoM(), rhoG(), omegaM(), omegaG(),
617  omegaW(), sigmaM(), sigmaG(), sigmaW() {}
618 
619  void initConstants();
620 
621  void initHadronicCurrent(vector<HelicityParticle>&);
622 
623 private:
624 
625  // Hadronic currents.
626  Wave4 Ja(Wave4 &q, Wave4 &q1, Wave4 &q2, Wave4 &q3, Wave4 &q4, Wave4 &q5);
627  Wave4 Jb(Wave4 &q, Wave4 &q1, Wave4 &q2, Wave4 &q3, Wave4 &q4, Wave4 &q5);
628 
629  // Simplified s-wave Breit-Wigner assuming massless products.
630  complex breitWigner(double s, double M, double G);
631 
632  // Masses and widths of intermediates.
633  double a1M, a1G, rhoM, rhoG, omegaM, omegaG, omegaW, sigmaM, sigmaG, sigmaW;
634 
635 };
636 
637 //==========================================================================
638 
639 // Helicity matrix element for a tau decay into flat phase space.
640 
641 class HMETau2PhaseSpace : public HMETauDecay {
642 
643 public:
644 
645  void initWaves(vector<HelicityParticle>&) {};
646 
647  complex calculateME(vector<int>) {return 1;}
648 
649  void calculateD(vector<HelicityParticle>&) {};
650 
651  void calculateRho(unsigned int, vector<HelicityParticle>&) {};
652 
653  double decayWeight(vector<HelicityParticle>&) {return 1.0;}
654 
655  double decayWeightMax(vector<HelicityParticle>&) {return 1.0;}
656 
657 };
658 
659 //==========================================================================
660 
661 } // end namespace Pythia8
662 
663 #endif // end Pythia8_HelicityMatrixElements_H