Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 74 of file G4ElectronIonPair.hh.

Constructor & Destructor Documentation

G4ElectronIonPair::G4ElectronIonPair ( G4int  verb)
explicit

Definition at line 58 of file G4ElectronIonPair.cc.

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

Definition at line 70 of file G4ElectronIonPair.cc.

71 {}

Member Function Documentation

void G4ElectronIonPair::DumpG4MeanEnergyPerIonPair ( ) const

Definition at line 190 of file G4ElectronIonPair.cc.

191 {
192  if(nMaterials > 0) {
193  G4cout << "### G4ElectronIonPair: mean energy per ion pair "
194  << " for Geant4 materials" << G4endl;
195  for(G4int i=0; i<nMaterials; ++i) {
196  G4cout << " " << g4MatNames[i] << " Epair= "
197  << g4MatData[i]/eV << " eV" << G4endl;
198  }
199  }
200 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
void G4ElectronIonPair::DumpMeanEnergyPerIonPair ( ) const

Definition at line 170 of file G4ElectronIonPair.cc.

171 {
174  if(nmat > 0) {
175  G4cout << "### G4ElectronIonPair: mean energy per ion pair avalable:"
176  << G4endl;
177  for(G4int i=0; i<nmat; ++i) {
178  const G4Material* mat = (*mtable)[i];
180  if(x > 0.0) {
181  G4cout << " " << mat->GetName() << " Epair= "
182  << x/eV << " eV" << G4endl;
183  }
184  }
185  }
186 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetMeanEnergyPerIonPair() const
const G4String & GetName() const
Definition: G4Material.hh:178
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ElectronIonPair::FindG4MeanEnergyPerIonPair ( const G4Material mat) const

Definition at line 147 of file G4ElectronIonPair.cc.

148 {
149  G4String name = mat->GetName();
150  G4double res = 0.0;
151 
152  // is this material in the vector?
153  for(G4int j=0; j<nMaterials; j++) {
154  if(name == g4MatNames[j]) {
155  res = g4MatData[j];
157  if(verbose > 0) {
158  G4cout << "### G4ElectronIonPair::FindG4MeanEnergyPerIonPair for "
159  << name << " Epair= " << res/eV << " eV is set"
160  << G4endl;
161  }
162  break;
163  }
164  }
165  return res;
166 }
const XML_Char * name
Definition: expat.h:151
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
const G4String & GetName() const
Definition: G4Material.hh:178
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetMeanEnergyPerIonPair(G4double value)

Here is the call graph for this function:

Here is the caller graph for this function:

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 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetMeanEnergyPerIonPair() const
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double FindG4MeanEnergyPerIonPair(const G4Material *) const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ElectronIonPair::MeanNumberOfIonsAlongStep ( const G4Step step)
inline

Definition at line 139 of file G4ElectronIonPair.hh.

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

Here is the call graph for this function:

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

Definition at line 134 of file G4ElectronIonPair.cc.

137 {
138  G4int nholes = 0;
139 
140  if(2 == subType || 12 == subType || 13 == subType) { nholes = 1; }
141  return nholes;
142 }
int G4int
Definition: G4Types.hh:78

Here is the caller graph for this function:

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

Definition at line 157 of file G4ElectronIonPair.hh.

158 {
159  G4int subtype = -1;
160  const G4VProcess* proc = step->GetPostStepPoint()->GetProcessDefinedStep();
161  if(proc) { subtype = proc->GetProcessSubType(); }
163  step->GetSecondary(),
164  subtype);
165 }
int G4int
Definition: G4Types.hh:78
G4int ResidualeChargePostStep(const G4ParticleDefinition *, const G4TrackVector *secondary=nullptr, G4int processSubType=-1) const
const G4ParticleDefinition * GetParticleDefinition() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4Track * GetTrack() const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
const G4TrackVector * GetSecondary() const

Here is the call graph for this function:

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

Definition at line 108 of file G4ElectronIonPair.cc.

109 {
110  std::vector<G4ThreeVector>* v = 0;
111 
112  G4int nion = SampleNumberOfIonsAlongStep(step);
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 }
int G4int
Definition: G4Types.hh:78
G4StepPoint * GetPreStepPoint() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPosition() const
G4int SampleNumberOfIonsAlongStep(const G4Step *)
G4StepPoint * GetPostStepPoint() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4int G4ElectronIonPair::SampleNumberOfIonsAlongStep ( const G4Step step)
inline

Definition at line 148 of file G4ElectronIonPair.hh.

149 {
150  // use gamma distribution with mean value n=meanion and
151  // dispersion D=meanion/invFanoFactor
152  G4double meanion = MeanNumberOfIonsAlongStep(step);
153  return G4lrint(G4RandGamma::shoot(meanion*invFanoFactor,invFanoFactor));
154 }
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double MeanNumberOfIonsAlongStep(const G4ParticleDefinition *, const G4Material *, G4double edepTotal, G4double edepNIEL=0.0)
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ElectronIonPair::SetVerbose ( G4int  val)
inline

Definition at line 167 of file G4ElectronIonPair.hh.

168 {
169  verbose = val;
170 }

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