Geant4 11.2.2
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)
 
G4ElectronIonPairoperator= (const G4ElectronIonPair &right)=delete
 
 G4ElectronIonPair (const G4ElectronIonPair &)=delete
 

Detailed Description

Definition at line 71 of file G4ElectronIonPair.hh.

Constructor & Destructor Documentation

◆ G4ElectronIonPair() [1/2]

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 ( )
virtualdefault

◆ G4ElectronIonPair() [2/2]

G4ElectronIonPair::G4ElectronIonPair ( const G4ElectronIonPair & )
delete

Member Function Documentation

◆ DumpG4MeanEnergyPerIonPair()

void G4ElectronIonPair::DumpG4MeanEnergyPerIonPair ( ) const

Definition at line 188 of file G4ElectronIonPair.cc.

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

◆ DumpMeanEnergyPerIonPair()

void G4ElectronIonPair::DumpMeanEnergyPerIonPair ( ) const

Definition at line 168 of file G4ElectronIonPair.cc.

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

◆ FindG4MeanEnergyPerIonPair()

G4double G4ElectronIonPair::FindG4MeanEnergyPerIonPair ( const G4Material * mat) const

Definition at line 145 of file G4ElectronIonPair.cc.

146{
147 G4String name = mat->GetName();
148 G4double res = 0.0;
149
150 // is this material in the vector?
151 for(G4int j=0; j<nMaterials; ++j) {
152 if(name == g4MatNames[j]) {
153 res = g4MatData[j];
155 if(verbose > 0) {
156 G4cout << "### G4ElectronIonPair::FindG4MeanEnergyPerIonPair for "
157 << name << " Epair= " << res/eV << " eV is set"
158 << G4endl;
159 }
160 break;
161 }
162 }
163 return res;
164}
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 73 of file G4ElectronIonPair.cc.

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

Referenced by MeanNumberOfIonsAlongStep(), and SampleNumberOfIonsAlongStep().

◆ MeanNumberOfIonsAlongStep() [2/2]

G4double G4ElectronIonPair::MeanNumberOfIonsAlongStep ( const G4Step * step)
inline

Definition at line 136 of file G4ElectronIonPair.hh.

137{
139 step->GetPreStepPoint()->GetMaterial(),
140 step->GetTotalEnergyDeposit(),
142}
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

◆ operator=()

G4ElectronIonPair & G4ElectronIonPair::operator= ( const G4ElectronIonPair & right)
delete

◆ ResidualeChargePostStep() [1/2]

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

Definition at line 132 of file G4ElectronIonPair.cc.

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

Referenced by ResidualeChargePostStep().

◆ ResidualeChargePostStep() [2/2]

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

Definition at line 154 of file G4ElectronIonPair.hh.

155{
156 G4int subtype = -1;
157 const G4VProcess* proc = step->GetPostStepPoint()->GetProcessDefinedStep();
158 if(proc) { subtype = proc->GetProcessSubType(); }
160 step->GetSecondary(),
161 subtype);
162}
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

◆ SampleIonsAlongStep()

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

Definition at line 106 of file G4ElectronIonPair.cc.

107{
108 std::vector<G4ThreeVector>* v = nullptr;
109
111
112 // sample ionisation along step
113 if(nion > 0) {
114
115 v = new std::vector<G4ThreeVector>;
116 G4ThreeVector prePos = step->GetPreStepPoint()->GetPosition();
117 G4ThreeVector deltaPos = step->GetPostStepPoint()->GetPosition() - prePos;
118 for(G4int i=0; i<nion; ++i) {
119 v->push_back( prePos + deltaPos*G4UniformRand() );
120 }
121 if(verbose > 1 ) {
122 G4cout << "### G4ElectronIonPair::SampleIonisationPoints: "
123 << v->size() << " ion pairs are added" << G4endl;
124 }
125 }
126 return v;
127}
#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 145 of file G4ElectronIonPair.hh.

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

Referenced by SampleIonsAlongStep().

◆ SetVerbose()

void G4ElectronIonPair::SetVerbose ( G4int val)
inline

Definition at line 164 of file G4ElectronIonPair.hh.

165{
166 verbose = val;
167}

Referenced by G4LossTableManager::ResetParameters().


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