Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4EmSaturation Class Reference

#include <G4EmSaturation.hh>

Public Member Functions

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

Detailed Description

Definition at line 69 of file G4EmSaturation.hh.

Constructor & Destructor Documentation

G4EmSaturation::G4EmSaturation ( )

Definition at line 58 of file G4EmSaturation.cc.

References G4NistManager::Instance().

59 {
60  verbose = 1;
61  manager = 0;
62 
63  curMaterial = 0;
64  curBirks = 0.0;
65  curRatio = 1.0;
66  curChargeSq = 1.0;
67  nMaterials = 0;
68 
69  electron = 0;
70  proton = 0;
71  nist = G4NistManager::Instance();
72 
73  Initialise();
74 }
static G4NistManager * Instance()
G4EmSaturation::~G4EmSaturation ( )
virtual

Definition at line 78 of file G4EmSaturation.cc.

79 {}

Member Function Documentation

void G4EmSaturation::DumpBirksCoefficients ( )

Definition at line 233 of file G4EmSaturation.cc.

References python.hepunit::cm2, g(), G4cout, G4endl, python.hepunit::MeV, and python.hepunit::mm.

234 {
235  if(nMaterials > 0) {
236  G4cout << "### Birks coeffitients used in run time" << G4endl;
237  for(G4int i=0; i<nMaterials; ++i) {
238  G4double br = matPointers[i]->GetIonisation()->GetBirksConstant();
239  G4cout << " " << matNames[i] << " "
240  << br*MeV/mm << " mm/MeV" << " "
241  << br*matPointers[i]->GetDensity()*MeV*cm2/g
242  << " g/cm^2/MeV"
243  << G4endl;
244  }
245  }
246 }
int G4int
Definition: G4Types.hh:78
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4EmSaturation::DumpG4BirksCoefficients ( )

Definition at line 250 of file G4EmSaturation.cc.

References G4cout, G4endl, python.hepunit::MeV, and python.hepunit::mm.

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

Definition at line 136 of file G4EmSaturation.cc.

References G4cout, G4endl, G4Material::GetName(), python.hepunit::MeV, and python.hepunit::mm.

137 {
138  G4String name = mat->GetName();
139  // is this material in the vector?
140 
141  for(G4int j=0; j<nG4Birks; ++j) {
142  if(name == g4MatNames[j]) {
143  if(verbose > 0)
144  G4cout << "### G4EmSaturation::FindG4BirksCoefficient for "
145  << name << " is " << g4MatData[j]*MeV/mm << " mm/MeV "
146  << G4endl;
147  return g4MatData[j];
148  }
149  }
150  return FindBirksCoefficient(mat);
151 }
const G4String & GetName() const
Definition: G4Material.hh:176
const XML_Char * name
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4EmSaturation::SetVerbose ( G4int  val)
inline

Definition at line 132 of file G4EmSaturation.hh.

Referenced by G4EmManager::SetVerbose(), and G4LossTableManager::SetVerbose().

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

Definition at line 83 of file G4EmSaturation.cc.

References G4MaterialCutsCouple::GetMaterial(), G4ParticleDefinition::GetPDGEncoding(), G4LossTableManager::GetRange(), and G4Proton::Proton().

Referenced by G4Scintillation::PostStepDoIt(), and VisibleEnergyDeposition().

89 {
90  if(edep <= 0.0) { return 0.0; }
91 
92  G4double evis = edep;
93  G4double bfactor = FindBirksCoefficient(couple->GetMaterial());
94 
95  if(bfactor > 0.0) {
96 
97  G4int pdgCode = p->GetPDGEncoding();
98  // atomic relaxations for gamma incident
99  if(22 == pdgCode) {
100  evis /= (1.0 + bfactor*edep/manager->GetRange(electron,edep,couple));
101 
102  // energy loss
103  } else {
104 
105  // protections
106  G4double nloss = niel;
107  if(nloss < 0.0) nloss = 0.0;
108  G4double eloss = edep - nloss;
109 
110  // neutrons
111  if(2112 == pdgCode || eloss < 0.0 || length <= 0.0) {
112  nloss = edep;
113  eloss = 0.0;
114  }
115 
116  // continues energy loss
117  if(eloss > 0.0) { eloss /= (1.0 + bfactor*eloss/length); }
118 
119  // non-ionizing energy loss
120  if(nloss > 0.0) {
121  if(!proton) { proton = G4Proton::Proton(); }
122  G4double escaled = nloss*curRatio;
123  G4double range = manager->GetRange(proton,escaled,couple)/curChargeSq;
124  nloss /= (1.0 + bfactor*nloss/range);
125  }
126 
127  evis = eloss + nloss;
128  }
129  }
130 
131  return evis;
132 }
int G4int
Definition: G4Types.hh:78
static G4Proton * Proton()
Definition: G4Proton.cc:93
double G4double
Definition: G4Types.hh:76
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const
G4double G4EmSaturation::VisibleEnergyDeposition ( const G4Step step)
inline

Definition at line 137 of file G4EmSaturation.hh.

References G4Track::GetMaterialCutsCouple(), G4Step::GetNonIonizingEnergyDeposit(), G4Track::GetParticleDefinition(), G4Step::GetStepLength(), G4Step::GetTotalEnergyDeposit(), G4Step::GetTrack(), and VisibleEnergyDeposition().

139 {
140  G4Track* track = step->GetTrack();
142  track->GetMaterialCutsCouple(),
143  step->GetStepLength(),
144  step->GetTotalEnergyDeposit(),
146 }
G4double GetStepLength() const
G4double VisibleEnergyDeposition(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0)
G4double GetNonIonizingEnergyDeposit() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetTotalEnergyDeposit() const
G4Track * GetTrack() const

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