#include <G4EmSaturation.hh>
Definition at line 69 of file G4EmSaturation.hh.
◆ G4EmSaturation() [1/2]
G4EmSaturation::G4EmSaturation |
( |
G4int |
verb | ) |
|
|
explicit |
Definition at line 62 of file G4EmSaturation.cc.
63{
64 verbose = verb;
65 nWarnings = nG4Birks = 0;
66
67 electron = nullptr;
68 proton = nullptr;
71}
void InitialiseG4Saturation()
static G4NistManager * Instance()
◆ ~G4EmSaturation()
G4EmSaturation::~G4EmSaturation |
( |
| ) |
|
|
virtualdefault |
◆ G4EmSaturation() [2/2]
◆ DumpBirksCoefficients()
void G4EmSaturation::DumpBirksCoefficients |
( |
| ) |
|
Definition at line 237 of file G4EmSaturation.cc.
238{
239 G4cout <<
"### Birks coefficients used in run time" <<
G4endl;
241 for(std::size_t i=0; i<nMaterials; ++i) {
244 if(br > 0.0) {
246 << br*MeV/mm << " mm/MeV" << " "
248 << " g/cm^2/MeV massFactor= " << massFactors[i]
249 <<
" effCharge= " << effCharges[i] <<
G4endl;
250 }
251 }
252}
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 256 of file G4EmSaturation.cc.
257{
258 if(nG4Birks > 0) {
259 G4cout <<
"### Birks coefficients for Geant4 materials" <<
G4endl;
260 for(
G4int i=0; i<nG4Birks; ++i) {
261 G4cout <<
" " << g4MatNames[i] <<
" "
262 << g4MatData[i]*MeV/mm <<
" mm/MeV" <<
G4endl;
263 }
264 }
265}
◆ FindG4BirksCoefficient()
Definition at line 160 of file G4EmSaturation.cc.
161{
162 if(0 == nG4Birks) { InitialiseG4materials(); }
163
165
166
167 for(
G4int j=0; j<nG4Birks; ++j) {
168 if(name == g4MatNames[j]) {
169 if(verbose > 0)
170 G4cout <<
"### G4EmSaturation::FindG4BirksCoefficient for "
171 <<
name <<
" is " << g4MatData[j]*MeV/mm <<
" mm/MeV "
173 return g4MatData[j];
174 }
175 }
176 return 0.0;
177}
const char * name(G4int ptype)
◆ InitialiseG4Saturation()
void G4EmSaturation::InitialiseG4Saturation |
( |
| ) |
|
Definition at line 143 of file G4EmSaturation.cc.
144{
147 massFactors.resize(nMaterials, 1.0);
148 effCharges.resize(nMaterials, 1.0);
149
150 if(0 == nG4Birks) { InitialiseG4materials(); }
151
152 for(std::size_t i=0; i<nMaterials; ++i) {
154 }
156}
void DumpBirksCoefficients()
static size_t GetNumberOfMaterials()
Referenced by G4EmSaturation().
◆ operator=()
◆ SetVerbose()
void G4EmSaturation::SetVerbose |
( |
G4int |
val | ) |
|
|
inline |
◆ VisibleEnergyDeposition()
Definition at line 79 of file G4EmSaturation.cc.
85{
86
87 if(edep <= 0.0) { return 0.0; }
88
89
90
91 if(length <= 0.0) { return edep; }
92
95
96 if(bfactor > 0.0) {
97
98
100
101 evis /= (1.0 + bfactor*edep/
103
104
105 } else {
106
107
108 G4double nloss = std::max(niel, 0.0);
110
111
113 nloss = edep;
114 eloss = 0.0;
115 } else {
116
117
118 eloss /= (1.0 + bfactor*eloss/length);
119 }
120
121 if(nloss > 0.0) {
123 G4double escaled = nloss*massFactors[idx];
124
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 139 of file G4EmSaturation.hh.
141{
147}
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: