Geant4 9.6.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 ()
 
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=0, G4int processSubType=-1)
 
G4int ResidualeChargePostStep (const G4Step *)
 
G4double FindG4MeanEnergyPerIonPair (const G4Material *)
 
void DumpMeanEnergyPerIonPair ()
 
void DumpG4MeanEnergyPerIonPair ()
 
void SetVerbose (G4int)
 

Detailed Description

Definition at line 73 of file G4ElectronIonPair.hh.

Constructor & Destructor Documentation

◆ G4ElectronIonPair()

G4ElectronIonPair::G4ElectronIonPair ( )

Definition at line 58 of file G4ElectronIonPair.cc.

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

◆ ~G4ElectronIonPair()

G4ElectronIonPair::~G4ElectronIonPair ( )
virtual

Definition at line 70 of file G4ElectronIonPair.cc.

71{}

Member Function Documentation

◆ DumpG4MeanEnergyPerIonPair()

void G4ElectronIonPair::DumpG4MeanEnergyPerIonPair ( )

Definition at line 187 of file G4ElectronIonPair.cc.

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

◆ DumpMeanEnergyPerIonPair()

void G4ElectronIonPair::DumpMeanEnergyPerIonPair ( )

Definition at line 168 of file G4ElectronIonPair.cc.

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

◆ FindG4MeanEnergyPerIonPair()

G4double G4ElectronIonPair::FindG4MeanEnergyPerIonPair ( const G4Material mat)

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)

Referenced by MeanNumberOfIonsAlongStep().

◆ MeanNumberOfIonsAlongStep() [1/2]

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

Definition at line 75 of file G4ElectronIonPair.cc.

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

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)
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=0, G4int processSubType=-1)
const G4VProcess * GetProcessDefinedStep() const
const G4TrackVector * GetSecondary() const
G4StepPoint * GetPostStepPoint() const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397

◆ SampleIonsAlongStep()

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

Definition at line 108 of file G4ElectronIonPair.cc.

109{
110 std::vector<G4ThreeVector>* v = 0;
111
113
114 // sample ionisation along step
115 if(nion > 0) {
116
117 v = new std::vector<G4ThreeVector>;
118 G4ThreeVector prePos = step->GetPreStepPoint()->GetPosition();
119 G4ThreeVector deltaPos = step->GetPostStepPoint()->GetPosition() - prePos;
120 for(G4int i=0; i<nion; ++i) {
121 v->push_back( prePos + deltaPos*G4UniformRand() );
122 }
123 if(verbose > 1 ) {
124 G4cout << "### G4ElectronIonPair::SampleIonisationPoints: "
125 << v->size() << " ion pairs are added" << G4endl;
126 }
127 }
128 return v;
129}
#define G4UniformRand()
Definition: Randomize.hh:53
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(CLHEP::RandGamma::shoot(meanion*invFanoFactor,invFanoFactor));
153}
static double shoot()
int G4lrint(double ad)
Definition: templates.hh:163

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::SetVerbose().


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