Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNARuddIonisationModel.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 //
26 // $Id: G4DNARuddIonisationModel.hh 96060 2016-03-11 12:58:04Z gcosmo $
27 //
28 
29 #ifndef G4DNARuddIonisationModel_h
30 #define G4DNARuddIonisationModel_h 1
31 
32 #include "G4VEmModel.hh"
34 #include "G4ProductionCutsTable.hh"
35 
38 #include "G4Electron.hh"
39 #include "G4Proton.hh"
40 #include "G4LogLogInterpolation.hh"
41 
43 #include "G4VAtomDeexcitation.hh"
44 #include "G4NistManager.hh"
45 
47 {
48 
49 public:
50 
52  const G4String& nam = "DNARuddIonisationModel");
53 
54  virtual ~G4DNARuddIonisationModel();
55 
56  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
57 
58  virtual G4double CrossSectionPerVolume( const G4Material* material,
59  const G4ParticleDefinition* p,
60  G4double ekin,
61  G4double emin,
62  G4double emax);
63 
64  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
65  const G4MaterialCutsCouple*,
66  const G4DynamicParticle*,
67  G4double tmin,
68  G4double maxEnergy);
69 
70  inline void SelectStationary(G4bool input);
71 
72 protected:
73 
75 
76 private:
77 
78  G4bool statCode;
79 
80  // Water density table
81  const std::vector<G4double>* fpWaterDensity;
82 
83  //deexcitation manager to produce fluo photns and e-
84  G4VAtomDeexcitation* fAtomDeexcitation;
85 
86  std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
87  std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
88 
89  G4double lowEnergyLimitForZ1;
90  G4double lowEnergyLimitForZ2;
91  G4double lowEnergyLimitOfModelForZ1;
92  G4double lowEnergyLimitOfModelForZ2;
93  G4double killBelowEnergyForZ1;
94  G4double killBelowEnergyForZ2;
95 
96  G4bool isInitialised;
97  G4int verboseLevel;
98 
99  // Cross section
100 
101  typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
102  MapFile tableFile;
103 
104  typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
105  MapData tableData;
106 
107  // Final state
108 
109  G4DNAWaterIonisationStructure waterStructure;
110 
111  G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
112  G4double incomingParticleEnergy,
113  G4int shell);
114 
115  G4double DifferentialCrossSection(G4ParticleDefinition* particleDefinition,
116  G4double k,
117  G4double energyTransfer,
118  G4int shell);
119 
120  G4double CorrectionFactor(G4ParticleDefinition* particleDefinition, G4double k);
121 
122  G4double S_1s(G4double t,
123  G4double energyTransferred,
124  G4double slaterEffectiveChg,
125  G4double shellNumber);
126 
127  G4double S_2s(G4double t,
128  G4double energyTransferred,
129  G4double slaterEffectiveChg,
130  G4double shellNumber);
131 
132 
133  G4double S_2p(G4double t,
134  G4double energyTransferred,
135  G4double slaterEffectiveChg,
136  G4double shellNumber);
137 
138  G4double R(G4double t,
139  G4double energyTransferred,
140  G4double slaterEffectiveChg,
141  G4double shellNumber) ;
142 
143  G4double slaterEffectiveCharge[3];
144  G4double sCoefficient[3];
145 
146  // Partial cross section
147 
148  G4double PartialCrossSection(const G4Track& track);
149 
150  G4double Sum(G4double energy, const G4String& particle);
151 
152  G4int RandomSelect(G4double energy,const G4String& particle );
153 
154  //
155 
156  G4DNARuddIonisationModel & operator=(const G4DNARuddIonisationModel &right);
158 
159 };
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
164 {
165  statCode = input;
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
170 #endif
G4ParticleChangeForGamma * fParticleChangeForGamma
const char * p
Definition: xmltok.h:285
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
int G4int
Definition: G4Types.hh:78
G4DNARuddIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNARuddIonisationModel")
bool G4bool
Definition: G4Types.hh:79
static const G4double emax
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)