#include <G4EmSaturation.hh>
Definition at line 71 of file G4EmSaturation.hh.
G4EmSaturation::G4EmSaturation |
( |
G4int |
verb | ) |
|
|
explicit |
Definition at line 63 of file G4EmSaturation.cc.
67 nWarnings = nG4Birks = 0;
static G4NistManager * Instance()
G4EmSaturation::~G4EmSaturation |
( |
| ) |
|
|
virtual |
void G4EmSaturation::DumpBirksCoefficients |
( |
| ) |
|
Definition at line 232 of file G4EmSaturation.cc.
234 G4cout <<
"### Birks coeffitients used in run time" <<
G4endl;
236 for(
G4int i=0; i<nMaterials; ++i) {
241 << br*
MeV/
mm <<
" mm/MeV" <<
" "
243 <<
" g/cm^2/MeV massFactor= " << massFactors[i]
244 <<
" effCharge= " << effCharges[i] <<
G4endl;
G4IonisParamMat * GetIonisation() const
static constexpr double mm
static constexpr double cm2
const G4String & GetName() const
static G4MaterialTable * GetMaterialTable()
std::vector< G4Material * > G4MaterialTable
G4double GetDensity() const
static constexpr double g
G4double GetBirksConstant() const
G4GLOB_DLL std::ostream G4cout
static constexpr double MeV
void G4EmSaturation::DumpG4BirksCoefficients |
( |
| ) |
|
Definition at line 251 of file G4EmSaturation.cc.
254 G4cout <<
"### Birks coeffitients for Geant4 materials" <<
G4endl;
255 for(
G4int i=0; i<nG4Birks; ++i) {
256 G4cout <<
" " << g4MatNames[i] <<
" "
static constexpr double mm
G4GLOB_DLL std::ostream G4cout
static constexpr double MeV
Definition at line 155 of file G4EmSaturation.cc.
157 if(0 == nG4Birks) { InitialiseG4materials(); }
162 for(
G4int j=0; j<nG4Birks; ++j) {
163 if(name == g4MatNames[j]) {
165 G4cout <<
"### G4EmSaturation::FindG4BirksCoefficient for "
166 << name <<
" is " << g4MatData[j]*
MeV/
mm <<
" mm/MeV "
static constexpr double mm
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
static constexpr double MeV
void G4EmSaturation::InitialiseG4Saturation |
( |
| ) |
|
Definition at line 139 of file G4EmSaturation.cc.
142 massFactors.resize(nMaterials, 1.0);
143 effCharges.resize(nMaterials, 1.0);
145 if(0 == nG4Birks) { InitialiseG4materials(); }
147 for(
G4int i=0; i<nMaterials; ++i) {
static G4MaterialTable * GetMaterialTable()
void DumpBirksCoefficients()
static size_t GetNumberOfMaterials()
Definition at line 81 of file G4EmSaturation.cc.
88 if(edep <= 0.0) {
return 0.0; }
98 evis /= (1.0 + bfactor*edep/
109 if(0.0 == p->
GetPDGCharge() || eloss < 0.0 || length <= 0.0) {
115 eloss /= (1.0 + bfactor*eloss/length);
120 G4double escaled = nloss*massFactors[idx];
128 ->
GetRange(proton,escaled,couple)/effCharges[idx];
129 nloss /= (1.0 + bfactor*nloss/range);
131 evis = eloss + nloss;
G4IonisParamMat * GetIonisation() const
static G4LossTableManager * Instance()
G4int GetPDGEncoding() const
G4double GetBirksConstant() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetPDGCharge() const
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const
G4double G4EmSaturation::VisibleEnergyDepositionAtAStep |
( |
const G4Step * |
step | ) |
const |
|
inline |
Definition at line 141 of file G4EmSaturation.hh.
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
The documentation for this class was generated from the following files: