source: trunk/examples/extended/electromagnetic/TestEm7/include/G4ScreenedNuclearRecoil.hh@ 1036

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

update

File size: 11.4 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//
27// $Id: G4ScreenedNuclearRecoil.hh,v 1.3 2007/12/07 17:51:10 vnivanch Exp $
28// GEANT4 tag $Name: $
29//
30//
31//
32// Class Description
33// Process for screened electromagnetic nuclear elastic scattering;
34// Physics comes from:
35// Marcus H. Mendenhall and Robert A. Weller,
36// "Algorithms for the rapid computation of classical cross sections
37// for screened Coulomb collisions "
38// Nuclear Instruments and Methods in Physics Research B58 (1991) 11-17
39// The only input required is a screening function phi(r/a) which is the ratio
40// of the actual interatomic potential for two atoms with atomic numbers Z1 and Z2,
41// to the unscreened potential Z1*Z2*e^2/r where e^2 is elm_coupling in Geant4 units
42// the actual screening tables are computed externally in a python module "screened_scattering.py"
43// to allow very specific screening functions to be added if desired, without messing
44// with the insides of this code.
45//
46// First version, April 2004, Marcus H. Mendenhall, Vanderbilt University
47//
48// Class Description - End
49
50
51#ifndef G4ScreenedNuclearRecoil_h
52#define G4ScreenedNuclearRecoil_h 1
53
54#include "globals.hh"
55#include "G4VDiscreteProcess.hh"
56#include "c2_function.hh"
57
58#include <map>
59#include <vector>
60
61class G4VParticleChange;
62
63typedef struct G4ScreeningTables {
64 G4double z1, z2, m1, m2, au, emin;
65 c2_function<G4double> *EMphiData;
66} G4ScreeningTables;
67
68// A class for loading ScreenedCoulombCrossSections
69class G4ScreenedCoulombCrossSectionInfo
70{
71public:
72 G4ScreenedCoulombCrossSectionInfo() { }
73 ~G4ScreenedCoulombCrossSectionInfo() { }
74
75 const char *CVSHeaderVers() { return
76 "";
77 }
78
79 const char *CVSFileVers() { return ""; }
80};
81
82// A class for loading ScreenedCoulombCrossSections
83class G4ScreenedCoulombCrossSection : public G4ScreenedCoulombCrossSectionInfo
84{
85public:
86
87 G4ScreenedCoulombCrossSection() : verbosity(1) { }
88 G4ScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src) :
89 G4ScreenedCoulombCrossSectionInfo(),verbosity(src.verbosity) { }
90 virtual ~G4ScreenedCoulombCrossSection();
91
92 typedef std::map<G4int, G4ScreeningTables> ScreeningMap;
93
94 // a local, fast-access mapping of a particle's Z to its full definition
95 typedef std::map<G4int, class G4ParticleDefinition *> ParticleCache;
96
97 // LoadData is called by G4ScreenedNuclearRecoil::GetMeanFreePath
98 // It loads the data tables, builds the elemental cross-section tables.
99 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff) = 0;
100
101 // BuildMFPTables is called by G4ScreenedNuclearRecoil::GetMeanFreePath to build the MFP tables for each material
102 void BuildMFPTables(void); // scan the MaterialsTable and construct MFP tables
103
104 virtual G4ScreenedCoulombCrossSection *create() = 0; // a 'virtual constructor' which clones the class
105 const G4ScreeningTables *GetScreening(G4int Z) { return &(screeningData[Z]); }
106 void SetVerbosity(G4int v) { verbosity=v; }
107
108 // this process needs element selection weighted only by number density
109 G4ParticleDefinition* SelectRandomUnweightedTarget(const G4MaterialCutsCouple* couple);
110
111 enum { nMassMapElements=116 };
112
113 G4double standardmass(G4int z1) { return z1 <= nMassMapElements ? massmap[z1] : 2.5*z1; }
114
115 // get the mean-free-path table for the indexed material
116 c2_function<G4double> * operator [] (G4int materialIndex) {
117 return MFPTables.find(materialIndex)!=MFPTables.end() ? MFPTables[materialIndex] : (c2_function<G4double> *)0;
118 }
119
120protected:
121 ScreeningMap screeningData; // screening tables for each element
122 ParticleCache targetMap;
123 G4int verbosity;
124 std::map<G4int, c2_function<G4double> *> sigmaMap; // total cross section for each element
125 std::map<G4int, c2_function<G4double> *> MFPTables; // MFP for each material
126
127private:
128 static const G4double massmap[nMassMapElements+1];
129
130};
131
132typedef struct G4CoulombKinematicsInfo {
133 G4double impactParameter;
134 G4ScreenedCoulombCrossSection *crossSection;
135 G4double a1, a2, sinTheta, cosTheta, sinZeta, cosZeta, eRecoil;
136 G4ParticleDefinition *recoilIon; } G4CoulombKinematicsInfo;
137
138class G4ScreenedCollisionStage {
139public:
140 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
141 const class G4Track& aTrack, const class G4Step& aStep)=0;
142 virtual ~G4ScreenedCollisionStage() {}
143};
144
145class G4ScreenedCoulombClassicalKinematics: public G4ScreenedCoulombCrossSectionInfo, public G4ScreenedCollisionStage {
146
147public:
148 G4ScreenedCoulombClassicalKinematics() { }
149 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
150 const class G4Track& aTrack, const class G4Step& aStep);
151
152 G4bool DoScreeningComputation(class G4ScreenedNuclearRecoil *master,
153 const G4ScreeningTables *screen,
154 G4double eps, G4double beta);
155 virtual ~G4ScreenedCoulombClassicalKinematics() {}
156};
157
158class G4SingleScatter: public G4ScreenedCoulombCrossSectionInfo, public G4ScreenedCollisionStage {
159
160public:
161 G4SingleScatter() { }
162 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
163 const class G4Track& aTrack, const class G4Step& aStep);
164 virtual ~G4SingleScatter() {}
165};
166
167class G4ScreenedNuclearRecoil : public G4ScreenedCoulombCrossSectionInfo, public G4VDiscreteProcess
168{
169public:
170
171 friend class G4ScreenedCollisionStage;
172
173 G4ScreenedNuclearRecoil(const G4String& processName = "ScreenedElastic",
174 const G4String &ScreeningKey="zbl", G4bool GenerateRecoils=1,
175 G4double RecoilCutoff=100.0*eV, G4double PhysicsCutoff=10.0*eV);
176
177 virtual ~G4ScreenedNuclearRecoil();
178
179 virtual G4double GetMeanFreePath(const G4Track&, G4double, G4ForceCondition* );
180
181 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
182
183 virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
184
185 virtual void BuildPhysicsTable(const G4ParticleDefinition&) { }
186
187 virtual void DumpPhysicsTable(const G4ParticleDefinition& aParticleType);
188
189 virtual G4bool CheckNuclearCollision(G4double A, G4double A1, G4double apsis); // return true if hard collision
190
191 virtual G4ScreenedCoulombCrossSection *GetNewCrossSectionHandler(void);
192
193 G4double GetNIEL() const { return NIEL; } // Get non-ionizing energy loss for last step
194
195 void ResetTables(); // clear all data tables to allow changing energy cutoff, materials, etc.
196
197 std::string GetScreeningKey() const { return screeningKey; }
198 void AllowEnergyDeposition(G4bool flag) { registerDepositedEnergy=flag; }
199 G4bool GetAllowEnergyDeposition() const { return registerDepositedEnergy; }
200 void EnableRecoils(G4bool flag) { generateRecoils=flag; }
201 G4bool GetEnableRecoils() const { return generateRecoils; }
202 void SetMFPScaling(G4double scale) { MFPScale=scale; }
203 G4double GetMFPScaling() const { return MFPScale; }
204 void AvoidNuclearReactions(G4bool flag) { avoidReactions=flag; }
205 G4bool GetAvoidNuclearReactions() const { return avoidReactions; }
206 void SetRecoilCutoff(G4double energy) { recoilCutoff=energy; }
207 G4double GetRecoilCutoff() const { return recoilCutoff; }
208 void SetPhysicsCutoff(G4double energy) { physicsCutoff=energy; ResetTables(); }
209 G4double GetPhysicsCutoff() const { return physicsCutoff; }
210 class G4ParticleChange &GetParticleChange() { return aParticleChange; }
211 void AddToNIEL(G4double energy) { NIEL+=energy; }
212 void SetCrossSectionHardening(G4double fraction, G4double HardeningFactor) {
213 hardeningFraction=fraction;
214 hardeningFactor=HardeningFactor;
215 }
216 G4double GetHardeningFraction() const { return hardeningFraction; }
217 G4double GetHardeningFactor() const { return hardeningFactor; }
218 G4double GetCurrentInteractionLength() const { return currentInteractionLength; }
219 void SetExternalCrossSectionHandler(G4ScreenedCoulombCrossSection *cs) {
220 externalCrossSectionConstructor=cs;
221 }
222 G4int GetVerboseLevel() const { return verboseLevel; }
223 std::map<G4int, G4ScreenedCoulombCrossSection*> &GetCrossSectionHandlers()
224 { return crossSectionHandlers; }
225 void ClearStages(void);
226 void AddStage(G4ScreenedCollisionStage *stage) { collisionStages.push_back(stage); }
227 G4CoulombKinematicsInfo &GetKinematics() { return kinematics; }
228 void SetValidCollision(G4bool flag) { validCollision=flag; }
229 G4bool GetValidCollision() const { return validCollision; }
230
231protected:
232 G4double highEnergyLimit, lowEnergyLimit;
233 G4String screeningKey;
234 G4bool generateRecoils, avoidReactions;
235 G4double recoilCutoff, physicsCutoff;
236 G4bool registerDepositedEnergy;
237 G4double NIEL;
238 G4double MFPScale;
239 G4double hardeningFraction, hardeningFactor;
240
241 G4ScreenedCoulombCrossSection *externalCrossSectionConstructor;
242 std::vector<G4ScreenedCollisionStage *> collisionStages;
243
244 std::map<G4int, G4ScreenedCoulombCrossSection*> crossSectionHandlers;
245 std::map<G4int, c2_function<G4double>*> meanFreePathTables;
246
247 G4bool validCollision;
248 G4CoulombKinematicsInfo kinematics;
249
250};
251
252// A customized G4CrossSectionHandler which gets its data from an external program
253class G4NativeScreenedCoulombCrossSection: public G4ScreenedCoulombCrossSection
254{
255public:
256 G4NativeScreenedCoulombCrossSection();
257
258 G4NativeScreenedCoulombCrossSection(const G4NativeScreenedCoulombCrossSection &src)
259 : G4ScreenedCoulombCrossSection(src), phiMap(src.phiMap) { }
260
261 G4NativeScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src) : G4ScreenedCoulombCrossSection(src) { }
262 virtual ~G4NativeScreenedCoulombCrossSection();
263
264 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff);
265 virtual G4ScreenedCoulombCrossSection *create()
266 { return new G4NativeScreenedCoulombCrossSection(*this); }
267 // get a list of available keys
268 std::vector<G4String> GetScreeningKeys() const;
269
270 typedef c2_function<G4double> &(*ScreeningFunc)(G4int z1, G4int z2, size_t nPoints, G4double rMax, G4double *au);
271
272 void AddScreeningFunction(G4String name, ScreeningFunc fn) {
273 phiMap[name]=fn;
274 }
275
276private:
277 // this is a map used to look up screening function generators
278 std::map<std::string, ScreeningFunc> phiMap;
279};
280
281
282
283#endif
Note: See TracBrowser for help on using the repository browser.