Geant4-11
G4EmSaturation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4EmSaturation
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 18.02.2008
37//
38// Modifications:
39//
40// -------------------------------------------------------------
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
45#include "G4EmSaturation.hh"
47#include "G4SystemOfUnits.hh"
48#include "G4LossTableManager.hh"
49#include "G4NistManager.hh"
50#include "G4Material.hh"
52#include "G4ParticleTable.hh"
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
57std::vector<G4double> G4EmSaturation::massFactors;
58std::vector<G4double> G4EmSaturation::effCharges;
59std::vector<G4double> G4EmSaturation::g4MatData;
60std::vector<G4String> G4EmSaturation::g4MatNames;
61
63{
64 verbose = verb;
65 nWarnings = nG4Birks = 0;
66
67 electron = nullptr;
68 proton = nullptr;
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
76{}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81 const G4ParticleDefinition* p,
82 const G4MaterialCutsCouple* couple,
83 G4double length,
84 G4double edep,
85 G4double niel) const
86{
87 // no energy deposition
88 if(edep <= 0.0) { return 0.0; }
89
90 // zero step length may happens only if step limiter process
91 // is applied, in that case saturation should not be applied
92 if(length <= 0.0) { return edep; }
93
94 G4double evis = edep;
95 G4double bfactor = couple->GetMaterial()->GetIonisation()->GetBirksConstant();
96
97 if(bfactor > 0.0) {
98
99 // atomic relaxations for gamma incident
100 if(22 == p->GetPDGEncoding()) {
101 //G4cout << "%% gamma edep= " << edep/keV << " keV " << G4endl;
102 evis /= (1.0 + bfactor*edep/
104
105 // energy loss
106 } else {
107
108 // protections
109 G4double nloss = std::max(niel, 0.0);
110 G4double eloss = edep - nloss;
111
112 // neutrons and neutral hadrons
113 if(0.0 == p->GetPDGCharge() || eloss < 0.0) {
114 nloss = edep;
115 eloss = 0.0;
116 } else {
117
118 // continues energy loss
119 eloss /= (1.0 + bfactor*eloss/length);
120 }
121 // non-ionizing energy loss
122 if(nloss > 0.0) {
123 G4int idx = couple->GetMaterial()->GetIndex();
124 G4double escaled = nloss*massFactors[idx];
125 /*
126 G4cout << "%% p edep= " << nloss/keV << " keV Escaled= "
127 << escaled << " MeV in " << couple->GetMaterial()->GetName()
128 << " " << p->GetParticleName()
129 << G4endl;
130 G4cout << proton->GetParticleName() << G4endl;
131 */
133 ->GetRange(proton,escaled,couple)/effCharges[idx];
134 nloss /= (1.0 + bfactor*nloss/range);
135 }
136 evis = eloss + nloss;
137 }
138 }
139 return evis;
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
143
145{
148 massFactors.resize(nMaterials, 1.0);
149 effCharges.resize(nMaterials, 1.0);
150
151 if(0 == nG4Birks) { InitialiseG4materials(); }
152
153 for(size_t i=0; i<nMaterials; ++i) {
155 }
156 if(verbose > 0) { DumpBirksCoefficients(); }
157}
158
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160
162{
163 if(0 == nG4Birks) { InitialiseG4materials(); }
164
165 G4String name = mat->GetName();
166 // is this material in the vector?
167
168 for(G4int j=0; j<nG4Birks; ++j) {
169 if(name == g4MatNames[j]) {
170 if(verbose > 0)
171 G4cout << "### G4EmSaturation::FindG4BirksCoefficient for "
172 << name << " is " << g4MatData[j]*MeV/mm << " mm/MeV "
173 << G4endl;
174 return g4MatData[j];
175 }
176 }
177 return 0.0;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
183{
184 // electron and proton should exist in any case
185 if(nullptr == electron) {
188 if(nullptr == electron) {
189 G4Exception("G4EmSaturation::InitialiseBirksCoefficient", "em0001",
190 FatalException, "electron should exist");
191 }
192 }
193
194 G4double curBirks = mat->GetIonisation()->GetBirksConstant();
195
196 G4String name = mat->GetName();
197
198 // material has no Birks coeffitient defined
199 // seach in the Geant4 list
200 if(curBirks == 0.0) {
201 for(G4int j=0; j<nG4Birks; ++j) {
202 if(name == g4MatNames[j]) {
204 curBirks = g4MatData[j];
205 break;
206 }
207 }
208 }
209
210 if(curBirks == 0.0) { return; }
211
212 // compute mean mass ratio
213 G4double curRatio = 0.0;
214 G4double curChargeSq = 0.0;
215 G4double norm = 0.0;
216 const G4ElementVector* theElementVector = mat->GetElementVector();
217 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
218 size_t nelm = mat->GetNumberOfElements();
219 for (size_t i=0; i<nelm; ++i) {
220 const G4Element* elm = (*theElementVector)[i];
221 G4double Z = elm->GetZ();
222 G4double w = Z*Z*theAtomNumDensityVector[i];
223 curRatio += w/nist->GetAtomicMassAmu(G4int(Z));
224 curChargeSq = Z*Z*w;
225 norm += w;
226 }
227 curRatio *= proton_mass_c2/norm;
228 curChargeSq /= norm;
229
230 // store results
231 G4int idx = mat->GetIndex();
232 massFactors[idx] = curRatio;
233 effCharges[idx] = curChargeSq;
234}
235
236//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
237
239{
240 G4cout << "### Birks coefficients used in run time" << G4endl;
242 for(size_t i=0; i<nMaterials; ++i) {
243 const G4Material* mat = (*mtable)[i];
245 if(br > 0.0) {
246 G4cout << " " << mat->GetName() << " "
247 << br*MeV/mm << " mm/MeV" << " "
248 << br*mat->GetDensity()*MeV*cm2/g
249 << " g/cm^2/MeV massFactor= " << massFactors[i]
250 << " effCharge= " << effCharges[i] << G4endl;
251 }
252 }
253}
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256
258{
259 if(nG4Birks > 0) {
260 G4cout << "### Birks coefficients for Geant4 materials" << G4endl;
261 for(G4int i=0; i<nG4Birks; ++i) {
262 G4cout << " " << g4MatNames[i] << " "
263 << g4MatData[i]*MeV/mm << " mm/MeV" << G4endl;
264 }
265 }
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
269
271{
272 nG4Birks = 4;
273 g4MatData.reserve(nG4Birks);
274
275 // M.Hirschberg et al., IEEE Trans. Nuc. Sci. 39 (1992) 511
276 // SCSN-38 kB = 0.00842 g/cm^2/MeV; rho = 1.06 g/cm^3
277 g4MatNames.push_back("G4_POLYSTYRENE");
278 g4MatData.push_back(0.07943*mm/MeV);
279
280 // C.Fabjan (private communication)
281 // kB = 0.006 g/cm^2/MeV; rho = 7.13 g/cm^3
282 g4MatNames.push_back("G4_BGO");
283 g4MatData.push_back(0.008415*mm/MeV);
284
285 // A.Ribon analysis of publications
286 // Scallettar et al., Phys. Rev. A25 (1982) 2419.
287 // NIM A 523 (2004) 275.
288 // kB = 0.022 g/cm^2/MeV; rho = 1.396 g/cm^3;
289 // ATLAS Efield = 10 kV/cm provide the strongest effect
290 // kB = 0.1576*mm/MeV
291 // A. Kiryunin and P.Strizenec "Geant4 hadronic
292 // working group meeting " kB = 0.041/9.13 g/cm^2/MeV
293 g4MatNames.push_back("G4_lAr");
294 g4MatData.push_back(0.032*mm/MeV);
295
296 //G4_BARIUM_FLUORIDE
297 //G4_CESIUM_IODIDE
298 //G4_GEL_PHOTO_EMULSION
299 //G4_PHOTO_EMULSION
300 //G4_PLASTIC_SC_VINYLTOLUENE
301 //G4_SODIUM_IODIDE
302 //G4_STILBENE
303 //G4_lAr
304
305 //G4_PbWO4 - CMS value
306 g4MatNames.push_back("G4_PbWO4");
307 g4MatData.push_back(0.0333333*mm/MeV);
308
309 //G4_Lucite
310
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
std::vector< const G4Element * > G4ElementVector
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4Material * > G4MaterialTable
static constexpr double cm2
Definition: G4SIunits.hh:100
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetZ() const
Definition: G4Element.hh:131
virtual ~G4EmSaturation()
void InitialiseG4materials()
const G4ParticleDefinition * electron
static std::vector< G4double > g4MatData
virtual G4double VisibleEnergyDeposition(const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double length, G4double edepTotal, G4double edepNIEL=0.0) const
G4double FindG4BirksCoefficient(const G4Material *)
void DumpBirksCoefficients()
G4NistManager * nist
static std::vector< G4String > g4MatNames
static std::vector< G4double > massFactors
void InitialiseG4Saturation()
void DumpG4BirksCoefficients()
G4EmSaturation(G4int verb)
const G4ParticleDefinition * proton
static std::vector< G4double > effCharges
void InitialiseBirksCoefficient(const G4Material *)
static size_t nMaterials
void SetBirksConstant(G4double value)
G4double GetBirksConstant() const
static G4LossTableManager * Instance()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:176
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:679
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:672
const G4String & GetName() const
Definition: G4Material.hh:173
size_t GetIndex() const
Definition: G4Material.hh:256
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
float proton_mass_c2
Definition: hepunit.py:274