Geant4  10.00.p01
G4ScreenedNuclearRecoil.hh
Go to the documentation of this file.
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 //
28 //
29 // $Id: G4ScreenedNuclearRecoil.hh 66854 2013-01-14 16:56:24Z vnivanch $
30 //
31 //
32 // G4ScreenedNuclearRecoil.hh,v 1.24 2008/05/01 19:58:59 marcus Exp
33 // GEANT4 tag
34 //
35 //
36 //
37 // Class Description
38 // Process for screened electromagnetic nuclear elastic scattering;
39 // Physics comes from:
40 // Marcus H. Mendenhall and Robert A. Weller,
41 // "Algorithms for the rapid computation of classical cross sections
42 // for screened Coulomb collisions "
43 // Nuclear Instruments and Methods in Physics Research B58 (1991) 11-17
44 // The only input required is a screening function phi(r/a) which is the ratio
45 // of the actual interatomic potential for two atoms with atomic numbers Z1 and Z2,
46 // to the unscreened potential Z1*Z2*e^2/r where e^2 is elm_coupling in Geant4 units
47 // the actual screening tables are computed externally in a python module "screened_scattering.py"
48 // to allow very specific screening functions to be added if desired, without messing
49 // with the insides of this code.
50 //
51 // First version, April 2004, Marcus H. Mendenhall, Vanderbilt University
52 // May 1, 2008 -- Added code to allow process to have zero cross section above max energy, to coordinate with G4MSC. -- mhm
53 //
54 // Class Description - End
55 
56 
57 #ifndef G4ScreenedNuclearRecoil_h
58 #define G4ScreenedNuclearRecoil_h 1
59 
60 #include "globals.hh"
61 #include "G4VDiscreteProcess.hh"
62 #include "G4ParticleChange.hh"
63 #include "c2_function.hh"
64 
65 #include "CLHEP/Units/SystemOfUnits.h"
66 
67 #include <map>
68 #include <vector>
69 
71 
75 
76 typedef struct G4ScreeningTables {
80 
81 // A class for loading ScreenedCoulombCrossSections
83 {
84 public:
87 
88  static const char* CVSHeaderVers() { return
89  "G4ScreenedNuclearRecoil.hh,v 1.24 2008/05/01 19:58:59 marcus Exp GEANT4 tag ";
90  }
91  static const char* CVSFileVers();
92 };
93 
94 // A class for loading ScreenedCoulombCrossSections
96 {
97 public:
98 
103 
104  typedef std::map<G4int, G4ScreeningTables> ScreeningMap;
105 
106  // a local, fast-access mapping of a particle's Z to its full definition
107  typedef std::map<G4int, class G4ParticleDefinition *> ParticleCache;
108 
109  // LoadData is called by G4ScreenedNuclearRecoil::GetMeanFreePath
110  // It loads the data tables, builds the elemental cross-section tables.
111  virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff) = 0;
112 
113  // BuildMFPTables is called by G4ScreenedNuclearRecoil::GetMeanFreePath to build the MFP tables for each material
114  void BuildMFPTables(void); // scan the MaterialsTable and construct MFP tables
115 
116  virtual G4ScreenedCoulombCrossSection *create() = 0; // a 'virtual constructor' which clones the class
117  const G4ScreeningTables *GetScreening(G4int Z) { return &(screeningData[Z]); }
118  void SetVerbosity(G4int v) { verbosity=v; }
119 
120  // this process needs element selection weighted only by number density
122 
123  enum { nMassMapElements=116 };
124 
125  G4double standardmass(G4int z1) { return z1 <= nMassMapElements ? massmap[z1] : 2.5*z1; }
126 
127  // get the mean-free-path table for the indexed material
128  const G4_c2_function * operator [] (G4int materialIndex) {
129  return MFPTables.find(materialIndex)!=MFPTables.end() ? &(MFPTables[materialIndex].get()) : (G4_c2_function *)0;
130  }
131 
132 protected:
133  ScreeningMap screeningData; // screening tables for each element
134  ParticleCache targetMap;
136  std::map<G4int, G4_c2_const_ptr > sigmaMap; // total cross section for each element
137  std::map<G4int, G4_c2_const_ptr > MFPTables; // MFP for each material
138 
139 private:
141 
142 };
143 
144 typedef struct G4CoulombKinematicsInfo {
151 
153 public:
154  virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
155  const class G4Track& aTrack, const class G4Step& aStep)=0;
157 };
158 
160 
161 public:
163  virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
164  const class G4Track& aTrack, const class G4Step& aStep);
165 
167  const G4ScreeningTables *screen,
168  G4double eps, G4double beta);
170 protected:
171  // the c2_functions we need to do the work.
175 
176 };
177 
179 
180 public:
182  virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master,
183  const class G4Track& aTrack, const class G4Step& aStep);
184  virtual ~G4SingleScatter() {}
185 };
186 
193 {
194 public:
195 
197 
212  G4ScreenedNuclearRecoil(const G4String& processName = "ScreenedElastic",
213  const G4String &ScreeningKey="zbl", G4bool GenerateRecoils=1,
214  G4double RecoilCutoff=100.0*CLHEP::eV, G4double PhysicsCutoff=10.0*CLHEP::eV);
216  virtual ~G4ScreenedNuclearRecoil();
220  virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
223  virtual G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
226  virtual void BuildPhysicsTable(const G4ParticleDefinition& aParticleType);
229  virtual void DumpPhysicsTable(const G4ParticleDefinition& aParticleType);
234  virtual G4bool CheckNuclearCollision(G4double A, G4double A1, G4double apsis); // return true if hard collision
235 
237 
239  G4double GetNIEL() const { return NIEL; }
240 
242  void ResetTables(); // clear all data tables to allow changing energy cutoff, materials, etc.
243 
252  std::string GetScreeningKey() const { return screeningKey; }
262  void EnableRecoils(G4bool flag) { generateRecoils=flag; }
269  void SetMFPScaling(G4double scale) { MFPScale=scale; }
271  G4double GetMFPScaling() const { return MFPScale; }
292  void SetNIELPartitionFunction(const G4VNIELPartition *part);
296  void SetCrossSectionHardening(G4double fraction, G4double HardeningFactor) {
297  hardeningFraction=fraction;
298  hardeningFactor=HardeningFactor;
299  }
310  }
312  G4int GetVerboseLevel() const { return verboseLevel; }
313 
314  std::map<G4int, G4ScreenedCoulombCrossSection*> &GetCrossSectionHandlers()
315  { return crossSectionHandlers; }
316  void ClearStages(void);
317  void AddStage(G4ScreenedCollisionStage *stage) { collisionStages.push_back(stage); }
321 
323  class G4ParticleChange &GetParticleChange() { return static_cast<G4ParticleChange &>(*pParticleChange); }
325  void DepositEnergy(G4int z1, G4double a1, const G4Material *material, G4double energy);
326 
327 protected:
341 
343  std::vector<G4ScreenedCollisionStage *> collisionStages;
344 
345  std::map<G4int, G4ScreenedCoulombCrossSection*> crossSectionHandlers;
346 
350 };
351 
352 // A customized G4CrossSectionHandler which gets its data from an external program
354 {
355 public:
357 
360 
363 
364  virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff);
366  { return new G4NativeScreenedCoulombCrossSection(*this); }
367  // get a list of available keys
368  std::vector<G4String> GetScreeningKeys() const;
369 
370  typedef G4_c2_function &(*ScreeningFunc)(G4int z1, G4int z2, size_t nPoints, G4double rMax, G4double *au);
371 
373  phiMap[name]=fn;
374  }
375 
376 private:
377  // this is a map used to look up screening function generators
378  std::map<std::string, ScreeningFunc> phiMap;
379 };
380 
381 G4_c2_function &ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
382 G4_c2_function &MoliereScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
383 G4_c2_function &LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
384 G4_c2_function &LJZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval);
385 
386 #endif
G4NativeScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src)
std::vector< G4String > GetScreeningKeys() const
const G4_c2_function * operator[](G4int materialIndex)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4CoulombKinematicsInfo kinematics
G4double GetCurrentInteractionLength() const
the the interaciton length used in the last scattering.
void AddStage(G4ScreenedCollisionStage *stage)
virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff)=0
struct G4ScreeningTables G4ScreeningTables
G4double highEnergyLimit
the energy per nucleon above which the MFP is constant
virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master, const class G4Track &aTrack, const class G4Step &aStep)
G4ParticleDefinition * recoilIon
G4bool GetAvoidNuclearReactions() const
get the flag indicating whether hadronic collisions are ignored.
G4CoulombKinematicsInfo & GetKinematics()
G4ScreenedCoulombCrossSection(const G4ScreenedCoulombCrossSection &src)
G4String name
Definition: TRTMaterials.hh:40
static const G4double a1
virtual void DumpPhysicsTable(const G4ParticleDefinition &aParticleType)
Export physics tables for persistency. Not Implemented.
G4double GetRecoilCutoff() const
get the recoil cutoff
G4ScreenedNuclearRecoil(const G4String &processName="ScreenedElastic", const G4String &ScreeningKey="zbl", G4bool GenerateRecoils=1, G4double RecoilCutoff=100.0 *CLHEP::eV, G4double PhysicsCutoff=10.0 *CLHEP::eV)
Construct the process and set some physics parameters for it.
const G4VNIELPartition * NIELPartitionFunction
static const G4double massmap[nMassMapElements+1]
static const G4double eps
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
used internally by Geant4 machinery
void DepositEnergy(G4int z1, G4double a1, const G4Material *material, G4double energy)
take the given energy, and use the material information to partition it into NIEL and ionizing energy...
struct G4CoulombKinematicsInfo G4CoulombKinematicsInfo
int G4int
Definition: G4Types.hh:78
virtual void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
Build physics tables in advance. Not Implemented.
std::map< G4int, G4ScreeningTables > ScreeningMap
Provides the headers for the general c2_function algebra which supports fast, flexible operations on ...
std::vector< G4ScreenedCollisionStage * > collisionStages
virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master, const class G4Track &aTrack, const class G4Step &aStep)=0
G4_c2_function & LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
std::map< std::string, ScreeningFunc > phiMap
void SetPhysicsCutoff(G4double energy)
set the energy to which screening tables are computed. Typically, this is 10 eV or so...
void AllowEnergyDeposition(G4bool flag)
enable or disable all energy deposition by this process
c2_const_plugin_function_p< G4double > & phifunc
G4ScreenedCoulombCrossSection * externalCrossSectionConstructor
virtual G4ScreenedCoulombCrossSection * create()
G4ScreenedCoulombCrossSection * crossSection
G4double GetHardeningFactor() const
get the boost factor in use.
void SetRecoilCutoff(G4double energy)
set the minimum energy (per nucleon) at which recoils can be generated, and the energy (per nucleon) ...
std::map< G4int, G4_c2_const_ptr > MFPTables
std::map< G4int, class G4ParticleDefinition * > ParticleCache
void SetMFPScaling(G4double scale)
set the mean free path scaling as specified
bool G4bool
Definition: G4Types.hh:79
std::map< G4int, G4ScreenedCoulombCrossSection * > crossSectionHandlers
G4double GetMFPScaling() const
get the MFPScaling parameter
G4double currentInteractionLength
Definition: G4VProcess.hh:297
const G4ScreeningTables * GetScreening(G4int Z)
G4ParticleDefinition * SelectRandomUnweightedTarget(const G4MaterialCutsCouple *couple)
G4bool GetAllowEnergyDeposition() const
get flag indicating whether deposition is enabled
G4double lowEnergyLimit
the energy per nucleon below which the MFP is zero
std::map< G4int, G4_c2_const_ptr > sigmaMap
G4_c2_function & ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
void SetCrossSectionHardening(G4double fraction, G4double HardeningFactor)
set the cross section boost to provide faster computation of backscattering
Definition: G4Step.hh:76
G4bool GetEnableRecoils() const
find out if generation of recoils is enabled.
static const G4double A[nN]
void SetMaxEnergyForScattering(G4double energy)
set the upper energy beyond which this process has no cross section
void AvoidNuclearReactions(G4bool flag)
enable or disable whether this process will skip collisions which are close enough they need hadronic...
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
test if a prticle of type aParticleType can use this process
void SetExternalCrossSectionHandler(G4ScreenedCoulombCrossSection *cs)
set a function to compute screening tables, if the user needs non-standard behavior.
void SetNIELPartitionFunction(const G4VNIELPartition *part)
set the pointer to a class for paritioning energy into NIEL
virtual G4ScreenedCoulombCrossSection * create()=0
static const double eV
Definition: G4SIunits.hh:194
c2_function< G4double > G4_c2_function
G4double energy(const ThreeVector &p, const G4double m)
G4_c2_function & LJZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
std::string GetScreeningKey() const
find out what screening funciton we are using
c2_const_ptr< G4double > G4_c2_const_ptr
virtual ~G4ScreenedNuclearRecoil()
destructor
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
used internally by Geant4 machinery
virtual G4bool CheckNuclearCollision(G4double A, G4double A1, G4double apsis)
deterine if the moving particle is within the strong force range of the selected nucleus ...
void AddScreeningFunction(G4String name, ScreeningFunc fn)
c2_ptr< G4double > G4_c2_ptr
virtual G4ScreenedCoulombCrossSection * GetNewCrossSectionHandler(void)
G4double GetHardeningFraction() const
get the fraction of particles which will have boosted scattering
G4double GetPhysicsCutoff() const
get the physics cutoff energy.
virtual void DoCollisionStep(class G4ScreenedNuclearRecoil *master, const class G4Track &aTrack, const class G4Step &aStep)
double G4double
Definition: G4Types.hh:76
G4int GetVerboseLevel() const
get the verbosity.
G4_c2_function &(* ScreeningFunc)(G4int z1, G4int z2, size_t nPoints, G4double rMax, G4double *au)
G4double GetNIEL() const
Get non-ionizing energy loss for last step.
G4ForceCondition
class G4ParticleChange & GetParticleChange()
get the pointer to our ParticleChange object. for internal use, primarily.
void EnableRecoils(G4bool flag)
enable or disable the generation of recoils. If recoils are disabled, the energy they would have rece...
G4bool DoScreeningComputation(class G4ScreenedNuclearRecoil *master, const G4ScreeningTables *screen, G4double eps, G4double beta)
G4_c2_function & MoliereScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
void ResetTables()
clear precomputed screening tables
A process which handles screened Coulomb collisions between nuclei.
virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff)
G4NativeScreenedCoulombCrossSection(const G4NativeScreenedCoulombCrossSection &src)
G4double processMaxEnergy
the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC ...
std::map< G4int, G4ScreenedCoulombCrossSection * > & GetCrossSectionHandlers()