Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4WentzelVIRelXSection Class Reference

#include <G4WentzelVIRelXSection.hh>

Public Member Functions

 G4WentzelVIRelXSection ()
 
virtual ~G4WentzelVIRelXSection ()
 
void Initialise (const G4ParticleDefinition *, G4double CosThetaLim)
 
void SetupParticle (const G4ParticleDefinition *)
 
G4double SetupTarget (G4int Z, G4double cut=DBL_MAX)
 
G4double ComputeTransportCrossSectionPerAtom (G4double CosThetaMax)
 
G4ThreeVector SampleSingleScattering (G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
 
G4double ComputeNuclearCrossSection (G4double CosThetaMin, G4double CosThetaMax)
 
G4double ComputeElectronCrossSection (G4double CosThetaMin, G4double CosThetaMax)
 
G4double SetupKinematic (G4double kinEnergy, const G4Material *mat)
 
void SetTargetMass (G4double value)
 
void SetRelativisticMass (G4double value)
 
G4double GetMomentumSquare () const
 
G4double GetCosThetaNuc () const
 
G4double GetCosThetaElec () const
 

Detailed Description

Definition at line 71 of file G4WentzelVIRelXSection.hh.

Constructor & Destructor Documentation

◆ G4WentzelVIRelXSection()

G4WentzelVIRelXSection::G4WentzelVIRelXSection ( )

Definition at line 65 of file G4WentzelVIRelXSection.cc.

65 :
66 numlimit(0.1),
67 nwarnings(0),
68 nwarnlimit(50),
69 alpha2(fine_structure_const*fine_structure_const)
70{
71 fNistManager = G4NistManager::Instance();
72 fG4pow = G4Pow::GetInstance();
73 theElectron = G4Electron::Electron();
74 thePositron = G4Positron::Positron();
75 theProton = G4Proton::Proton();
76 lowEnergyLimit = 1.0*eV;
77 G4double p0 = electron_mass_c2*classic_electr_radius;
78 coeff = twopi*p0*p0;
79 particle = 0;
80
81 // Thomas-Fermi screening radii
82 // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282
83
84 if(0.0 == ScreenRSquare[0]) {
85 G4double a0 = electron_mass_c2/0.88534;
86 G4double constn = 6.937e-6/(MeV*MeV);
87
88 ScreenRSquare[0] = alpha2*a0*a0;
89 for(G4int j=1; j<100; ++j) {
90 G4double x = a0*fG4pow->Z13(j);
91 //ScreenRSquare[j] = 0.5*(1 + exp(-j*j*0.001))*alpha2*x*x;
92 ScreenRSquare[j] = 0.5*alpha2*x*x;
93 x = fNistManager->GetA27(j);
94 FormFactor[j] = constn*x*x;
95 }
96 }
97 currentMaterial = 0;
98 elecXSRatio = factB = factD = formfactA = screenZ = 0.0;
99 cosTetMaxElec = cosTetMaxNuc = invbeta2 = kinFactor = gam0pcmp = pcmp2 = 1.0;
100
101 factB1= 0.5*CLHEP::pi*fine_structure_const;
102
103 Initialise(theElectron, 1.0);
104 targetMass = proton_mass_c2;
105}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetA27(G4int Z)
static G4NistManager * Instance()
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
static G4Proton * Proton()
Definition: G4Proton.cc:93
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)

◆ ~G4WentzelVIRelXSection()

G4WentzelVIRelXSection::~G4WentzelVIRelXSection ( )
virtual

Definition at line 109 of file G4WentzelVIRelXSection.cc.

110{}

Member Function Documentation

◆ ComputeElectronCrossSection()

G4double G4WentzelVIRelXSection::ComputeElectronCrossSection ( G4double  CosThetaMin,
G4double  CosThetaMax 
)
inline

Definition at line 246 of file G4WentzelVIRelXSection.hh.

248{
249 G4double xsec = 0.0;
250 G4double cost1 = std::max(cosTMin,cosTetMaxElec);
251 G4double cost2 = std::max(cosTMax,cosTetMaxElec);
252 if(cost1 > cost2) {
253 xsec = kinFactor*(cost1 - cost2)/((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
254 }
255 return xsec;
256}

Referenced by G4hCoulombScatteringModel::ComputeCrossSectionPerAtom().

◆ ComputeNuclearCrossSection()

G4double G4WentzelVIRelXSection::ComputeNuclearCrossSection ( G4double  CosThetaMin,
G4double  CosThetaMax 
)
inline

Definition at line 232 of file G4WentzelVIRelXSection.hh.

234{
235 G4double xsec = 0.0;
236 if(cosTMax < cosTMin) {
237 xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
238 ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
239 }
240 return xsec;
241}

Referenced by G4hCoulombScatteringModel::ComputeCrossSectionPerAtom().

◆ ComputeTransportCrossSectionPerAtom()

G4double G4WentzelVIRelXSection::ComputeTransportCrossSectionPerAtom ( G4double  CosThetaMax)

Definition at line 192 of file G4WentzelVIRelXSection.cc.

193{
194 G4double xsec = 0.0;
195 if(cosTMax >= 1.0) { return xsec; }
196
197 G4double xSection = 0.0;
198 G4double x = 0;
199 G4double y = 0;
200 G4double x1= 0;
201 G4double x2= 0;
202 G4double xlog = 0.0;
203
204 G4double costm = std::max(cosTMax,cosTetMaxElec);
205 G4double fb = screenZ*factB;
206
207 // scattering off electrons
208 if(costm < 1.0) {
209 x = (1.0 - costm)/screenZ;
210 if(x < numlimit) {
211 x2 = 0.5*x*x;
212 y = x2*(1.0 - 1.3333333*x + 3*x2);
213 if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
214 } else {
215 x1= x/(1 + x);
216 xlog = log(1.0 + x);
217 y = xlog - x1;
218 if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
219 }
220
221 if(y < 0.0) {
222 ++nwarnings;
223 if(nwarnings < nwarnlimit) {
224 G4cout << "G4WentzelVIRelXSection::ComputeTransportCrossSectionPerAtom scattering on e- <0"
225 << G4endl;
226 G4cout << "y= " << y
227 << " e(MeV)= " << tkin << " p(MeV/c)= " << sqrt(mom2)
228 << " Z= " << targetZ << " "
229 << particle->GetParticleName() << G4endl;
230 G4cout << " 1-costm= " << 1.0-costm << " screenZ= " << screenZ
231 << " x= " << x << G4endl;
232 }
233 y = 0.0;
234 }
235 xSection = y;
236 }
237 /*
238 G4cout << "G4WentzelVI:XS per A " << " Z= " << targetZ
239 << " e(MeV)= " << tkin/MeV << " XSel= " << xSection
240 << " cut(MeV)= " << ecut/MeV
241 << " zmaxE= " << (1.0 - cosTetMaxElec)/screenZ
242 << " zmaxN= " << (1.0 - cosThetaMax)/screenZ
243 << " 1-costm= " << 1.0 - cosThetaMax << G4endl;
244 */
245 // scattering off nucleus
246 if(cosTMax < 1.0) {
247 x = (1.0 - cosTMax)/screenZ;
248 if(x < numlimit) {
249 x2 = 0.5*x*x;
250 y = x2*(1.0 - 1.3333333*x + 3*x2);
251 if(0.0 < factB) { y -= fb*x2*x*(0.6666667 - x); }
252 } else {
253 x1= x/(1 + x);
254 xlog = log(1.0 + x);
255 y = xlog - x1;
256 if(0.0 < factB) { y -= fb*(x + x1 - 2*xlog); }
257 }
258
259 if(y < 0.0) {
260 ++nwarnings;
261 if(nwarnings < nwarnlimit) {
262 G4cout << "G4WentzelVIRelXSection::ComputeTransportCrossSectionPerAtom scattering on e- <0"
263 << G4endl;
264 G4cout << "y= " << y
265 << " e(MeV)= " << tkin << " Z= " << targetZ << " "
266 << particle->GetParticleName() << G4endl;
267 G4cout << " formfactA= " << formfactA << " screenZ= " << screenZ
268 << " x= " << " x1= " << x1 <<G4endl;
269 }
270 y = 0.0;
271 }
272 xSection += y*targetZ;
273 }
274 xSection *= kinFactor;
275 /*
276 G4cout << "Z= " << targetZ << " XStot= " << xSection/barn
277 << " screenZ= " << screenZ << " formF= " << formfactA
278 << " for " << particle->GetParticleName()
279 << " m= " << mass << " 1/v= " << sqrt(invbeta2) << " p= " << sqrt(mom2)
280 << " x= " << x
281 << G4endl;
282 */
283 return xSection;
284}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const G4String & GetParticleName() const

Referenced by G4WentzelVIRelModel::ComputeCrossSectionPerAtom().

◆ GetCosThetaElec()

G4double G4WentzelVIRelXSection::GetCosThetaElec ( ) const
inline

Definition at line 224 of file G4WentzelVIRelXSection.hh.

225{
226 return cosTetMaxElec;
227}

◆ GetCosThetaNuc()

G4double G4WentzelVIRelXSection::GetCosThetaNuc ( ) const
inline

Definition at line 217 of file G4WentzelVIRelXSection.hh.

218{
219 return cosTetMaxNuc;
220}

◆ GetMomentumSquare()

G4double G4WentzelVIRelXSection::GetMomentumSquare ( ) const
inline

Definition at line 210 of file G4WentzelVIRelXSection.hh.

211{
212 return mom2;
213}

Referenced by G4hCoulombScatteringModel::SampleSecondaries().

◆ Initialise()

void G4WentzelVIRelXSection::Initialise ( const G4ParticleDefinition p,
G4double  CosThetaLim 
)

Definition at line 114 of file G4WentzelVIRelXSection.cc.

116{
117 SetupParticle(p);
118 tkin = mom2 = momCM2 = 0.0;
119 ecut = etag = DBL_MAX;
120 targetZ = 0;
121 cosThetaMax = CosThetaLim;
122 G4double a =
123 G4LossTableManager::Instance()->FactorForAngleLimit()*CLHEP::hbarc/CLHEP::fermi;
124 factorA2 = 0.5*a*a;
125 currentMaterial = 0;
126}
static G4LossTableManager * Instance()
G4double FactorForAngleLimit() const
void SetupParticle(const G4ParticleDefinition *)
#define DBL_MAX
Definition: templates.hh:83

Referenced by G4WentzelVIRelXSection(), G4hCoulombScatteringModel::Initialise(), and G4WentzelVIRelModel::Initialise().

◆ SampleSingleScattering()

G4ThreeVector G4WentzelVIRelXSection::SampleSingleScattering ( G4double  CosThetaMin,
G4double  CosThetaMax,
G4double  elecRatio = 0.0 
)

Definition at line 289 of file G4WentzelVIRelXSection.cc.

292{
293 G4ThreeVector v(0.0,0.0,1.0);
294
295 G4double formf = formfactA;
296 G4double cost1 = cosTMin;
297 G4double cost2 = cosTMax;
298 if(elecRatio > 0.0) {
299 if(G4UniformRand() <= elecRatio) {
300 formf = 0.0;
301 cost1 = std::max(cost1,cosTetMaxElec);
302 cost2 = std::max(cost2,cosTetMaxElec);
303 }
304 }
305 if(cost1 < cost2) { return v; }
306
307 G4double w1 = 1. - cost1 + screenZ;
308 G4double w2 = 1. - cost2 + screenZ;
309 G4double z1 = w1*w2/(w1 + G4UniformRand()*(w2 - w1)) - screenZ;
310
311 //G4double fm = 1.0 + formf*z1/(1.0 + (mass + tkin)*z1/targetMass);
312 G4double fm = 1.0 + formf*z1;
313 //G4double grej = (1. - z1*factB)/( (1.0 + z1*factD)*fm*fm );
314 G4double grej = (1. - z1*factB + factB1*targetZ*sqrt(z1*factB)*(2 - z1))/( (1.0 + z1*factD)*fm*fm );
315 // "false" scattering
316 if( G4UniformRand() > grej ) { return v; }
317 // }
318 G4double cost = 1.0 - z1;
319
320 if(cost > 1.0) { cost = 1.0; }
321 else if(cost < -1.0) { cost =-1.0; }
322 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
323 //G4cout << "sint= " << sint << G4endl;
324 G4double phi = twopi*G4UniformRand();
325 G4double vx1 = sint*cos(phi);
326 G4double vy1 = sint*sin(phi);
327
328 // only direction is changed
329 v.set(vx1,vy1,cost);
330 return v;
331}
#define G4UniformRand()
Definition: Randomize.hh:53

Referenced by G4WentzelVIRelModel::SampleScattering(), and G4hCoulombScatteringModel::SampleSecondaries().

◆ SetRelativisticMass()

void G4WentzelVIRelXSection::SetRelativisticMass ( G4double  value)
inline

Definition at line 203 of file G4WentzelVIRelXSection.hh.

204{
205 mass = value;
206}

◆ SetTargetMass()

void G4WentzelVIRelXSection::SetTargetMass ( G4double  value)
inline

Definition at line 195 of file G4WentzelVIRelXSection.hh.

196{
197 targetMass = value;
198 factD = std::sqrt(mom2)/value;
199}

Referenced by G4hCoulombScatteringModel::SampleSecondaries(), and SetupTarget().

◆ SetupKinematic()

G4double G4WentzelVIRelXSection::SetupKinematic ( G4double  kinEnergy,
const G4Material mat 
)
inline

Definition at line 179 of file G4WentzelVIRelXSection.hh.

180{
181 if(ekin != tkin || mat != currentMaterial) {
182 currentMaterial = mat;
183 tkin = ekin;
184 mom2 = tkin*(tkin + 2.0*mass);
185 invbeta2 = 1.0 + mass*mass/mom2;
186 factB = spin/invbeta2;
187 cosTetMaxNuc =
188 std::max(cosThetaMax,1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
189 }
190 return cosTetMaxNuc;
191}
G4double GetInvA23() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225

Referenced by G4hCoulombScatteringModel::ComputeCrossSectionPerAtom(), G4WentzelVIRelModel::ComputeCrossSectionPerAtom(), G4WentzelVIRelModel::ComputeGeomPathLength(), G4WentzelVIRelModel::ComputeTruePathLengthLimit(), and G4WentzelVIRelModel::ComputeTrueStepLength().

◆ SetupParticle()

void G4WentzelVIRelXSection::SetupParticle ( const G4ParticleDefinition p)

Definition at line 130 of file G4WentzelVIRelXSection.cc.

131{
132 particle = p;
133 mass = particle->GetPDGMass();
134 spin = particle->GetPDGSpin();
135 if(0.0 != spin) { spin = 0.5; }
136 G4double q = std::fabs(particle->GetPDGCharge()/eplus);
137 chargeSquare = q*q;
138 charge3 = chargeSquare*q;
139 tkin = 0.0;
140 currentMaterial = 0;
141 targetZ = 0;
142}
G4double GetPDGCharge() const

Referenced by Initialise(), and G4hCoulombScatteringModel::SetupParticle().

◆ SetupTarget()

G4double G4WentzelVIRelXSection::SetupTarget ( G4int  Z,
G4double  cut = DBL_MAX 
)

Definition at line 147 of file G4WentzelVIRelXSection.cc.

148{
149 G4double cosTetMaxNuc2 = cosTetMaxNuc;
150 if(Z != targetZ || tkin != etag) {
151 etag = tkin;
152 targetZ = Z;
153 if(targetZ > 99) { targetZ = 99; }
154 SetTargetMass(fNistManager->GetAtomicMassAmu(targetZ)*CLHEP::amu_c2);
155 //G4double tmass2 = targetMass*targetMass;
156 //G4double etot = tkin + mass;
157 //G4double invmass2 = mass*mass + tmass2 + 2*etot*targetMass;
158 //momCM2 = mom2*tmass2/invmass2;
159 //gam0pcmp = (etot + targetMass)*targetMass/invmass2;
160 //pcmp2 = tmass2/invmass2;
161
162 kinFactor = coeff*targetZ*chargeSquare*invbeta2/mom2;
163
164 screenZ = ScreenRSquare[targetZ]/mom2;
165 if(Z > 1) {
166 screenZ *= std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare);
167 /*
168 if(mass > MeV) {
169 screenZ *= std::min(Z*1.13,1.13 +3.76*Z*Z*invbeta2*alpha2*chargeSquare);
170 } else {
171 G4double tau = tkin/mass;
172 // screenZ *= std::min(Z*invbeta2,
173 screenZ *= std::min(Z*1.13,
174 (1.13 +3.76*Z*Z*invbeta2*alpha2*std::sqrt(tau/(tau + fG4pow->Z23(Z)))));
175 }
176 */
177 }
178 if(targetZ == 1 && cosTetMaxNuc2 < 0.0 && particle == theProton) {
179 cosTetMaxNuc2 = 0.0;
180 }
181 formfactA = FormFactor[targetZ]*mom2;
182
183 cosTetMaxElec = 1.0;
184 ComputeMaxElectronScattering(cut);
185 }
186 return cosTetMaxNuc2;
187}
G4double GetAtomicMassAmu(const G4String &symb) const
void SetTargetMass(G4double value)

Referenced by G4hCoulombScatteringModel::ComputeCrossSectionPerAtom(), G4WentzelVIRelModel::ComputeCrossSectionPerAtom(), and G4WentzelVIRelModel::SampleScattering().


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