source: trunk/source/processes/electromagnetic/standard/include/G4eCoulombScatteringModel.hh@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 8.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4eCoulombScatteringModel.hh,v 1.36 2008/08/04 08:49:09 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4eCoulombScatteringModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 19.02.2006
39//
40// Modifications:
41// 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
42// logic of building - only elements from G4ElementTable
43// 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
44// 19.08.06 V.Ivanchenko add inline function ScreeningParameter and
45// make some members protected
46// 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e-
47// 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion
48//
49// Class Description:
50//
51// Implementation of eCoulombScattering of pointlike charge particle
52// on Atomic Nucleus for interval of scattering anles in Lab system
53// thetaMin - ThetaMax, nucleus recoil is neglected.
54// The model based on analysis of J.M.Fernandez-Varea et al.
55// NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with
56// screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256.
57//
58
59// -------------------------------------------------------------------
60//
61
62#ifndef G4eCoulombScatteringModel_h
63#define G4eCoulombScatteringModel_h 1
64
65#include "G4VEmModel.hh"
66#include "G4PhysicsTable.hh"
67#include "G4NistManager.hh"
68#include "globals.hh"
69
70class G4ParticleChangeForGamma;
71class G4ParticleDefinition;
72
73class G4eCoulombScatteringModel : public G4VEmModel
74{
75
76public:
77
78 G4eCoulombScatteringModel(const G4String& nam = "eCoulombScattering");
79
80 virtual ~G4eCoulombScatteringModel();
81
82 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83
84 virtual G4double ComputeCrossSectionPerAtom(
85 const G4ParticleDefinition*,
86 G4double kinEnergy,
87 G4double Z,
88 G4double A,
89 G4double cut,
90 G4double emax);
91
92 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
93 const G4MaterialCutsCouple*,
94 const G4DynamicParticle*,
95 G4double tmin,
96 G4double maxEnergy);
97
98 inline void SetRecoilThreshold(G4double eth);
99
100protected:
101
102 G4double CrossSectionPerAtom();
103
104 G4double SampleCosineTheta();
105
106 inline void DefineMaterial(const G4MaterialCutsCouple*);
107
108 inline void SetupParticle(const G4ParticleDefinition*);
109
110 inline void SetupKinematic(G4double kinEnergy, G4double cut);
111
112 inline void SetupTarget(G4double Z, G4double kinEnergy);
113
114private:
115
116 void ComputeMaxElectronScattering(G4double cut);
117
118 // hide assignment operator
119 G4eCoulombScatteringModel & operator=(const G4eCoulombScatteringModel &right);
120 G4eCoulombScatteringModel(const G4eCoulombScatteringModel&);
121
122protected:
123
124 const G4ParticleDefinition* theProton;
125 const G4ParticleDefinition* theElectron;
126 const G4ParticleDefinition* thePositron;
127
128 G4ParticleTable* theParticleTable;
129 G4ParticleChangeForGamma* fParticleChange;
130 G4NistManager* fNistManager;
131 const G4DataVector* currentCuts;
132
133 const G4MaterialCutsCouple* currentCouple;
134 const G4Material* currentMaterial;
135 const G4Element* currentElement;
136 G4int currentMaterialIndex;
137
138 G4double coeff;
139 G4double constn;
140 G4double cosThetaMin;
141 G4double cosThetaMax;
142 G4double cosTetMinNuc;
143 G4double cosTetMaxNuc;
144 G4double cosTetMaxNuc2;
145 G4double cosTetMaxElec;
146 G4double cosTetMaxElec2;
147 G4double q2Limit;
148 G4double recoilThreshold;
149 G4double elecXSection;
150 G4double nucXSection;
151 G4double ecut;
152
153 // projectile
154 const G4ParticleDefinition* particle;
155
156 G4double chargeSquare;
157 G4double spin;
158 G4double mass;
159 G4double tkin;
160 G4double mom2;
161 G4double invbeta2;
162 G4double etag;
163 G4double lowEnergyLimit;
164
165 // target
166 G4double targetZ;
167 G4double screenZ;
168 G4double formfactA;
169 G4int idxelm;
170
171private:
172
173 G4double a0;
174 G4double alpha2;
175 G4double faclim;
176 G4double FF[100];
177
178 G4bool isInitialised;
179};
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182
183inline
184void G4eCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
185{
186 if(cup != currentCouple) {
187 currentCouple = cup;
188 currentMaterial = cup->GetMaterial();
189 currentMaterialIndex = currentCouple->GetIndex();
190 }
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
195inline
196void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
197{
198 // Initialise mass and charge
199 if(p != particle) {
200 particle = p;
201 mass = particle->GetPDGMass();
202 spin = particle->GetPDGSpin();
203 G4double q = particle->GetPDGCharge()/eplus;
204 chargeSquare = q*q;
205 tkin = 0.0;
206 lowEnergyLimit = keV*mass/electron_mass_c2;
207 }
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
212inline void G4eCoulombScatteringModel::SetupKinematic(G4double ekin,
213 G4double cut)
214{
215 if(ekin != tkin || ecut != cut) {
216 tkin = ekin;
217 mom2 = tkin*(tkin + 2.0*mass);
218 invbeta2 = 1.0 + mass*mass/mom2;
219 cosTetMinNuc = cosThetaMin;
220 cosTetMaxNuc = cosThetaMax;
221 if(ekin <= 10.*cut && mass < MeV && cosThetaMin < 1.0) {
222 cosTetMinNuc = ekin*(cosThetaMin + 1.0)/(10.*cut) - 1.0;
223 }
224 ComputeMaxElectronScattering(cut);
225 }
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
230inline void G4eCoulombScatteringModel::SetupTarget(G4double Z, G4double e)
231{
232 if(Z != targetZ || e != etag) {
233 etag = e;
234 targetZ = Z;
235 G4int iz= G4int(Z);
236 if(iz > 99) iz = 99;
237 G4double x = fNistManager->GetZ13(iz);
238 screenZ = a0*x*x/mom2;
239 if(iz > 1) screenZ *=(1.13 + 3.76*invbeta2*Z*Z*chargeSquare*alpha2);
240 //screenZ = a0*x*x*(1.13 + 3.76*Z*Z*chargeSquare*alpha2)/mom2;
241 // A.V. Butkevich et al., NIM A 488 (2002) 282
242 formfactA = FF[iz];
243 if(formfactA == 0.0) {
244 x = fNistManager->GetA27(iz);
245 formfactA = constn*x*x;
246 FF[iz] = formfactA;
247 }
248 formfactA *= mom2;
249 cosTetMaxNuc2 = cosTetMaxNuc;
250 if(particle == theProton && 1 == iz && cosTetMaxNuc2 < 0.0) {
251 cosTetMaxNuc2 = 0.0;
252 }
253 /*
254 G4double ee = 10.*eV*Z;
255 if(1 == iz) ee *= 2.0;
256 G4double z = std::min(cosTetMaxElec, 1.0 - std::max(ecut,ee)*amu_c2
257 *fNistManager->GetAtomicMassAmu(iz)/mom2);
258 cosTetMaxElec2 = std::max(cosTetMaxNuc2, z);
259 */
260 cosTetMaxElec2 = cosTetMaxElec;
261 }
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
265
266inline void G4eCoulombScatteringModel::SetRecoilThreshold(G4double eth)
267{
268 recoilThreshold = eth;
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272
273#endif
Note: See TracBrowser for help on using the repository browser.