Geant4-11
Public Member Functions | Private Attributes
G4AngularDistribution Class Reference

#include <G4AngularDistribution.hh>

Inheritance diagram for G4AngularDistribution:
G4VAngularDistribution

Public Member Functions

virtual G4double CosTheta (G4double s, G4double m1, G4double m2) const
 
G4double Cross (G4double tpPion, G4double tpSigma, G4double tpOmega, G4double tmPion, G4double tmSigma, G4double tmOmega, G4double bMix_o1, G4double bMix_s1, G4double bMix_Omega, G4double bMix_sm, G4double bMix_oL, G4double bMix_sL, G4double bOmega_0, G4double bOmega_1, G4double bOmega_2, G4double bOmega_3, G4double bOmega_m, G4double bOmega_L) const
 
G4double DifferentialCrossSection (G4double sIn, G4double m1, G4double m2, G4double cosTheta) const
 
 G4AngularDistribution (G4bool symmetrize)
 
virtual G4double Phi () const
 
virtual ~G4AngularDistribution ()
 

Private Attributes

G4double cm2go
 
G4double cm2gs
 
G4double cm6gp
 
G4double cMix_o1
 
G4double cMix_oLc
 
G4double cMix_oLs
 
G4double cMix_Omega
 
G4double cMix_s1
 
G4double cMix_sLc
 
G4double cMix_sLs
 
G4double cMix_sm
 
G4double cmOmega
 
G4double cmOmega2
 
G4double cmOmega4
 
G4double cmOmega6
 
G4double cmPion
 
G4double cmPion2
 
G4double cmSigma
 
G4double cmSigma2
 
G4double cmSigma4
 
G4double cmSigma6
 
G4double cOmega_1
 
G4double cOmega_2
 
G4double cOmega_3
 
G4double cOmega_L
 
G4double cOmega_m
 
G4double cPion_0
 
G4double cPion_1
 
G4double cPion_2
 
G4double cPion_3
 
G4double cPion_L
 
G4double cPion_m
 
G4double cSigma_0
 
G4double cSigma_1
 
G4double cSigma_2
 
G4double cSigma_3
 
G4double cSigma_L
 
G4double cSigma_m
 
G4double dMix1
 
G4double dMix2
 
G4double dMix3
 
G4double dOmega1
 
G4double dOmega2
 
G4double dOmega3
 
G4double dPion1
 
G4double dPion2
 
G4double dSigma1
 
G4double dSigma2
 
G4double dSigma3
 
G4double fac1
 
G4double fac2
 
G4double fac3
 
G4double gOmega
 
G4double gPion
 
G4double gSigma
 
G4double m42
 
G4double mNucleon
 
G4double mOmega
 
G4double mOmega2
 
G4double mPion
 
G4double mPion2
 
G4double mSigma
 
G4double mSigma2
 
G4double sOmega1
 
G4bool sym
 

Detailed Description

Definition at line 42 of file G4AngularDistribution.hh.

Constructor & Destructor Documentation

◆ G4AngularDistribution()

G4AngularDistribution::G4AngularDistribution ( G4bool  symmetrize)

Definition at line 34 of file G4AngularDistribution.cc.

35 : sym(symmetrize)
36{
37 // The following are parameters of the model - not to be confused with the PDG values!
38
39 mSigma = 0.55;
40 cmSigma = 1.20;
41 gSigma = 9.4;
42
43 mOmega = 0.783;
44 cmOmega = 0.808;
45 gOmega = 10.95;
46
47 mPion = 0.138;
48 cmPion = 0.51;
49 gPion = 7.27;
50
51 mNucleon = 0.938;
52
53 // Definition of constants for pion-Term (no s-dependence)
54
55 m42 = 4. * mNucleon * mNucleon;
56 mPion2 = mPion * mPion;
61
62 cPion_3 = -(cm6gp/3.);
64 cPion_1 = -(cm6gp * mPion2 * (2. * cmPion2 + mPion2) / dPion2);
66 cPion_L = -(cm6gp * 2. * cmPion2 * mPion2 * (cmPion2 + mPion2) / dPion2 / dPion1);
68
69 // Definition of constants for sigma-Term (no s-dependence)
70
71 G4double gSigmaSq = gSigma * gSigma;
72
80
81 G4double dSigma1Sq = dSigma1 * dSigma1;
82 G4double dSigma2Sq = dSigma2 * dSigma2;
83 G4double dSigma3Sq = dSigma3 * dSigma3;
84
85 cm2gs = 0.5 * cmSigma2 * gSigmaSq*gSigmaSq / dSigma3Sq;
86
87
88 cSigma_3 = -(cm2gs * dSigma1Sq / 3.);
90 cSigma_1 = -(cm2gs * cmSigma4 * (2. * dSigma1 + dSigma2) * dSigma2 / dSigma3Sq);
91 cSigma_m = -(cm2gs * cmSigma6 * dSigma2Sq / mSigma2 / dSigma3Sq);
92 cSigma_L = -(cm2gs * cmSigma6 * dSigma2 * (dSigma1 + dSigma2) * 2. / (dSigma3 * dSigma3Sq));
94
95 // Definition of constants for omega-Term
96
97 G4double gOmegaSq = gOmega * gOmega;
98
103 dOmega1 = m42 - cmOmega2;
104 dOmega2 = m42 - mOmega2;
107
108 G4double dOmega3Sq = dOmega3 * dOmega3;
109
110 cm2go = 0.5 * cmOmega2 * gOmegaSq * gOmegaSq / dOmega3Sq;
111
112 cOmega_3 = cm2go / 3.;
114 cOmega_1 = cm2go * cmOmega4 / dOmega3Sq;
115 cOmega_m = cm2go * cmOmega6 / (dOmega3Sq * mOmega2);
116 cOmega_L = -(cm2go * cmOmega6 * 4. / (dOmega3 * dOmega3Sq));
117
118 // Definition of constants for mix-Term
119
120 G4double fac1Tmp = (gSigma * gOmega * cmSigma2 * cmOmega2);
121 fac1 = -(fac1Tmp * fac1Tmp * m42);
125
126 G4double dMix1Sq = dMix1 * dMix1;
127 G4double dMix2Sq = dMix2 * dMix2;
128 G4double dMix3Sq = dMix3 * dMix3;
129
130 cMix_o1 = fac1 / (cmOmega2 * dMix1Sq * dMix2 * dOmega3);
131 cMix_s1 = fac1 / (cmSigma2 * dMix1Sq * dMix3 * dSigma3);
132 cMix_Omega = fac1 / (dOmega3Sq * dMix3Sq * (mOmega2 - mSigma2));
133 cMix_sm = fac1 / (dSigma3Sq * dMix2Sq * (mSigma2 - mOmega2));
134 fac2 = (-fac1) / (dMix1*dMix1Sq * dOmega3Sq * dMix2Sq);
135 fac3 = (-fac1) / (dMix1*dMix1Sq * dSigma3Sq * dMix3Sq);
136
138 - 2. * cmOmega4 * mOmega2 - 2. * cmOmega4 * mSigma2
140 - 4. * cmOmega4 * m42 + 2. * cmOmega2 * cmSigma2 * m42
141 + 3. * cmOmega2 * mOmega2 * m42 - cmSigma2 * mOmega2 * m42
142 + 3. * cmOmega2 * mSigma2 * m42 - cmSigma2 * mSigma2 * m42
143 - 2. * mOmega2 * mSigma2 * m42);
144
145 cMix_oLs = fac2 * (8. * cmOmega4 - 4. * cmOmega2 * cmSigma2
146 - 6. * cmOmega2 * mOmega2 + 2. * cmSigma2 * mOmega2
147 - 6. * cmOmega2 * mSigma2 + 2. * cmSigma2 * mSigma2
148 + 4. * mOmega2 * mSigma2);
149
150 cMix_sLc = fac3 * (cmOmega2 * cmSigma4 - 3. * cmSigma6
151 + 2. * cmSigma4 * mOmega2 + 2. * cmSigma4 * mSigma2
153 - 2. * cmOmega2 * cmSigma2 * m42 + 4. * cmSigma4 * m42
154 + cmOmega2 * mOmega2 * m42 - 3. * cmSigma2 * mOmega2 * m42
155 + cmOmega2 * mSigma2 * m42 - 3. * cmSigma2 * mSigma2 * m42
156 + 2. * mOmega2 * mSigma2 * m42);
157
158 cMix_sLs = fac3 * (4. * cmOmega2 * cmSigma2 - 8. * cmSigma4
159 - 2. * cmOmega2 * mOmega2 + 6. * cmSigma2 * mOmega2
160 - 2. * cmOmega2 * mSigma2 + 6. * cmSigma2 * mSigma2
161 - 4. * mOmega2 * mSigma2);
162}
double G4double
Definition: G4Types.hh:83

References cm2go, cm2gs, cm6gp, cMix_o1, cMix_oLc, cMix_oLs, cMix_Omega, cMix_s1, cMix_sLc, cMix_sLs, cMix_sm, cmOmega, cmOmega2, cmOmega4, cmOmega6, cmPion, cmPion2, cmSigma, cmSigma2, cmSigma4, cmSigma6, cOmega_1, cOmega_2, cOmega_3, cOmega_L, cOmega_m, cPion_0, cPion_1, cPion_2, cPion_3, cPion_L, cPion_m, cSigma_0, cSigma_1, cSigma_2, cSigma_3, cSigma_L, cSigma_m, dMix1, dMix2, dMix3, dOmega1, dOmega2, dOmega3, dPion1, dPion2, dSigma1, dSigma2, dSigma3, fac1, fac2, fac3, gOmega, gPion, gSigma, m42, mNucleon, mOmega, mOmega2, mPion, mPion2, mSigma, mSigma2, and sOmega1.

◆ ~G4AngularDistribution()

G4AngularDistribution::~G4AngularDistribution ( )
virtual

Definition at line 165 of file G4AngularDistribution.cc.

167{ }

Member Function Documentation

◆ CosTheta()

G4double G4AngularDistribution::CosTheta ( G4double  s,
G4double  m1,
G4double  m2 
) const
virtual

Implements G4VAngularDistribution.

Definition at line 170 of file G4AngularDistribution.cc.

171{
172 G4double random = G4UniformRand();
173 G4double dCosTheta = 2.;
174 G4double cosTheta = -1.;
175
176 // For jmax=12 the accuracy is better than 0.1 degree
177 G4int jMax = 12;
178
179 for (G4int j = 1; j <= jMax; ++j)
180 {
181 // Accuracy is 2^-jmax
182 dCosTheta *= 0.5;
183 G4double cosTh = cosTheta + dCosTheta;
184 if(DifferentialCrossSection(S, m_1, m_2, cosTh) <= random) cosTheta = cosTh;
185 }
186
187 // Randomize in final interval in order to avoid discrete angles
188 cosTheta += G4UniformRand() * dCosTheta;
189
190
191 if (cosTheta > 1. || cosTheta < -1.)
192 throw G4HadronicException(__FILE__, __LINE__, "G4AngularDistribution::CosTheta - std::cos(theta) outside allowed range");
193
194 return cosTheta;
195}
G4double S(G4double temp)
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
G4double DifferentialCrossSection(G4double sIn, G4double m1, G4double m2, G4double cosTheta) const

References DifferentialCrossSection(), G4UniformRand, and S().

◆ Cross()

G4double G4AngularDistribution::Cross ( G4double  tpPion,
G4double  tpSigma,
G4double  tpOmega,
G4double  tmPion,
G4double  tmSigma,
G4double  tmOmega,
G4double  bMix_o1,
G4double  bMix_s1,
G4double  bMix_Omega,
G4double  bMix_sm,
G4double  bMix_oL,
G4double  bMix_sL,
G4double  bOmega_0,
G4double  bOmega_1,
G4double  bOmega_2,
G4double  bOmega_3,
G4double  bOmega_m,
G4double  bOmega_L 
) const

Definition at line 295 of file G4AngularDistribution.cc.

313{
314 G4double cross = 0;
315 // Pion
316 cross += ((cPion_3 * tpPion + cPion_2) * tpPion + cPion_1) * tpPion + cPion_m/tmPion + cPion_0 + cPion_L * G4Log(tpPion*tmPion);
317// G4cout << "cross1 "<< cross<<G4endl;
318 // Sigma
319 cross += ((cSigma_3 * tpSigma + cSigma_2) * tpSigma + cSigma_1) * tpSigma + cSigma_m/tmSigma + cSigma_0 + cSigma_L * G4Log(tpSigma*tmSigma);
320// G4cout << "cross2 "<< cross<<G4endl;
321 // Omega
322 cross += ((bOmega_3 * tpOmega + bOmega_2) * tpOmega + bOmega_1) * tpOmega + bOmega_m/tmOmega + bOmega_0 + bOmega_L * G4Log(tpOmega*tmOmega)
323 // Mix
324 + bMix_o1 * (tpOmega - 1.)
325 + bMix_s1 * (tpSigma - 1.)
326 + bMix_Omega * G4Log(tmOmega)
327 + bMix_sm * G4Log(tmSigma)
328 + bMix_oL * G4Log(tpOmega)
329 + bMix_sL * G4Log(tpSigma);
330/* G4cout << "cross3 "<< cross<<" "
331 <<bMix_o1<<" "
332 <<bMix_s1<<" "
333 <<bMix_Omega<<" "
334 <<bMix_sm<<" "
335 <<bMix_oL<<" "
336 <<bMix_sL<<" "
337 <<tpOmega<<" "
338 <<tpSigma<<" "
339 <<tmOmega<<" "
340 <<tmSigma<<" "
341 <<tpOmega<<" "
342 <<tpSigma
343 <<G4endl;
344*/
345 return cross;
346
347}
G4double G4Log(G4double x)
Definition: G4Log.hh:226

References cPion_0, cPion_1, cPion_2, cPion_3, cPion_L, cPion_m, cSigma_0, cSigma_1, cSigma_2, cSigma_3, cSigma_L, cSigma_m, and G4Log().

Referenced by DifferentialCrossSection().

◆ DifferentialCrossSection()

G4double G4AngularDistribution::DifferentialCrossSection ( G4double  sIn,
G4double  m1,
G4double  m2,
G4double  cosTheta 
) const

Definition at line 198 of file G4AngularDistribution.cc.

200{
201// local calculus is in GeV, ie. normalize input
202 sIn = sIn/sqr(GeV)+m42/2.;
203 m_1 = m_1/GeV;
204 m_2 = m_2/GeV;
205// G4cout << "Here we go"<<sIn << " "<<m1 << " " << m2 <<" " m42<< G4endl;
206// scaling from masses other than p,p.
207 G4double S = sIn - (m_1+m_2) * (m_1+m_2) + m42;
208 G4double tMax = S - m42;
209 G4double tp = 0.5 * (cosTheta + 1.) * tMax;
210 G4double twoS = 2. * S;
211
212 // Define s-dependent stuff for omega-Term
213 G4double brak1 = (twoS-m42) * (twoS-m42);
214 G4double bOmega_3 = cOmega_3 * (-2. * cmOmega4 - 2. * cmOmega2 * twoS - brak1);
215 G4double bOmega_2 = cOmega_2 * ( 2. * cmOmega2 * mOmega2 + sOmega1 * twoS + brak1);
216 G4double bOmega_1 = cOmega_1 * (-4. * cmOmega2 * mOmega2
217 - 2. * mOmega2*mOmega2
218 - 2. * (cmOmega2 + 2 * mOmega2) * twoS
219 - 3. * brak1);
220 G4double bOmega_m = cOmega_m * (-2. * mOmega2*mOmega2 - 2. * mOmega2 * twoS - brak1);
221 G4double bOmega_L = cOmega_L * (sOmega1 * mOmega2 + (cmOmega2 + 3. * mOmega2) * S + brak1);
222 G4double bOmega_0 = -(bOmega_3 + bOmega_2 + bOmega_1 + bOmega_m);
223
224 // Define s-dependent stuff for mix-Term
225 G4double bMix_o1 = cMix_o1 * (dOmega1 - twoS);
226 G4double bMix_s1 = cMix_s1 * (dSigma1 - twoS);
227 G4double bMix_Omega = cMix_Omega * (dOmega2 - twoS);
228 G4double bMix_sm = cMix_sm * (dSigma2 - twoS);
229 G4double bMix_oL = cMix_oLc + cMix_oLs * S;
230 G4double bMix_sL = cMix_sLc + cMix_sLs * S;
231
232 G4double t1_Pion = 1. / (1. + tMax / cmPion2);
233 G4double t2_Pion = 1. + tMax / mPion2;
234 G4double t1_Sigma = 1. / (1. + tMax / cmSigma2);
235 G4double t2_Sigma = 1. + tMax / mSigma2;
236 G4double t1_Omega = 1. / (1. + tMax / cmOmega2);
237 G4double t2_Omega = 1. + tMax / mOmega2;
238
239 G4double norm = Cross(t1_Pion, t1_Sigma, t1_Omega,
240 t2_Pion, t2_Sigma, t2_Omega,
241 bMix_o1, bMix_s1, bMix_Omega,
242 bMix_sm, bMix_oL, bMix_sL,
243 bOmega_0, bOmega_1, bOmega_2,
244 bOmega_3, bOmega_m, bOmega_L);
245
246 t1_Pion = 1. / (1. + tp / cmPion2);
247 t2_Pion = 1. + tp / mPion2;
248 t1_Sigma = 1. / (1. + tp / cmSigma2);
249 t2_Sigma = 1. + tp / mSigma2;
250 t1_Omega = 1. / (1. + tp / cmOmega2);
251 t2_Omega = 1. + tp / mOmega2;
252
253 G4double dSigma;
254 if (sym)
255 {
256 G4double to;
257 norm = 2. * norm;
258 to = tMax - tp;
259 G4double t3_Pion = 1. / (1. + to / cmPion2);
260 G4double t4_Pion = 1. + to / mPion2;
261 G4double t3_Sigma = 1. / (1. + to / cmSigma2);
262 G4double t4_Sigma = 1. + to / mSigma2;
263 G4double t3_Omega = 1. / (1. + to / cmOmega2);
264 G4double t4_Omega = 1. + to / mOmega2;
265
266 dSigma = ( Cross(t1_Pion, t1_Sigma, t1_Omega,
267 t2_Pion,t2_Sigma, t2_Omega,
268 bMix_o1, bMix_s1, bMix_Omega,
269 bMix_sm, bMix_oL, bMix_sL,
270 bOmega_0, bOmega_1, bOmega_2,
271 bOmega_3, bOmega_m, bOmega_L) -
272 Cross(t3_Pion,t3_Sigma, t3_Omega,
273 t4_Pion, t4_Sigma, t4_Omega,
274 bMix_o1, bMix_s1, bMix_Omega,
275 bMix_sm, bMix_oL, bMix_sL,
276 bOmega_0, bOmega_1, bOmega_2,
277 bOmega_3, bOmega_m, bOmega_L) )
278 / norm + 0.5;
279 }
280 else
281 {
282 dSigma = Cross(t1_Pion, t1_Sigma, t1_Omega,
283 t2_Pion, t2_Sigma, t2_Omega,
284 bMix_o1, bMix_s1, bMix_Omega,
285 bMix_sm, bMix_oL, bMix_sL,
286 bOmega_0, bOmega_1, bOmega_2,
287 bOmega_3, bOmega_m, bOmega_L)
288 / norm;
289 }
290
291 return dSigma;
292}
static constexpr double GeV
Definition: G4SIunits.hh:203
G4double Cross(G4double tpPion, G4double tpSigma, G4double tpOmega, G4double tmPion, G4double tmSigma, G4double tmOmega, G4double bMix_o1, G4double bMix_s1, G4double bMix_Omega, G4double bMix_sm, G4double bMix_oL, G4double bMix_sL, G4double bOmega_0, G4double bOmega_1, G4double bOmega_2, G4double bOmega_3, G4double bOmega_m, G4double bOmega_L) const
T sqr(const T &x)
Definition: templates.hh:128

References cMix_o1, cMix_oLc, cMix_oLs, cMix_Omega, cMix_s1, cMix_sLc, cMix_sLs, cMix_sm, cmOmega2, cmOmega4, cmPion2, cmSigma2, cOmega_1, cOmega_2, cOmega_3, cOmega_L, cOmega_m, Cross(), dOmega1, dOmega2, dSigma1, dSigma2, GeV, m42, mOmega2, mPion2, mSigma2, S(), sOmega1, sqr(), sym, and G4InuclParticleNames::tp.

Referenced by CosTheta().

◆ Phi()

virtual G4double G4VAngularDistribution::Phi ( ) const
inlinevirtualinherited

Reimplemented in G4AngularDistributionNP, and G4AngularDistributionPP.

Definition at line 66 of file G4VAngularDistribution.hh.

66{ return 2.*CLHEP::pi*G4UniformRand(); }
static constexpr double pi
Definition: SystemOfUnits.h:55

References G4UniformRand, and CLHEP::pi.

Referenced by G4VElasticCollision::FinalState(), and G4VScatteringCollision::FinalState().

Field Documentation

◆ cm2go

G4double G4AngularDistribution::cm2go
private

Definition at line 131 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cm2gs

G4double G4AngularDistribution::cm2gs
private

Definition at line 111 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cm6gp

G4double G4AngularDistribution::cm6gp
private

Definition at line 93 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cMix_o1

G4double G4AngularDistribution::cMix_o1
private

Definition at line 145 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_oLc

G4double G4AngularDistribution::cMix_oLc
private

Definition at line 152 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_oLs

G4double G4AngularDistribution::cMix_oLs
private

Definition at line 153 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_Omega

G4double G4AngularDistribution::cMix_Omega
private

Definition at line 147 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_s1

G4double G4AngularDistribution::cMix_s1
private

Definition at line 146 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_sLc

G4double G4AngularDistribution::cMix_sLc
private

Definition at line 154 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_sLs

G4double G4AngularDistribution::cMix_sLs
private

Definition at line 155 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cMix_sm

G4double G4AngularDistribution::cMix_sm
private

Definition at line 148 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cmOmega

G4double G4AngularDistribution::cmOmega
private

Definition at line 78 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cmOmega2

G4double G4AngularDistribution::cmOmega2
private

Definition at line 123 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cmOmega4

G4double G4AngularDistribution::cmOmega4
private

Definition at line 124 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cmOmega6

G4double G4AngularDistribution::cmOmega6
private

Definition at line 125 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cmPion

G4double G4AngularDistribution::cmPion
private

Definition at line 76 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cmPion2

G4double G4AngularDistribution::cmPion2
private

Definition at line 90 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cmSigma

G4double G4AngularDistribution::cmSigma
private

Definition at line 77 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cmSigma2

G4double G4AngularDistribution::cmSigma2
private

Definition at line 105 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cmSigma4

G4double G4AngularDistribution::cmSigma4
private

Definition at line 106 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cmSigma6

G4double G4AngularDistribution::cmSigma6
private

Definition at line 107 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ cOmega_1

G4double G4AngularDistribution::cOmega_1
private

Definition at line 135 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cOmega_2

G4double G4AngularDistribution::cOmega_2
private

Definition at line 134 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cOmega_3

G4double G4AngularDistribution::cOmega_3
private

Definition at line 133 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cOmega_L

G4double G4AngularDistribution::cOmega_L
private

Definition at line 137 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cOmega_m

G4double G4AngularDistribution::cOmega_m
private

Definition at line 136 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ cPion_0

G4double G4AngularDistribution::cPion_0
private

Definition at line 100 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cPion_1

G4double G4AngularDistribution::cPion_1
private

Definition at line 97 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cPion_2

G4double G4AngularDistribution::cPion_2
private

Definition at line 96 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cPion_3

G4double G4AngularDistribution::cPion_3
private

Definition at line 95 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cPion_L

G4double G4AngularDistribution::cPion_L
private

Definition at line 99 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cPion_m

G4double G4AngularDistribution::cPion_m
private

Definition at line 98 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_0

G4double G4AngularDistribution::cSigma_0
private

Definition at line 118 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_1

G4double G4AngularDistribution::cSigma_1
private

Definition at line 115 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_2

G4double G4AngularDistribution::cSigma_2
private

Definition at line 114 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_3

G4double G4AngularDistribution::cSigma_3
private

Definition at line 113 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_L

G4double G4AngularDistribution::cSigma_L
private

Definition at line 117 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ cSigma_m

G4double G4AngularDistribution::cSigma_m
private

Definition at line 116 of file G4AngularDistribution.hh.

Referenced by Cross(), and G4AngularDistribution().

◆ dMix1

G4double G4AngularDistribution::dMix1
private

Definition at line 142 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dMix2

G4double G4AngularDistribution::dMix2
private

Definition at line 143 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dMix3

G4double G4AngularDistribution::dMix3
private

Definition at line 144 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dOmega1

G4double G4AngularDistribution::dOmega1
private

Definition at line 126 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ dOmega2

G4double G4AngularDistribution::dOmega2
private

Definition at line 127 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ dOmega3

G4double G4AngularDistribution::dOmega3
private

Definition at line 128 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dPion1

G4double G4AngularDistribution::dPion1
private

Definition at line 91 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dPion2

G4double G4AngularDistribution::dPion2
private

Definition at line 92 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ dSigma1

G4double G4AngularDistribution::dSigma1
private

Definition at line 108 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ dSigma2

G4double G4AngularDistribution::dSigma2
private

Definition at line 109 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ dSigma3

G4double G4AngularDistribution::dSigma3
private

Definition at line 110 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ fac1

G4double G4AngularDistribution::fac1
private

Definition at line 141 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ fac2

G4double G4AngularDistribution::fac2
private

Definition at line 149 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ fac3

G4double G4AngularDistribution::fac3
private

Definition at line 150 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ gOmega

G4double G4AngularDistribution::gOmega
private

Definition at line 82 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ gPion

G4double G4AngularDistribution::gPion
private

Definition at line 80 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ gSigma

G4double G4AngularDistribution::gSigma
private

Definition at line 81 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ m42

G4double G4AngularDistribution::m42
private

Definition at line 88 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ mNucleon

G4double G4AngularDistribution::mNucleon
private

Definition at line 84 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ mOmega

G4double G4AngularDistribution::mOmega
private

Definition at line 74 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ mOmega2

G4double G4AngularDistribution::mOmega2
private

Definition at line 122 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ mPion

G4double G4AngularDistribution::mPion
private

Definition at line 72 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ mPion2

G4double G4AngularDistribution::mPion2
private

Definition at line 89 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ mSigma

G4double G4AngularDistribution::mSigma
private

Definition at line 73 of file G4AngularDistribution.hh.

Referenced by G4AngularDistribution().

◆ mSigma2

G4double G4AngularDistribution::mSigma2
private

Definition at line 104 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ sOmega1

G4double G4AngularDistribution::sOmega1
private

Definition at line 129 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection(), and G4AngularDistribution().

◆ sym

G4bool G4AngularDistribution::sym
private

Definition at line 68 of file G4AngularDistribution.hh.

Referenced by DifferentialCrossSection().


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