File indexing completed on 2025-01-18 09:58:11
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026 #ifndef G4ElectronIonPair_h
0027 #define G4ElectronIonPair_h 1
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059 #include "globals.hh"
0060 #include "G4Step.hh"
0061 #include "G4ParticleDefinition.hh"
0062 #include "G4ThreeVector.hh"
0063 #include "G4VProcess.hh"
0064 #include "Randomize.hh"
0065 #include <vector>
0066
0067
0068
0069 class G4Material;
0070
0071 class G4ElectronIonPair
0072 {
0073 public:
0074
0075 explicit G4ElectronIonPair(G4int verb);
0076
0077 virtual ~G4ElectronIonPair();
0078
0079
0080 G4double MeanNumberOfIonsAlongStep(const G4ParticleDefinition*,
0081 const G4Material*,
0082 G4double edepTotal,
0083 G4double edepNIEL = 0.0);
0084
0085 inline G4double MeanNumberOfIonsAlongStep(const G4Step*);
0086
0087 inline G4int SampleNumberOfIonsAlongStep(const G4Step*);
0088
0089
0090
0091 std::vector<G4ThreeVector>* SampleIonsAlongStep(const G4Step*);
0092
0093
0094 G4int ResidualeChargePostStep(const G4ParticleDefinition*,
0095 const G4TrackVector* secondary = nullptr,
0096 G4int processSubType = -1) const;
0097
0098 inline G4int ResidualeChargePostStep(const G4Step*) const;
0099
0100
0101 G4double FindG4MeanEnergyPerIonPair(const G4Material*) const;
0102
0103
0104 void DumpMeanEnergyPerIonPair() const;
0105
0106
0107 void DumpG4MeanEnergyPerIonPair() const;
0108
0109 inline void SetVerbose(G4int);
0110
0111
0112 G4ElectronIonPair & operator=(const G4ElectronIonPair &right) = delete;
0113 G4ElectronIonPair(const G4ElectronIonPair&) = delete;
0114
0115 private:
0116
0117 void Initialise();
0118
0119 G4double FindMeanEnergyPerIonPair(const G4Material*) const;
0120
0121
0122 const G4Material* curMaterial;
0123 G4double curMeanEnergy;
0124
0125 G4double invFanoFactor;
0126
0127 G4int verbose;
0128 G4int nMaterials;
0129
0130
0131 std::vector<G4double> g4MatData;
0132 std::vector<G4String> g4MatNames;
0133 };
0134
0135 inline G4double
0136 G4ElectronIonPair::MeanNumberOfIonsAlongStep(const G4Step* step)
0137 {
0138 return MeanNumberOfIonsAlongStep(step->GetTrack()->GetParticleDefinition(),
0139 step->GetPreStepPoint()->GetMaterial(),
0140 step->GetTotalEnergyDeposit(),
0141 step->GetNonIonizingEnergyDeposit());
0142 }
0143
0144 inline G4int
0145 G4ElectronIonPair::SampleNumberOfIonsAlongStep(const G4Step* step)
0146 {
0147
0148
0149 G4double meanion = MeanNumberOfIonsAlongStep(step);
0150 return G4lrint(G4RandGamma::shoot(meanion*invFanoFactor,invFanoFactor));
0151 }
0152
0153 inline G4int
0154 G4ElectronIonPair::ResidualeChargePostStep(const G4Step* step) const
0155 {
0156 G4int subtype = -1;
0157 const G4VProcess* proc = step->GetPostStepPoint()->GetProcessDefinedStep();
0158 if(proc) { subtype = proc->GetProcessSubType(); }
0159 return ResidualeChargePostStep(step->GetTrack()->GetParticleDefinition(),
0160 step->GetSecondary(),
0161 subtype);
0162 }
0163
0164 inline void G4ElectronIonPair::SetVerbose(G4int val)
0165 {
0166 verbose = val;
0167 }
0168
0169 #endif
0170