#include <G4EmSaturation.hh>
Definition at line 70 of file G4EmSaturation.hh.
◆ G4EmSaturation()
G4EmSaturation::G4EmSaturation |
( |
G4int |
verb | ) |
|
|
explicit |
Definition at line 62 of file G4EmSaturation.cc.
63{
64 verbose = verb;
65
66 nWarnings = nG4Birks = 0;
67
68 electron = nullptr;
69 proton = nullptr;
71}
static G4NistManager * Instance()
◆ ~G4EmSaturation()
G4EmSaturation::~G4EmSaturation |
( |
| ) |
|
|
virtual |
◆ DumpBirksCoefficients()
void G4EmSaturation::DumpBirksCoefficients |
( |
| ) |
|
Definition at line 236 of file G4EmSaturation.cc.
237{
238 G4cout <<
"### Birks coefficients used in run time" <<
G4endl;
240 for(
G4int i=0; i<nMaterials; ++i) {
243 if(br > 0.0) {
245 << br*MeV/mm << " mm/MeV" << " "
247 << " g/cm^2/MeV massFactor= " << massFactors[i]
248 <<
" effCharge= " << effCharges[i] <<
G4endl;
249 }
250 }
251}
std::vector< G4Material * > G4MaterialTable
G4GLOB_DLL std::ostream G4cout
G4double GetBirksConstant() const
G4double GetDensity() const
G4IonisParamMat * GetIonisation() const
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const
Referenced by InitialiseG4Saturation().
◆ DumpG4BirksCoefficients()
void G4EmSaturation::DumpG4BirksCoefficients |
( |
| ) |
|
Definition at line 255 of file G4EmSaturation.cc.
256{
257 if(nG4Birks > 0) {
258 G4cout <<
"### Birks coefficients for Geant4 materials" <<
G4endl;
259 for(
G4int i=0; i<nG4Birks; ++i) {
260 G4cout <<
" " << g4MatNames[i] <<
" "
261 << g4MatData[i]*MeV/mm <<
" mm/MeV" <<
G4endl;
262 }
263 }
264}
◆ FindG4BirksCoefficient()
Definition at line 159 of file G4EmSaturation.cc.
160{
161 if(0 == nG4Birks) { InitialiseG4materials(); }
162
164
165
166 for(
G4int j=0; j<nG4Birks; ++j) {
167 if(name == g4MatNames[j]) {
168 if(verbose > 0)
169 G4cout <<
"### G4EmSaturation::FindG4BirksCoefficient for "
170 <<
name <<
" is " << g4MatData[j]*MeV/mm <<
" mm/MeV "
172 return g4MatData[j];
173 }
174 }
175 return 0.0;
176}
const char * name(G4int ptype)
◆ InitialiseG4Saturation()
void G4EmSaturation::InitialiseG4Saturation |
( |
| ) |
|
Definition at line 143 of file G4EmSaturation.cc.
144{
146 massFactors.resize(nMaterials, 1.0);
147 effCharges.resize(nMaterials, 1.0);
148
149 if(0 == nG4Birks) { InitialiseG4materials(); }
150
151 for(
G4int i=0; i<nMaterials; ++i) {
153 }
155}
void DumpBirksCoefficients()
static size_t GetNumberOfMaterials()
Referenced by G4EmParameters::SetBirksActive().
◆ SetVerbose()
void G4EmSaturation::SetVerbose |
( |
G4int |
val | ) |
|
|
inline |
◆ VisibleEnergyDeposition()
Definition at line 80 of file G4EmSaturation.cc.
86{
87
88 if(edep <= 0.0) { return 0.0; }
89
90
91
92 if(length <= 0.0) { return edep; }
93
96
97 if(bfactor > 0.0) {
98
99
101
102 evis /= (1.0 + bfactor*edep/
104
105
106 } else {
107
108
109 G4double nloss = std::max(niel, 0.0);
111
112
114 nloss = edep;
115 eloss = 0.0;
116 } else {
117
118
119 eloss /= (1.0 + bfactor*eloss/length);
120 }
121
122 if(nloss > 0.0) {
124 G4double escaled = nloss*massFactors[idx];
125
126
127
128
129
130
132 ->
GetRange(proton,escaled,couple)/effCharges[idx];
133 nloss /= (1.0 + bfactor*nloss/range);
134 }
135 evis = eloss + nloss;
136 }
137 }
138 return evis;
139}
static G4LossTableManager * Instance()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const
G4int GetPDGEncoding() const
G4double GetPDGCharge() const
Referenced by VisibleEnergyDepositionAtAStep().
◆ VisibleEnergyDepositionAtAStep()
G4double G4EmSaturation::VisibleEnergyDepositionAtAStep |
( |
const G4Step * |
step | ) |
const |
|
inline |
Definition at line 140 of file G4EmSaturation.hh.
142{
148}
virtual G4double VisibleEnergyDeposition(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0) const
G4Track * GetTrack() const
G4double GetNonIonizingEnergyDeposit() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
Referenced by G4Scintillation::PostStepDoIt().
The documentation for this class was generated from the following files: