00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 #ifndef G4DNARuddIonisationModel_h
00030 #define G4DNARuddIonisationModel_h 1
00031
00032 #include "G4VEmModel.hh"
00033 #include "G4ParticleChangeForGamma.hh"
00034 #include "G4ProductionCutsTable.hh"
00035
00036 #include "G4DNAGenericIonsManager.hh"
00037 #include "G4DNACrossSectionDataSet.hh"
00038 #include "G4Electron.hh"
00039 #include "G4Proton.hh"
00040 #include "G4LogLogInterpolation.hh"
00041
00042 #include "G4DNAWaterIonisationStructure.hh"
00043 #include "G4VAtomDeexcitation.hh"
00044 #include "G4NistManager.hh"
00045
00046 class G4DNARuddIonisationModel : public G4VEmModel
00047 {
00048
00049 public:
00050
00051 G4DNARuddIonisationModel(const G4ParticleDefinition* p = 0,
00052 const G4String& nam = "DNARuddIonisationModel");
00053
00054 virtual ~G4DNARuddIonisationModel();
00055
00056 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00057
00058 virtual G4double CrossSectionPerVolume( const G4Material* material,
00059 const G4ParticleDefinition* p,
00060 G4double ekin,
00061 G4double emin,
00062 G4double emax);
00063
00064 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00065 const G4MaterialCutsCouple*,
00066 const G4DynamicParticle*,
00067 G4double tmin,
00068 G4double maxEnergy);
00069
00070 protected:
00071
00072 G4ParticleChangeForGamma* fParticleChangeForGamma;
00073
00074 private:
00075
00076 const std::vector<G4double>* fpWaterDensity;
00077
00078
00079 G4VAtomDeexcitation* fAtomDeexcitation;
00080
00081 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
00082 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
00083
00084 G4double lowEnergyLimitForZ1;
00085 G4double lowEnergyLimitForZ2;
00086 G4double lowEnergyLimitOfModelForZ1;
00087 G4double lowEnergyLimitOfModelForZ2;
00088 G4double killBelowEnergyForZ1;
00089 G4double killBelowEnergyForZ2;
00090
00091 G4bool isInitialised;
00092 G4int verboseLevel;
00093
00094
00095
00096 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
00097 MapFile tableFile;
00098
00099 typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
00100 MapData tableData;
00101
00102
00103
00104 G4DNAWaterIonisationStructure waterStructure;
00105
00106 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition* particleDefinition,
00107 G4double incomingParticleEnergy,
00108 G4int shell);
00109
00110 void RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
00111 G4double incomingParticleEnergy,
00112 G4double outgoingParticleEnergy,
00113 G4double & cosTheta,
00114 G4double & phi);
00115
00116 G4double DifferentialCrossSection(G4ParticleDefinition* particleDefinition,
00117 G4double k,
00118 G4double energyTransfer,
00119 G4int shell);
00120
00121 G4double CorrectionFactor(G4ParticleDefinition* particleDefinition, G4double k);
00122
00123 G4double S_1s(G4double t,
00124 G4double energyTransferred,
00125 G4double slaterEffectiveChg,
00126 G4double shellNumber);
00127
00128 G4double S_2s(G4double t,
00129 G4double energyTransferred,
00130 G4double slaterEffectiveChg,
00131 G4double shellNumber);
00132
00133
00134 G4double S_2p(G4double t,
00135 G4double energyTransferred,
00136 G4double slaterEffectiveChg,
00137 G4double shellNumber);
00138
00139 G4double R(G4double t,
00140 G4double energyTransferred,
00141 G4double slaterEffectiveChg,
00142 G4double shellNumber) ;
00143
00144 G4double slaterEffectiveCharge[3];
00145 G4double sCoefficient[3];
00146
00147
00148
00149 G4double PartialCrossSection(const G4Track& track);
00150
00151 G4double Sum(G4double energy, const G4String& particle);
00152
00153 G4int RandomSelect(G4double energy,const G4String& particle );
00154
00155
00156
00157 G4DNARuddIonisationModel & operator=(const G4DNARuddIonisationModel &right);
00158 G4DNARuddIonisationModel(const G4DNARuddIonisationModel&);
00159
00160 };
00161
00162
00163
00164 #endif