source: trunk/source/processes/electromagnetic/lowenergy/include/G4DNARuddIonisationModel.hh@ 992

Last change on this file since 992 was 966, checked in by garnier, 17 years ago

fichier ajoutes

File size: 5.2 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: G4DNARuddIonisationModel.hh,v 1.1 2009/01/12 14:26:03 sincerti Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29
30#ifndef G4DNARuddIonisationModel_h
31#define G4DNARuddIonisationModel_h 1
32
33#include "G4VEmModel.hh"
34#include "G4ParticleChangeForGamma.hh"
35#include "G4ProductionCutsTable.hh"
36
37#include "G4DNAGenericIonsManager.hh"
38#include "G4DNACrossSectionDataSet.hh"
39#include "G4Electron.hh"
40#include "G4Proton.hh"
41#include "G4LogLogInterpolation.hh"
42
43#include "G4WaterIonisationStructure.hh"
44
45class G4DNARuddIonisationModel : public G4VEmModel
46{
47
48public:
49
50 G4DNARuddIonisationModel(const G4ParticleDefinition* p = 0,
51 const G4String& nam = "DNARuddIonisationModel");
52
53 virtual ~G4DNARuddIonisationModel();
54
55 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
56
57 virtual G4double CrossSectionPerVolume( const G4Material* material,
58 const G4ParticleDefinition* p,
59 G4double ekin,
60 G4double emin,
61 G4double emax);
62
63 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
64 const G4MaterialCutsCouple*,
65 const G4DynamicParticle*,
66 G4double tmin,
67 G4double maxEnergy);
68
69protected:
70
71 G4ParticleChangeForGamma* fParticleChangeForGamma;
72
73private:
74
75 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
76 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
77
78 G4double lowEnergyLimitForZ1;
79 G4double lowEnergyLimitForZ2;
80 G4double lowEnergyLimitOfModelForZ1;
81 G4double lowEnergyLimitOfModelForZ2;
82 G4double killBelowEnergyForZ1;
83 G4double killBelowEnergyForZ2;
84
85 G4bool isInitialised;
86 G4int verboseLevel;
87
88 // Cross section
89
90 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
91 MapFile tableFile;
92
93 typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
94 MapData tableData;
95
96 // Final state
97
98 G4WaterIonisationStructure waterStructure;
99
100 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
101 G4double incomingParticleEnergy,
102 G4int shell);
103
104 void RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
105 G4double incomingParticleEnergy,
106 G4double outgoingParticleEnergy,
107 G4double & cosTheta,
108 G4double & phi);
109
110 G4double DifferentialCrossSection(G4ParticleDefinition* particleDefinition,
111 G4double k,
112 G4double energyTransfer,
113 G4int shell);
114
115 G4double CorrectionFactor(G4ParticleDefinition* particleDefinition, G4double k);
116
117 G4double S_1s(G4double t,
118 G4double energyTransferred,
119 G4double slaterEffectiveChg,
120 G4double shellNumber);
121
122 G4double S_2s(G4double t,
123 G4double energyTransferred,
124 G4double slaterEffectiveChg,
125 G4double shellNumber);
126
127
128 G4double S_2p(G4double t,
129 G4double energyTransferred,
130 G4double slaterEffectiveChg,
131 G4double shellNumber);
132
133 G4double R(G4double t,
134 G4double energyTransferred,
135 G4double slaterEffectiveChg,
136 G4double shellNumber) ;
137
138 G4double slaterEffectiveCharge[3];
139 G4double sCoefficient[3];
140
141 // Partial cross section
142
143 G4double PartialCrossSection(const G4Track& track);
144
145 G4double Sum(G4double energy, const G4String& particle);
146
147 G4int RandomSelect(G4double energy,const G4String& particle );
148
149 // Test water material
150
151 G4bool flagMaterialIsWater;
152 G4double densityWater;
153
154 //
155
156 G4DNARuddIonisationModel & operator=(const G4DNARuddIonisationModel &right);
157 G4DNARuddIonisationModel(const G4DNARuddIonisationModel&);
158
159};
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162
163#endif
Note: See TracBrowser for help on using the repository browser.