Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmSaturation Class Reference

#include <G4EmSaturation.hh>

Public Member Functions

 G4EmSaturation (G4int verb)
 
virtual ~G4EmSaturation ()
 
virtual G4double VisibleEnergyDeposition (const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0) const
 
void InitialiseG4Saturation ()
 
G4double FindG4BirksCoefficient (const G4Material *)
 
void DumpBirksCoefficients ()
 
void DumpG4BirksCoefficients ()
 
G4double VisibleEnergyDepositionAtAStep (const G4Step *) const
 
void SetVerbose (G4int)
 

Detailed Description

Definition at line 71 of file G4EmSaturation.hh.

Constructor & Destructor Documentation

G4EmSaturation::G4EmSaturation ( G4int  verb)
explicit

Definition at line 63 of file G4EmSaturation.cc.

64 {
65  verbose = verb;
66 
67  nWarnings = nG4Birks = 0;
68 
69  electron = nullptr;
70  proton = nullptr;
71  nist = G4NistManager::Instance();
72 }
static G4NistManager * Instance()

Here is the call graph for this function:

G4EmSaturation::~G4EmSaturation ( )
virtual

Definition at line 76 of file G4EmSaturation.cc.

77 {}

Member Function Documentation

void G4EmSaturation::DumpBirksCoefficients ( )

Definition at line 232 of file G4EmSaturation.cc.

233 {
234  G4cout << "### Birks coeffitients used in run time" << G4endl;
236  for(G4int i=0; i<nMaterials; ++i) {
237  const G4Material* mat = (*mtable)[i];
238  G4double br = mat->GetIonisation()->GetBirksConstant();
239  if(br > 0.0) {
240  G4cout << " " << mat->GetName() << " "
241  << br*MeV/mm << " mm/MeV" << " "
242  << br*mat->GetDensity()*MeV*cm2/g
243  << " g/cm^2/MeV massFactor= " << massFactors[i]
244  << " effCharge= " << effCharges[i] << G4endl;
245  }
246  }
247 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static constexpr double mm
Definition: G4SIunits.hh:115
static constexpr double cm2
Definition: G4SIunits.hh:120
const G4String & GetName() const
Definition: G4Material.hh:178
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4double GetDensity() const
Definition: G4Material.hh:180
static constexpr double g
Definition: G4SIunits.hh:183
int G4int
Definition: G4Types.hh:78
G4double GetBirksConstant() const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmSaturation::DumpG4BirksCoefficients ( )

Definition at line 251 of file G4EmSaturation.cc.

252 {
253  if(nG4Birks > 0) {
254  G4cout << "### Birks coeffitients for Geant4 materials" << G4endl;
255  for(G4int i=0; i<nG4Birks; ++i) {
256  G4cout << " " << g4MatNames[i] << " "
257  << g4MatData[i]*MeV/mm << " mm/MeV" << G4endl;
258  }
259  }
260 }
static constexpr double mm
Definition: G4SIunits.hh:115
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double G4EmSaturation::FindG4BirksCoefficient ( const G4Material mat)

Definition at line 155 of file G4EmSaturation.cc.

156 {
157  if(0 == nG4Birks) { InitialiseG4materials(); }
158 
159  G4String name = mat->GetName();
160  // is this material in the vector?
161 
162  for(G4int j=0; j<nG4Birks; ++j) {
163  if(name == g4MatNames[j]) {
164  if(verbose > 0)
165  G4cout << "### G4EmSaturation::FindG4BirksCoefficient for "
166  << name << " is " << g4MatData[j]*MeV/mm << " mm/MeV "
167  << G4endl;
168  return g4MatData[j];
169  }
170  }
171  return 0.0;
172 }
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
const G4String & GetName() const
Definition: G4Material.hh:178
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

void G4EmSaturation::InitialiseG4Saturation ( )

Definition at line 139 of file G4EmSaturation.cc.

140 {
141  nMaterials = G4Material::GetNumberOfMaterials();
142  massFactors.resize(nMaterials, 1.0);
143  effCharges.resize(nMaterials, 1.0);
144 
145  if(0 == nG4Birks) { InitialiseG4materials(); }
146 
147  for(G4int i=0; i<nMaterials; ++i) {
148  InitialiseBirksCoefficient((*G4Material::GetMaterialTable())[i]);
149  }
150  if(verbose > 0) { DumpBirksCoefficients(); }
151 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
void DumpBirksCoefficients()
int G4int
Definition: G4Types.hh:78
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmSaturation::SetVerbose ( G4int  val)
inline

Definition at line 134 of file G4EmSaturation.hh.

135 {
136  verbose = val;
137 }
G4double G4EmSaturation::VisibleEnergyDeposition ( const G4ParticleDefinition p,
const G4MaterialCutsCouple couple,
G4double  length,
G4double  edepTotal,
G4double  edepNIEL = 0.0 
) const
virtual

Definition at line 81 of file G4EmSaturation.cc.

87 {
88  if(edep <= 0.0) { return 0.0; }
89 
90  G4double evis = edep;
91  G4double bfactor = couple->GetMaterial()->GetIonisation()->GetBirksConstant();
92 
93  if(bfactor > 0.0) {
94 
95  // atomic relaxations for gamma incident
96  if(22 == p->GetPDGEncoding()) {
97  //G4cout << "%% gamma edep= " << edep/keV << " keV " <<manager << G4endl;
98  evis /= (1.0 + bfactor*edep/
99  G4LossTableManager::Instance()->GetRange(electron,edep,couple));
100 
101  // energy loss
102  } else {
103 
104  // protections
105  G4double nloss = std::max(niel, 0.0);
106  G4double eloss = edep - nloss;
107 
108  // neutrons and neutral hadrons
109  if(0.0 == p->GetPDGCharge() || eloss < 0.0 || length <= 0.0) {
110  nloss = edep;
111  eloss = 0.0;
112  } else {
113 
114  // continues energy loss
115  eloss /= (1.0 + bfactor*eloss/length);
116  }
117  // non-ionizing energy loss
118  if(nloss > 0.0) {
119  G4int idx = couple->GetMaterial()->GetIndex();
120  G4double escaled = nloss*massFactors[idx];
121  /*
122  G4cout << "%% p edep= " << nloss/keV << " keV Escaled= "
123  << escaled << " MeV in " << couple->GetMaterial()->GetName()
124  << " " << p->GetParticleName()
125  << G4endl;
126  */
128  ->GetRange(proton,escaled,couple)/effCharges[idx];
129  nloss /= (1.0 + bfactor*nloss/range);
130  }
131  evis = eloss + nloss;
132  }
133  }
134  return evis;
135 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static G4LossTableManager * Instance()
size_t GetIndex() const
Definition: G4Material.hh:262
int G4int
Definition: G4Types.hh:78
G4double GetBirksConstant() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmSaturation::VisibleEnergyDepositionAtAStep ( const G4Step step) const
inline

Definition at line 141 of file G4EmSaturation.hh.

143 {
145  step->GetTrack()->GetMaterialCutsCouple(),
146  step->GetStepLength(),
147  step->GetTotalEnergyDeposit(),
149 }
G4double GetStepLength() const
G4double GetNonIonizingEnergyDeposit() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double VisibleEnergyDeposition(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0) const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetTotalEnergyDeposit() const
G4Track * GetTrack() const

Here is the call graph for this function:

Here is the caller graph for this function:


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