a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
LEP2TwoFermions Class Reference

A class for cross sections and forward-backward asymmetries of \(e^+e^-\to f \bar{f}\) above the \(Z\) pole. More...

#include <LEP2TwoFermions.h>

Detailed Description

A class for cross sections and forward-backward asymmetries of \(e^+e^-\to f \bar{f}\) above the \(Z\) pole.

Author
HEPfit Collaboration

Definition at line 21 of file LEP2TwoFermions.h.

Public Member Functions

double AFB_l (const QCD::lepton l, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 
double AFB_q (const QCD::quark q, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 
double dsigma_l (const QCD::lepton l, const double mf, const double s, const double cosTheta, const double Mw, const double GammaZ, const bool bWeak) const
 
double dsigma_l_box (const QCD::lepton l, const double mf, const double s, const double cosTheta, const double Mw, const double GammaZ) const
 
double dsigma_q (const QCD::quark q, const double mf, const double s, const double cosTheta, const double Mw, const double GammaZ, const bool bWeak) const
 
double dsigma_q_box (const QCD::quark q, const double mf, const double s, const double cosTheta, const double Mw, const double GammaZ) const
 
double G_3prime_l (const QCD::lepton l, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 
double G_3prime_q (const QCD::quark q, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 
double H_ISR (const double x, const double s) const
 
double H_ISR_FB (const double x, const double s) const
 
 LEP2TwoFermions (const StandardModel &SM_i)
 LEP2TwoFermions constructor. More...
 
double QCD_FSR_forAFB (const QCD::quark q, const double mf, const double s) const
 
double QCD_FSR_forSigma (const double s) const
 
double QED_FSR_forSigma (const double s, const double Qf) const
 
double sigma_l (const QCD::lepton l, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 
double sigma_q (const QCD::quark q, const double mf, const double s, const double Mw, const double GammaZ, const bool bWeak) const
 

Private Member Functions

double AFB (const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWeak) const
 
double alpha_at_s (const double s) const
 
double dsigma (const double s, const double cosTheta, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf, const bool bWeak) const
 
double dsigma_box (const double s, const double cosTheta, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf) const
 
double sigma (const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf, const bool bWeak) const
 

Private Attributes

const StandardModelSM
 

Constructor & Destructor Documentation

◆ LEP2TwoFermions()

LEP2TwoFermions::LEP2TwoFermions ( const StandardModel SM_i)

LEP2TwoFermions constructor.

Parameters
[in]SM_ian object of StandardModel class

Definition at line 12 of file LEP2TwoFermions.cpp.

13: SM(SM_i)
14{
15}
const StandardModel & SM

Member Function Documentation

◆ AFB()

double LEP2TwoFermions::AFB ( const double  s,
const double  Mw,
const double  GammaZ,
const double  I3f,
const double  Qf,
const double  mf,
const double  mfp,
const bool  bWeak 
) const
private

Definition at line 292 of file LEP2TwoFermions.cpp.

296{
297 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
298 double G1 = SM.getMyTwoFermionsLEP2()->G_1_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
299 double G2 = SM.getMyTwoFermionsLEP2()->G_2_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
300 double G3 = SM.getMyTwoFermionsLEP2()->G_3_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
301
302 return ( 3.0/4.0*betaf*G3/(G1 + 2.0*mf*mf/s*G2) );
303}
double G_2_noBox(const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWeak) const
double G_1_noBox(const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWeak) const
double G_3_noBox(const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWeak) const
An observable class for the total decay width of the boson.
Definition: GammaZ.h:32
An observable class for the -boson mass.
Definition: Mw.h:22
EWSMTwoFermionsLEP2 * getMyTwoFermionsLEP2() const
A get method to retrieve the member pointer of type EWSMTwoFermionsLEP2.
Test Observable.

◆ AFB_l()

double LEP2TwoFermions::AFB_l ( const QCD::lepton  l,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]llepton in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the forward-backward asymmetry for e^+ e^- -> l lbar

Definition at line 109 of file LEP2TwoFermions.cpp.

112{
113 double I3f = SM.getLeptons(l).getIsospin();
114 double Qf = SM.getLeptons(l).getCharge();
115
116 return AFB(s, Mw, GammaZ, I3f, Qf, mf, 0.0, bWeak);
117}
double AFB(const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWeak) const
double getIsospin() const
A get method to access the particle isospin.
Definition: Particle.h:115
double getCharge() const
A get method to access the particle charge.
Definition: Particle.h:97
const Particle & getLeptons(const QCD::lepton p) const
A get method to retrieve the member object of a lepton.

◆ AFB_q()

double LEP2TwoFermions::AFB_q ( const QCD::quark  q,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the forward-backward asymmetry for e^+ e^- -> q qbar

Definition at line 120 of file LEP2TwoFermions.cpp.

123{
124 double I3f = SM.getQuarks(q).getIsospin();
125 double Qf = SM.getQuarks(q).getCharge();
126 double mfp;
127 if (q==SM.TOP)
128 throw std::runtime_error("Error in LEP2TwoFermions::AFB_q()");
129 else if (q==SM.BOTTOM)
130 mfp = SM.getMtpole();
131 else
132 mfp = 0.0;
133
134 return ( AFB(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak) );
135}
@ BOTTOM
Definition: QCD.h:329
@ TOP
Definition: QCD.h:328
const double getMtpole() const
A get method to access the pole mass of the top quark.
Definition: QCD.h:600
const Particle & getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:536

◆ alpha_at_s()

double LEP2TwoFermions::alpha_at_s ( const double  s) const
private
Parameters
sinvariant mass squared
Returns
the electromagnetic coupling at s

Definition at line 227 of file LEP2TwoFermions.cpp.

228{
229 double alpha;
230 //alpha = SM.getAle()/complex(1.0715119759, -0.0186242179, false).real(); // for debug, s=(200GeV)^2
231 alpha = SM.ale_OS(sqrt(s), FULLNLO);
232
233 return alpha;
234}
@ FULLNLO
Definition: OrderScheme.h:38
const double ale_OS(const double mu, orders order=FULLNLO) const
The running electromagnetic coupling in the on-shell scheme.

◆ dsigma()

double LEP2TwoFermions::dsigma ( const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ,
const double  I3f,
const double  Qf,
const double  mf,
const double  mfp,
const double  Ncf,
const bool  bWeak 
) const
private

Definition at line 237 of file LEP2TwoFermions.cpp.

242{
243 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
244 double G1 = SM.getMyTwoFermionsLEP2()->G_1_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
245 double G2 = SM.getMyTwoFermionsLEP2()->G_2_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
246 double G3 = SM.getMyTwoFermionsLEP2()->G_3_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
247
248 return ( M_PI*SM.getAle()*SM.getAle()/s*Ncf*betaf
249 *( G1*(1.0 + cosTheta*cosTheta)/2.0
250 + 2.0*mf*mf/s*G2*(1.0 - cosTheta*cosTheta)
251 + betaf*G3*cosTheta ) );
252}
const double getAle() const
A get method to retrieve the fine-structure constant .

◆ dsigma_box()

double LEP2TwoFermions::dsigma_box ( const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ,
const double  I3f,
const double  Qf,
const double  mf,
const double  mfp,
const double  Ncf 
) const
private

Definition at line 255 of file LEP2TwoFermions.cpp.

260{
261 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
262
263 //double t = mf*mf - s/2.0*(1.0 - betaf*cosTheta);
264 double t = - s/2.0*(1.0 - betaf*cosTheta);
265
266 double G1 = SM.getMyTwoFermionsLEP2()->G_1_box(s, t, Mw, GammaZ, I3f, Qf, mf, mfp);
267 double G2 = SM.getMyTwoFermionsLEP2()->G_2_box(s, t, Mw, GammaZ, I3f, Qf, mf, mfp);
268 double G3 = SM.getMyTwoFermionsLEP2()->G_3_box(s, t, Mw, GammaZ, I3f, Qf, mf, mfp);
269
270 return ( M_PI*SM.getAle()*SM.getAle()/s*Ncf*betaf
271 *( G1*(1.0 + cosTheta*cosTheta)/2.0
272 + 2.0*mf*mf/s*G2*(1.0 - cosTheta*cosTheta)
273 + betaf*G3*cosTheta ) );
274}
double G_3_box(const double s, const double t, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWWbox=true, const bool bZZbox=true) const
double G_1_box(const double s, const double t, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWWbox=true, const bool bZZbox=true) const
double G_2_box(const double s, const double t, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const bool bWWbox=true, const bool bZZbox=true) const
Test Observable.

◆ dsigma_l()

double LEP2TwoFermions::dsigma_l ( const QCD::lepton  l,
const double  mf,
const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]llepton in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]cosThetacosine of the scattering angle theta
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the differential cross section d sigma(e^+ e^- -> l lbar)/d cosTheta in GeV^{-2}

Definition at line 20 of file LEP2TwoFermions.cpp.

24{
25 double I3f = SM.getLeptons(l).getIsospin();
26 double Qf = SM.getLeptons(l).getCharge();
27
28 return ( dsigma(s, cosTheta, Mw, GammaZ, I3f, Qf, mf, 0.0, 1.0, bWeak) );
29}
double dsigma(const double s, const double cosTheta, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf, const bool bWeak) const

◆ dsigma_l_box()

double LEP2TwoFermions::dsigma_l_box ( const QCD::lepton  l,
const double  mf,
const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ 
) const
Parameters
[in]llepton in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]cosThetacosine of the scattering angle theta
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
Returns
the box contribution to the differential cross section d sigma(e^+ e^- -> l lbar)/d cosTheta in GeV^{-2}

Definition at line 51 of file LEP2TwoFermions.cpp.

54{
55 double I3f = SM.getLeptons(l).getIsospin();
56 double Qf = SM.getLeptons(l).getCharge();
57
58 return ( dsigma_box(s, cosTheta, Mw, GammaZ, I3f, Qf, mf, 0.0, 1.0) );
59}
double dsigma_box(const double s, const double cosTheta, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf) const

◆ dsigma_q()

double LEP2TwoFermions::dsigma_q ( const QCD::quark  q,
const double  mf,
const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]cosThetacosine of the scattering angle theta
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the differential cross section d sigma(e^+ e^- -> q qbar)/d cosTheta in GeV^{-2}

Definition at line 32 of file LEP2TwoFermions.cpp.

36{
37 double I3f = SM.getQuarks(q).getIsospin();
38 double Qf = SM.getQuarks(q).getCharge();
39 double mfp;
40 if (q==SM.TOP)
41 throw std::runtime_error("Error in LEP2TwoFermions::sigma_q()");
42 else if (q==SM.BOTTOM)
43 mfp = SM.getMtpole();
44 else
45 mfp = 0.0;
46
47 return ( dsigma(s, cosTheta, Mw, GammaZ, I3f, Qf, mf, mfp, 3.0, bWeak) );
48}

◆ dsigma_q_box()

double LEP2TwoFermions::dsigma_q_box ( const QCD::quark  q,
const double  mf,
const double  s,
const double  cosTheta,
const double  Mw,
const double  GammaZ 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]cosThetacosine of the scattering angle theta
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
Returns
the box contribution to the differential cross section d sigma(e^+ e^- -> q qbar)/d cosTheta in GeV^{-2}

Definition at line 62 of file LEP2TwoFermions.cpp.

65{
66 double I3f = SM.getQuarks(q).getIsospin();
67 double Qf = SM.getQuarks(q).getCharge();
68 double mfp;
69 if (q==SM.TOP)
70 throw std::runtime_error("Error in LEP2TwoFermions::sigma_q()");
71 else if (q==SM.BOTTOM)
72 mfp = SM.getMtpole();
73 else
74 mfp = 0.0;
75
76 return ( dsigma_box(s, cosTheta, Mw, GammaZ, I3f, Qf, mf, mfp, 3.0) );
77}

◆ G_3prime_l()

double LEP2TwoFermions::G_3prime_l ( const QCD::lepton  l,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]llepton in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the form factor beta_f^2*G_3(s) for e^+ e^- -> l lbar

Definition at line 190 of file LEP2TwoFermions.cpp.

194{
195 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
196 double I3f = SM.getLeptons(l).getIsospin();
197 double Qf = SM.getLeptons(l).getCharge();
198 double G3 = SM.getMyTwoFermionsLEP2()->G_3_noBox(s, Mw, GammaZ, I3f, Qf, mf, 0.0, bWeak);
199
200 return ( betaf*betaf*G3 );
201}

◆ G_3prime_q()

double LEP2TwoFermions::G_3prime_q ( const QCD::quark  q,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the form factor beta_f^2*G_3(s) for e^+ e^- -> q qbar

Definition at line 204 of file LEP2TwoFermions.cpp.

208{
209 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
210 double I3f = SM.getQuarks(q).getIsospin();
211 double Qf = SM.getQuarks(q).getCharge();
212 double mfp;
213 if (q==SM.TOP)
214 throw std::runtime_error("Error in LEP2TwoFermions::G_3prime_q()");
215 else if (q==SM.BOTTOM)
216 mfp = SM.getMtpole();
217 else
218 mfp = 0.0;
219 double G3 = SM.getMyTwoFermionsLEP2()->G_3_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
220
221 return ( betaf*betaf*G3 );
222}

◆ H_ISR()

double LEP2TwoFermions::H_ISR ( const double  x,
const double  s 
) const
Parameters
[in]xs'=(1-x)s
[in]sthe invariant mass squared of the initial-state e^+ e^- pair
Returns
the additive radiator of initial-state radiations in cross sections

Definition at line 160 of file LEP2TwoFermions.cpp.

161{
162 double me = SM.getLeptons(SM.ELECTRON).getMass();
163 double alphaOverPi = SM.getAle()/M_PI; // alpha(0)/Pi
164 double L = log(s/(me*me));
165 double beta = 2.0*alphaOverPi*(L - 1.0);
166 double deltaVS_1 = 3.0/2.0*L + M_PI*M_PI/3.0 - 2.0;
167 double deltaH_1 = - (2.0 - x)*(L - 1.0);
168
169 return ( beta*pow(x, beta-1.0)*(1.0 + alphaOverPi*deltaVS_1)
170 + alphaOverPi*deltaH_1 );
171}
const double & getMass() const
A get method to access the particle mass.
Definition: Particle.h:61
@ ELECTRON
Definition: QCD.h:312

◆ H_ISR_FB()

double LEP2TwoFermions::H_ISR_FB ( const double  x,
const double  s 
) const
Parameters
[in]xs'=(1-x)s
[in]sthe invariant mass squared of the initial-state e^+ e^- pair
Returns
the additive radiator of initial-state radiations in forward-backward asysmmetries

Definition at line 174 of file LEP2TwoFermions.cpp.

175{
176 double me = SM.getLeptons(SM.ELECTRON).getMass();
177 double alphaOverPi = SM.getAle()/M_PI; // alpha(0)/Pi
178 double L = log(s/(me*me));
179 double beta = 2.0*alphaOverPi*(L - 1.0);
180 double deltaVS_1 = 3.0/2.0*L + M_PI*M_PI/3.0 - 2.0;
181 double tmp = (1.0-x)/(1.0-x/2.0)/(1.0-x/2.0);
182 double deltaH_FB_1 = (1.0+(1.0-x)*(1.0-x))/x*tmp*(L - 1.0 - log(tmp))
183 - 2.0/x*(L - 1.0);
184
185 return ( beta*pow(x, beta-1.0)*(1.0 + alphaOverPi*deltaVS_1)
186 + alphaOverPi*deltaH_FB_1 );
187}

◆ QCD_FSR_forAFB()

double LEP2TwoFermions::QCD_FSR_forAFB ( const QCD::quark  q,
const double  mf,
const double  s 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sthe invariant mass squared of the initial-state e^+ e^- pair
Returns
the final-state QCD corrections to forward-backward asymmetries

Definition at line 144 of file LEP2TwoFermions.cpp.

146{
147 return ( 1.0 - SM.Als(sqrt(s), FULLNLO)/M_PI*(1.0 - 16.0/3.0*mf/sqrt(s)) );
148}
const double Als(const double mu, const orders order, const bool Nf_thr, const bool qed_flag) const
The running QCD coupling in the scheme including QED corrections.

◆ QCD_FSR_forSigma()

double LEP2TwoFermions::QCD_FSR_forSigma ( const double  s) const
Parameters
[in]sthe invariant mass squared of the initial-state e^+ e^- pair
Returns
the final-state QCD corrections to cross sections

Definition at line 138 of file LEP2TwoFermions.cpp.

139{
140 return ( 1.0 + SM.Als(sqrt(s), FULLNLO)/M_PI );
141}

◆ QED_FSR_forSigma()

double LEP2TwoFermions::QED_FSR_forSigma ( const double  s,
const double  Qf 
) const
Parameters
[in]sthe invariant mass squared of the initial-state e^+ e^- pair
[in]Qfthe electromagnetic charge of the final-state fermion
Returns
the final-state QED corrections to cross sections

Definition at line 151 of file LEP2TwoFermions.cpp.

152{
153 //double alpha = SM.getAle();
154 double alpha = alpha_at_s(s);
155
156 return ( 1.0 + 3.0*alpha/(4.0*M_PI)*Qf*Qf );
157}
double alpha_at_s(const double s) const

◆ sigma()

double LEP2TwoFermions::sigma ( const double  s,
const double  Mw,
const double  GammaZ,
const double  I3f,
const double  Qf,
const double  mf,
const double  mfp,
const double  Ncf,
const bool  bWeak 
) const
private

Definition at line 277 of file LEP2TwoFermions.cpp.

282{
283 double betaf = sqrt(1.0 - 4.0*mf*mf/s);
284 double G1 = SM.getMyTwoFermionsLEP2()->G_1_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
285 double G2 = SM.getMyTwoFermionsLEP2()->G_2_noBox(s, Mw, GammaZ, I3f, Qf, mf, mfp, bWeak);
286
287 return ( 4.0*M_PI*SM.getAle()*SM.getAle()/(3.0*s)*Ncf*betaf
288 *(G1 + 2.0*mf*mf/s*G2) );
289}

◆ sigma_l()

double LEP2TwoFermions::sigma_l ( const QCD::lepton  l,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]llepton in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the total cross section for e^+ e^- -> l lbar in GeV^{-2}

Definition at line 80 of file LEP2TwoFermions.cpp.

83{
84 double I3f = SM.getLeptons(l).getIsospin();
85 double Qf = SM.getLeptons(l).getCharge();
86
87 return ( sigma(s, Mw, GammaZ, I3f, Qf, mf, 0.0, 1.0, bWeak) );
88}
double sigma(const double s, const double Mw, const double GammaZ, const double I3f, const double Qf, const double mf, const double mfp, const double Ncf, const bool bWeak) const

◆ sigma_q()

double LEP2TwoFermions::sigma_q ( const QCD::quark  q,
const double  mf,
const double  s,
const double  Mw,
const double  GammaZ,
const bool  bWeak 
) const
Parameters
[in]qquark in the final state
[in]mfthe mass of the final-state fermion
[in]sinvariant mass squared of the initial-state e^+ e^- pair
[in]Mwthe W-boson mass
[in]GammaZthe Z-boson decay width
[in]bWeakflag to control weak corrections (not including box diagrams)
Returns
the total cross section for e^+ e^- -> q qbar in GeV^{-2}

Definition at line 91 of file LEP2TwoFermions.cpp.

94{
95 double I3f = SM.getQuarks(q).getIsospin();
96 double Qf = SM.getQuarks(q).getCharge();
97 double mfp;
98 if (q==SM.TOP)
99 throw std::runtime_error("Error in LEP2TwoFermions::sigma_q()");
100 else if (q==SM.BOTTOM)
101 mfp = SM.getMtpole();
102 else
103 mfp = 0.0;
104
105 return ( sigma(s, Mw, GammaZ, I3f, Qf, mf, mfp, 3.0, bWeak) );
106}

Member Data Documentation

◆ SM

const StandardModel& LEP2TwoFermions::SM
private

Definition at line 225 of file LEP2TwoFermions.h.


The documentation for this class was generated from the following files: