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

Last change on this file since 1268 was 1230, checked in by garnier, 16 years ago

update to geant4.9.3

File size: 17.7 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// G4ScreenedNuclearRecoil.hh,v 1.24 2008/05/01 19:58:59 marcus Exp
28// GEANT4 tag
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// May 1, 2008 -- Added code to allow process to have zero cross section above max energy, to coordinate with G4MSC. -- mhm
48//
49// Class Description - End
50
51
52#ifndef G4ScreenedNuclearRecoil_h
53#define G4ScreenedNuclearRecoil_h 1
54
55#include "globals.hh"
56#include "G4VDiscreteProcess.hh"
57#include "G4ParticleChange.hh"
58#include "c2_function.hh"
59
60#include <map>
61#include <vector>
62
63class G4VNIELPartition;
64
65typedef c2_const_ptr<G4double> G4_c2_const_ptr;
66typedef c2_ptr<G4double> G4_c2_ptr;
67typedef c2_function<G4double> G4_c2_function;
68
69typedef struct G4ScreeningTables {
70 G4double z1, z2, m1, m2, au, emin;
71 G4_c2_const_ptr EMphiData;
72} G4ScreeningTables;
73
74// A class for loading ScreenedCoulombCrossSections
75class G4ScreenedCoulombCrossSectionInfo
76{
77public:
78 G4ScreenedCoulombCrossSectionInfo() { }
79 ~G4ScreenedCoulombCrossSectionInfo() { }
80
81 static const char* CVSHeaderVers() { return
82 "G4ScreenedNuclearRecoil.hh,v 1.24 2008/05/01 19:58:59 marcus Exp GEANT4 tag ";
83 }
84 static const char* CVSFileVers();
85};
86
87// A class for loading ScreenedCoulombCrossSections
88class G4ScreenedCoulombCrossSection : public G4ScreenedCoulombCrossSectionInfo
89{
90public:
91
92 G4ScreenedCoulombCrossSection() : verbosity(1) { }
93 G4ScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src) :
94 G4ScreenedCoulombCrossSectionInfo(),verbosity(src.verbosity) { }
95 virtual ~G4ScreenedCoulombCrossSection();
96
97 typedef std::map<G4int, G4ScreeningTables> ScreeningMap;
98
99 // a local, fast-access mapping of a particle's Z to its full definition
100 typedef std::map<G4int, class G4ParticleDefinition *> ParticleCache;
101
102 // LoadData is called by G4ScreenedNuclearRecoil::GetMeanFreePath
103 // It loads the data tables, builds the elemental cross-section tables.
104 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff) = 0;
105
106 // BuildMFPTables is called by G4ScreenedNuclearRecoil::GetMeanFreePath to build the MFP tables for each material
107 void BuildMFPTables(void); // scan the MaterialsTable and construct MFP tables
108
109 virtual G4ScreenedCoulombCrossSection *create() = 0; // a 'virtual constructor' which clones the class
110 const G4ScreeningTables *GetScreening(G4int Z) { return &(screeningData[Z]); }
111 void SetVerbosity(G4int v) { verbosity=v; }
112
113 // this process needs element selection weighted only by number density
114 G4ParticleDefinition* SelectRandomUnweightedTarget(const G4MaterialCutsCouple* couple);
115
116 enum { nMassMapElements=116 };
117
118 G4double standardmass(G4int z1) { return z1 <= nMassMapElements ? massmap[z1] : 2.5*z1; }
119
120 // get the mean-free-path table for the indexed material
121 const G4_c2_function * operator [] (G4int materialIndex) {
122 return MFPTables.find(materialIndex)!=MFPTables.end() ? &(MFPTables[materialIndex].get()) : (G4_c2_function *)0;
123 }
124
125protected:
126 ScreeningMap screeningData; // screening tables for each element
127 ParticleCache targetMap;
128 G4int verbosity;
129 std::map<G4int, G4_c2_const_ptr > sigmaMap; // total cross section for each element
130 std::map<G4int, G4_c2_const_ptr > MFPTables; // MFP for each material
131
132private:
133 static const G4double massmap[nMassMapElements+1];
134
135};
136
137typedef struct G4CoulombKinematicsInfo {
138 G4double impactParameter;
139 G4ScreenedCoulombCrossSection *crossSection;
140 G4double a1, a2, sinTheta, cosTheta, sinZeta, cosZeta, eRecoil;
141 G4ParticleDefinition *recoilIon;
142 const G4Material *targetMaterial;
143} G4CoulombKinematicsInfo;
144
145class G4ScreenedCollisionStage {
146public:
147 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
148 const class G4Track& aTrack, const class G4Step& aStep)=0;
149 virtual ~G4ScreenedCollisionStage() {}
150};
151
152class G4ScreenedCoulombClassicalKinematics: public G4ScreenedCoulombCrossSectionInfo, public G4ScreenedCollisionStage {
153
154public:
155 G4ScreenedCoulombClassicalKinematics();
156 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
157 const class G4Track& aTrack, const class G4Step& aStep);
158
159 G4bool DoScreeningComputation(class G4ScreenedNuclearRecoil *master,
160 const G4ScreeningTables *screen,
161 G4double eps, G4double beta);
162 virtual ~G4ScreenedCoulombClassicalKinematics() { }
163protected:
164 // the c2_functions we need to do the work.
165 c2_const_plugin_function_p<G4double> &phifunc;
166 c2_linear_p<G4double> &xovereps;
167 G4_c2_ptr diff;
168
169};
170
171class G4SingleScatter: public G4ScreenedCoulombCrossSectionInfo, public G4ScreenedCollisionStage {
172
173public:
174 G4SingleScatter() { }
175 virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
176 const class G4Track& aTrack, const class G4Step& aStep);
177 virtual ~G4SingleScatter() {}
178};
179
180/**
181 \brief A process which handles screened Coulomb collisions between nuclei
182
183*/
184
185class G4ScreenedNuclearRecoil : public G4ScreenedCoulombCrossSectionInfo, public G4VDiscreteProcess
186{
187public:
188
189 friend class G4ScreenedCollisionStage;
190
191 /// \brief Construct the process and set some physics parameters for it.
192 /// \param processName the name to assign the process
193 /// \param ScreeningKey the name of a screening function to use.
194 /// The default functions are "zbl" (recommended for soft scattering),
195 /// "lj" (recommended for backscattering) and "mol" (Moliere potential)
196 /// \param GenerateRecoils if frue, ions struck by primary are converted into new moving particles.
197 /// If false, energy is deposited, but no new moving ions are created.
198 /// \param RecoilCutoff energy below which no new moving particles will be created,
199 /// even if \a GenerateRecoils is true.
200 /// Also, a moving primary particle will be stopped if its energy falls below this limit.
201 /// \param PhysicsCutoff the energy transfer to which screening tables are calucalted.
202 /// There is no really
203 /// compelling reason to change it from the 10.0 eV default. However, see the paper on running this
204 /// in thin targets for further discussion, and its interaction with SetMFPScaling()
205 G4ScreenedNuclearRecoil(const G4String& processName = "ScreenedElastic",
206 const G4String &ScreeningKey="zbl", G4bool GenerateRecoils=1,
207 G4double RecoilCutoff=100.0*eV, G4double PhysicsCutoff=10.0*eV);
208 /// \brief destructor
209 virtual ~G4ScreenedNuclearRecoil();
210 /// \brief used internally by Geant4 machinery
211 virtual G4double GetMeanFreePath(const G4Track&, G4double, G4ForceCondition* );
212 /// \brief used internally by Geant4 machinery
213 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
214 /// \brief test if a prticle of type \a aParticleType can use this process
215 /// \param aParticleType the particle to test
216 virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
217 /// \brief Build physics tables in advance. Not Implemented.
218 /// \param aParticleType the type of particle to build tables for
219 virtual void BuildPhysicsTable(const G4ParticleDefinition&) { }
220 /// \brief Export physics tables for persistency. Not Implemented.
221 /// \param aParticleType the type of particle to build tables for
222 virtual void DumpPhysicsTable(const G4ParticleDefinition& aParticleType);
223 /// \brief deterine if the moving particle is within the strong force range of the selected nucleus
224 /// \param A the nucleon number of the beam
225 /// \param A1 the nucleon number of the target
226 /// \param apsis the distance of closest approach
227 virtual G4bool CheckNuclearCollision(G4double A, G4double A1, G4double apsis); // return true if hard collision
228
229 virtual G4ScreenedCoulombCrossSection *GetNewCrossSectionHandler(void);
230
231 /// \brief Get non-ionizing energy loss for last step
232 G4double GetNIEL() const { return NIEL; }
233
234 /// \brief clear precomputed screening tables
235 void ResetTables(); // clear all data tables to allow changing energy cutoff, materials, etc.
236
237 /// \brief set the upper energy beyond which this process has no cross section
238 ///
239 /// This funciton is used to coordinate this process with G4MSC. Typically, G4MSC should
240 /// not be allowed to operate in a range which overlaps that of this process. The criterion which is most reasonable
241 /// is that the transition should be somewhere in the modestly relativistic regime (500 MeV/u for example).
242 /// \param energy energy per nucleon for the cutoff
243 void SetMaxEnergyForScattering(G4double energy) { processMaxEnergy=energy; }
244 /// \brief find out what screening funciton we are using
245 std::string GetScreeningKey() const { return screeningKey; }
246 /// \brief enable or disable all energy deposition by this process
247 /// \param flag if true, enable deposition of energy (the default). If false, disable deposition.
248 void AllowEnergyDeposition(G4bool flag) { registerDepositedEnergy=flag; }
249 /// \brief get flag indicating whether deposition is enabled
250 G4bool GetAllowEnergyDeposition() const { return registerDepositedEnergy; }
251 /// \brief enable or disable the generation of recoils.
252 /// If recoils are disabled, the energy they would have received is just deposited.
253 /// \param flag if true, create recoil ions in cases in which the energy is above the recoilCutoff.
254 /// If false, just deposit the energy.
255 void EnableRecoils(G4bool flag) { generateRecoils=flag; }
256 /// \brief find out if generation of recoils is enabled.
257 G4bool GetEnableRecoils() const { return generateRecoils; }
258 /// \brief set the mean free path scaling as specified
259 /// \param scale the factor by which the default MFP will be scaled.
260 /// Set to less than 1 for very thin films, typically, to sample multiple scattering,
261 /// or to greater than 1 for quick simulaitons with a very long flight path.
262 void SetMFPScaling(G4double scale) { MFPScale=scale; }
263 /// \brief get the MFPScaling parameter
264 G4double GetMFPScaling() const { return MFPScale; }
265 /// \brief enable or disable whether this process will skip collisions
266 /// which are close enough they need hadronic phsyics. Default is true (skip close collisions).
267 /// Disabling this results in excess nuclear stopping power.
268 /// \param flag true results in hard collisions being skipped. false allows hard collisions.
269 void AvoidNuclearReactions(G4bool flag) { avoidReactions=flag; }
270 /// \brief get the flag indicating whether hadronic collisions are ignored.
271 G4bool GetAvoidNuclearReactions() const { return avoidReactions; }
272 /// \brief set the minimum energy (per nucleon) at which recoils can be generated,
273 /// and the energy (per nucleon) below which all ions are stopped.
274 /// \param energy energy per nucleon
275 void SetRecoilCutoff(G4double energy) { recoilCutoff=energy; }
276 /// \brief get the recoil cutoff
277 G4double GetRecoilCutoff() const { return recoilCutoff; }
278 /// \brief set the energy to which screening tables are computed. Typically, this is 10 eV or so, and not often changed.
279 /// \param energy the cutoff energy
280 void SetPhysicsCutoff(G4double energy) { physicsCutoff=energy; ResetTables(); }
281 /// \brief get the physics cutoff energy.
282 G4double GetPhysicsCutoff() const { return physicsCutoff; }
283 /// \brief set the pointer to a class for paritioning energy into NIEL
284 /// \brief part the pointer to the class.
285 void SetNIELPartitionFunction(const G4VNIELPartition *part);
286 /// \brief set the cross section boost to provide faster computation of backscattering
287 /// \param fraction the fraction of particles to have their cross section boosted.
288 /// \param HardeningFactor the factor by which to boost the scattering cross section.
289 void SetCrossSectionHardening(G4double fraction, G4double HardeningFactor) {
290 hardeningFraction=fraction;
291 hardeningFactor=HardeningFactor;
292 }
293 /// \brief get the fraction of particles which will have boosted scattering
294 G4double GetHardeningFraction() const { return hardeningFraction; }
295 /// \brief get the boost factor in use.
296 G4double GetHardeningFactor() const { return hardeningFactor; }
297 /// \brief the the interaciton length used in the last scattering.
298 G4double GetCurrentInteractionLength() const { return currentInteractionLength; }
299 /// \brief set a function to compute screening tables, if the user needs non-standard behavior.
300 /// \param cs a class which constructs the screening tables.
301 void SetExternalCrossSectionHandler(G4ScreenedCoulombCrossSection *cs) {
302 externalCrossSectionConstructor=cs;
303 }
304 /// \brief get the verbosity.
305 G4int GetVerboseLevel() const { return verboseLevel; }
306
307 std::map<G4int, G4ScreenedCoulombCrossSection*> &GetCrossSectionHandlers()
308 { return crossSectionHandlers; }
309 void ClearStages(void);
310 void AddStage(G4ScreenedCollisionStage *stage) { collisionStages.push_back(stage); }
311 G4CoulombKinematicsInfo &GetKinematics() { return kinematics; }
312 void SetValidCollision(G4bool flag) { validCollision=flag; }
313 G4bool GetValidCollision() const { return validCollision; }
314
315 /// \brief get the pointer to our ParticleChange object. for internal use, primarily.
316 class G4ParticleChange &GetParticleChange() { return static_cast<G4ParticleChange &>(*pParticleChange); }
317 /// \brief take the given energy, and use the material information to partition it into NIEL and ionizing energy.
318 void DepositEnergy(G4int z1, G4double a1, const G4Material *material, G4double energy);
319
320protected:
321 /// \brief the energy per nucleon above which the MFP is constant
322 G4double highEnergyLimit;
323 /// \brief the energy per nucleon below which the MFP is zero
324 G4double lowEnergyLimit;
325 /// \brief the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC
326 G4double processMaxEnergy;
327 G4String screeningKey;
328 G4bool generateRecoils, avoidReactions;
329 G4double recoilCutoff, physicsCutoff;
330 G4bool registerDepositedEnergy;
331 G4double IonizingLoss, NIEL;
332 G4double MFPScale;
333 G4double hardeningFraction, hardeningFactor;
334
335 G4ScreenedCoulombCrossSection *externalCrossSectionConstructor;
336 std::vector<G4ScreenedCollisionStage *> collisionStages;
337
338 std::map<G4int, G4ScreenedCoulombCrossSection*> crossSectionHandlers;
339
340 G4bool validCollision;
341 G4CoulombKinematicsInfo kinematics;
342 const G4VNIELPartition *NIELPartitionFunction;
343};
344
345// A customized G4CrossSectionHandler which gets its data from an external program
346class G4NativeScreenedCoulombCrossSection: public G4ScreenedCoulombCrossSection
347{
348public:
349 G4NativeScreenedCoulombCrossSection();
350
351 G4NativeScreenedCoulombCrossSection(const G4NativeScreenedCoulombCrossSection &src)
352 : G4ScreenedCoulombCrossSection(src), phiMap(src.phiMap) { }
353
354 G4NativeScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src) : G4ScreenedCoulombCrossSection(src) { }
355 virtual ~G4NativeScreenedCoulombCrossSection();
356
357 virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff);
358 virtual G4ScreenedCoulombCrossSection *create()
359 { return new G4NativeScreenedCoulombCrossSection(*this); }
360 // get a list of available keys
361 std::vector<G4String> GetScreeningKeys() const;
362
363 typedef G4_c2_function &(*ScreeningFunc)(G4int z1, G4int z2, size_t nPoints, G4double rMax, G4double *au);
364
365 void AddScreeningFunction(G4String name, ScreeningFunc fn) {
366 phiMap[name]=fn;
367 }
368
369private:
370 // this is a map used to look up screening function generators
371 std::map<std::string, ScreeningFunc> phiMap;
372};
373
374G4_c2_function &ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
375G4_c2_function &MoliereScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
376G4_c2_function &LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
377G4_c2_function &LJZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
378
379#endif
Note: See TracBrowser for help on using the repository browser.