Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ElectronIonPair Class Reference

#include <G4ElectronIonPair.hh>

Public Member Functions

 G4ElectronIonPair (G4int verb)
 
virtual ~G4ElectronIonPair ()
 
G4double MeanNumberOfIonsAlongStep (const G4ParticleDefinition *, const G4Material *, G4double edepTotal, G4double edepNIEL=0.0)
 
G4double MeanNumberOfIonsAlongStep (const G4Step *)
 
G4int SampleNumberOfIonsAlongStep (const G4Step *)
 
std::vector< G4ThreeVector > * SampleIonsAlongStep (const G4Step *)
 
G4int ResidualeChargePostStep (const G4ParticleDefinition *, const G4TrackVector *secondary=nullptr, G4int processSubType=-1) const
 
G4int ResidualeChargePostStep (const G4Step *) const
 
G4double FindG4MeanEnergyPerIonPair (const G4Material *) const
 
void DumpMeanEnergyPerIonPair () const
 
void DumpG4MeanEnergyPerIonPair () const
 
void SetVerbose (G4int)
 

Detailed Description

Definition at line 73 of file G4ElectronIonPair.hh.

Constructor & Destructor Documentation

◆ G4ElectronIonPair()

G4ElectronIonPair::G4ElectronIonPair ( G4int  verb)
explicit

Definition at line 57 of file G4ElectronIonPair.cc.

58{
59 verbose = verb;
60 curMaterial = nullptr;
61 curMeanEnergy = 0.0;
62 nMaterials = 0;
63 invFanoFactor = 1.0/0.2;
64 Initialise();
65}

◆ ~G4ElectronIonPair()

G4ElectronIonPair::~G4ElectronIonPair ( )
virtual

Definition at line 69 of file G4ElectronIonPair.cc.

70{}

Member Function Documentation

◆ DumpG4MeanEnergyPerIonPair()

void G4ElectronIonPair::DumpG4MeanEnergyPerIonPair ( ) const

Definition at line 189 of file G4ElectronIonPair.cc.

190{
191 if(nMaterials > 0) {
192 G4cout << "### G4ElectronIonPair: mean energy per ion pair "
193 << " for Geant4 materials" << G4endl;
194 for(G4int i=0; i<nMaterials; ++i) {
195 G4cout << " " << g4MatNames[i] << " Epair= "
196 << g4MatData[i]/eV << " eV" << G4endl;
197 }
198 }
199}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ DumpMeanEnergyPerIonPair()

void G4ElectronIonPair::DumpMeanEnergyPerIonPair ( ) const

Definition at line 169 of file G4ElectronIonPair.cc.

170{
173 if(nmat > 0) {
174 G4cout << "### G4ElectronIonPair: mean energy per ion pair available:"
175 << G4endl;
176 for(G4int i=0; i<nmat; ++i) {
177 const G4Material* mat = (*mtable)[i];
179 if(x > 0.0) {
180 G4cout << " " << mat->GetName() << " Epair= "
181 << x/eV << " eV" << G4endl;
182 }
183 }
184 }
185}
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
G4double GetMeanEnergyPerIonPair() const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
const G4String & GetName() const
Definition: G4Material.hh:175

◆ FindG4MeanEnergyPerIonPair()

G4double G4ElectronIonPair::FindG4MeanEnergyPerIonPair ( const G4Material mat) const

Definition at line 146 of file G4ElectronIonPair.cc.

147{
148 G4String name = mat->GetName();
149 G4double res = 0.0;
150
151 // is this material in the vector?
152 for(G4int j=0; j<nMaterials; j++) {
153 if(name == g4MatNames[j]) {
154 res = g4MatData[j];
156 if(verbose > 0) {
157 G4cout << "### G4ElectronIonPair::FindG4MeanEnergyPerIonPair for "
158 << name << " Epair= " << res/eV << " eV is set"
159 << G4endl;
160 }
161 break;
162 }
163 }
164 return res;
165}
void SetMeanEnergyPerIonPair(G4double value)
const char * name(G4int ptype)

Referenced by MeanNumberOfIonsAlongStep().

◆ MeanNumberOfIonsAlongStep() [1/2]

G4double G4ElectronIonPair::MeanNumberOfIonsAlongStep ( const G4ParticleDefinition part,
const G4Material material,
G4double  edepTotal,
G4double  edepNIEL = 0.0 
)

Definition at line 74 of file G4ElectronIonPair.cc.

79{
80 G4double nion = 0.0;
81
82 // NIEL does not provide ionisation clusters
83 if(edep > niel) {
84
85 // neutral particles do not produce ionisation along step
86 if(part->GetPDGCharge() != 0.0) {
87
88 // select material
89 if(material != curMaterial) {
90 curMaterial = material;
91 curMeanEnergy = material->GetIonisation()->GetMeanEnergyPerIonPair();
92
93 // if mean energy is not defined then look into G4 DB
94 if(0.0 == curMeanEnergy) {
95 curMeanEnergy = FindG4MeanEnergyPerIonPair(material);
96 }
97 }
98 if(curMeanEnergy > 0.0) { nion = (edep - niel)/curMeanEnergy; }
99 }
100 }
101 return nion;
102}
G4double FindG4MeanEnergyPerIonPair(const G4Material *) const
G4double GetPDGCharge() const

Referenced by MeanNumberOfIonsAlongStep(), and SampleNumberOfIonsAlongStep().

◆ MeanNumberOfIonsAlongStep() [2/2]

G4double G4ElectronIonPair::MeanNumberOfIonsAlongStep ( const G4Step step)
inline

Definition at line 138 of file G4ElectronIonPair.hh.

139{
141 step->GetPreStepPoint()->GetMaterial(),
142 step->GetTotalEnergyDeposit(),
144}
G4double MeanNumberOfIonsAlongStep(const G4ParticleDefinition *, const G4Material *, G4double edepTotal, G4double edepNIEL=0.0)
G4Material * GetMaterial() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetNonIonizingEnergyDeposit() const
G4double GetTotalEnergyDeposit() const
const G4ParticleDefinition * GetParticleDefinition() const

◆ ResidualeChargePostStep() [1/2]

G4int G4ElectronIonPair::ResidualeChargePostStep ( const G4ParticleDefinition ,
const G4TrackVector secondary = nullptr,
G4int  processSubType = -1 
) const

Definition at line 133 of file G4ElectronIonPair.cc.

136{
137 G4int nholes = 0;
138
139 if(2 == subType || 12 == subType || 13 == subType) { nholes = 1; }
140 return nholes;
141}

Referenced by ResidualeChargePostStep().

◆ ResidualeChargePostStep() [2/2]

G4int G4ElectronIonPair::ResidualeChargePostStep ( const G4Step step) const
inline

Definition at line 156 of file G4ElectronIonPair.hh.

157{
158 G4int subtype = -1;
159 const G4VProcess* proc = step->GetPostStepPoint()->GetProcessDefinedStep();
160 if(proc) { subtype = proc->GetProcessSubType(); }
162 step->GetSecondary(),
163 subtype);
164}
G4int ResidualeChargePostStep(const G4ParticleDefinition *, const G4TrackVector *secondary=nullptr, G4int processSubType=-1) const
const G4VProcess * GetProcessDefinedStep() const
const G4TrackVector * GetSecondary() const
G4StepPoint * GetPostStepPoint() const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400

◆ SampleIonsAlongStep()

std::vector< G4ThreeVector > * G4ElectronIonPair::SampleIonsAlongStep ( const G4Step step)

Definition at line 107 of file G4ElectronIonPair.cc.

108{
109 std::vector<G4ThreeVector>* v = 0;
110
112
113 // sample ionisation along step
114 if(nion > 0) {
115
116 v = new std::vector<G4ThreeVector>;
117 G4ThreeVector prePos = step->GetPreStepPoint()->GetPosition();
118 G4ThreeVector deltaPos = step->GetPostStepPoint()->GetPosition() - prePos;
119 for(G4int i=0; i<nion; ++i) {
120 v->push_back( prePos + deltaPos*G4UniformRand() );
121 }
122 if(verbose > 1 ) {
123 G4cout << "### G4ElectronIonPair::SampleIonisationPoints: "
124 << v->size() << " ion pairs are added" << G4endl;
125 }
126 }
127 return v;
128}
#define G4UniformRand()
Definition: Randomize.hh:52
G4int SampleNumberOfIonsAlongStep(const G4Step *)
const G4ThreeVector & GetPosition() const

◆ SampleNumberOfIonsAlongStep()

G4int G4ElectronIonPair::SampleNumberOfIonsAlongStep ( const G4Step step)
inline

Definition at line 147 of file G4ElectronIonPair.hh.

148{
149 // use gamma distribution with mean value n=meanion and
150 // dispersion D=meanion/invFanoFactor
151 G4double meanion = MeanNumberOfIonsAlongStep(step);
152 return G4lrint(G4RandGamma::shoot(meanion*invFanoFactor,invFanoFactor));
153}
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by SampleIonsAlongStep().

◆ SetVerbose()

void G4ElectronIonPair::SetVerbose ( G4int  val)
inline

Definition at line 166 of file G4ElectronIonPair.hh.

167{
168 verbose = val;
169}

Referenced by G4LossTableManager::ResetParameters().


The documentation for this class was generated from the following files: