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 G4ElectronIonPair_h
00030 #define G4ElectronIonPair_h 1
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062 #include "globals.hh"
00063 #include "G4Step.hh"
00064 #include "G4ParticleDefinition.hh"
00065 #include "G4ThreeVector.hh"
00066 #include "G4VProcess.hh"
00067 #include <vector>
00068
00069
00070
00071 class G4Material;
00072
00073 class G4ElectronIonPair
00074 {
00075 public:
00076
00077 G4ElectronIonPair();
00078
00079 virtual ~G4ElectronIonPair();
00080
00081
00082 G4double MeanNumberOfIonsAlongStep(const G4ParticleDefinition*,
00083 const G4Material*,
00084 G4double edepTotal,
00085 G4double edepNIEL = 0.0);
00086
00087 inline G4double MeanNumberOfIonsAlongStep(const G4Step*);
00088
00089 inline G4int SampleNumberOfIonsAlongStep(const G4Step*);
00090
00091
00092
00093 std::vector<G4ThreeVector>* SampleIonsAlongStep(const G4Step*);
00094
00095
00096 G4int ResidualeChargePostStep(const G4ParticleDefinition*,
00097 const G4TrackVector* secondary = 0,
00098 G4int processSubType = -1);
00099
00100 inline G4int ResidualeChargePostStep(const G4Step*);
00101
00102
00103 G4double FindG4MeanEnergyPerIonPair(const G4Material*);
00104
00105
00106 void DumpMeanEnergyPerIonPair();
00107
00108
00109 void DumpG4MeanEnergyPerIonPair();
00110
00111 inline void SetVerbose(G4int);
00112
00113 private:
00114
00115 void Initialise();
00116
00117 G4double FindMeanEnergyPerIonPair(const G4Material*);
00118
00119
00120 G4ElectronIonPair & operator=(const G4ElectronIonPair &right);
00121 G4ElectronIonPair(const G4ElectronIonPair&);
00122
00123
00124 const G4Material* curMaterial;
00125 G4double curMeanEnergy;
00126
00127 G4double invFanoFactor;
00128
00129 G4int verbose;
00130 G4int nMaterials;
00131
00132
00133 std::vector<G4double> g4MatData;
00134 std::vector<G4String> g4MatNames;
00135 };
00136
00137 inline G4double
00138 G4ElectronIonPair::MeanNumberOfIonsAlongStep(const G4Step* step)
00139 {
00140 return MeanNumberOfIonsAlongStep(step->GetTrack()->GetParticleDefinition(),
00141 step->GetPreStepPoint()->GetMaterial(),
00142 step->GetTotalEnergyDeposit(),
00143 step->GetNonIonizingEnergyDeposit());
00144 }
00145
00146 inline
00147 G4int G4ElectronIonPair::SampleNumberOfIonsAlongStep(const G4Step* step)
00148 {
00149
00150
00151 G4double meanion = MeanNumberOfIonsAlongStep(step);
00152 return G4lrint(CLHEP::RandGamma::shoot(meanion*invFanoFactor,invFanoFactor));
00153 }
00154
00155 inline
00156 G4int G4ElectronIonPair::ResidualeChargePostStep(const G4Step* step)
00157 {
00158 G4int subtype = -1;
00159 const G4VProcess* proc = step->GetPostStepPoint()->GetProcessDefinedStep();
00160 if(proc) { subtype = proc->GetProcessSubType(); }
00161 return ResidualeChargePostStep(step->GetTrack()->GetParticleDefinition(),
00162 step->GetSecondary(),
00163 subtype);
00164 }
00165
00166 inline void G4ElectronIonPair::SetVerbose(G4int val)
00167 {
00168 verbose = val;
00169 }
00170
00171 #endif
00172