Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4hNuclearStoppingModel.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: G4hNuclearStoppingModel
33//
34// Author: V.Ivanchenko ([email protected])
35//
36// Creation date: 20 July 2000
37//
38// Modifications:
39// 20/07/2000 V.Ivanchenko First implementation
40// 22/08/2000 V.Ivanchenko Bug fixed in call of a model
41// 03/10/2000 V.Ivanchenko CodeWizard clean up
42//
43// Class Description:
44//
45// Low energy protons/ions nuclear stopping parametrisation
46//
47// Class Description: End
48//
49// -------------------------------------------------------------------
50//
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
55
56#include "globals.hh"
58#include "G4SystemOfUnits.hh"
59#include "G4UnitsTable.hh"
62#include "G4hICRU49Nuclear.hh"
63#include "G4DynamicParticle.hh"
65#include "G4ElementVector.hh"
66#include "G4Material.hh"
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
71 :G4VLowEnergyModel(name), modelName(name)
72{
73 InitializeMe() ;
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78void G4hNuclearStoppingModel::InitializeMe()
79{
80 // Constants
81 highEnergyLimit = 100.*MeV ;
82 lowEnergyLimit = 1.*eV ;
83 factorPDG2AMU = 1.007276/proton_mass_c2 ;
84 theZieglerFactor= eV*cm2*1.0e-15 ;
85
86 // Registration of parametrisation models of nuclear energy losses
87 G4String blank = G4String(" ") ;
88 G4String zi77 = G4String("Ziegler1977") ;
89 G4String ir49 = G4String("ICRU_R49") ;
90 G4String zi85 = G4String("Ziegler1985") ;
91 if(zi77 == modelName) {
92 nStopingPowerTable = new G4hZiegler1977Nuclear();
93
94 } else if(ir49 == modelName || blank == modelName) {
95 nStopingPowerTable = new G4hICRU49Nuclear();
96
97 } else if(zi85 == modelName) {
98 nStopingPowerTable = new G4hZiegler1985Nuclear();
99
100 } else {
101 G4cout <<
102 "G4hLowEnergyIonisation warning: There is no table with the modelName <"
103 << modelName << ">"
104 << " for nuclear stopping, <ICRU_R49> is applied "
105 << G4endl;
106 nStopingPowerTable = new G4hICRU49Nuclear();
107 }
108
109 // Default is nuclear stopping fluctuations On
110 // nStopingPowerTable->SetNuclearStoppingFluctuationsOn();
111 nStopingPowerTable->SetNuclearStoppingFluctuationsOff();
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115
117{
118 delete nStopingPowerTable;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
124 const G4DynamicParticle* particle,
125 const G4Material* material)
126{
127 // Projectile nucleus
128 G4double energy = particle->GetKineticEnergy() ;
129 G4double z1 = std::abs((particle->GetCharge())/eplus) ;
130 G4double m1 = (particle->GetMass())*factorPDG2AMU ;
131
132 G4double nloss = StoppingPower(material, energy, z1, m1) * theZieglerFactor;
133
134 return nloss;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
140 const G4ParticleDefinition* aParticle,
141 const G4Material* material,
142 G4double kineticEnergy)
143{
144 // Projectile nucleus
145 G4double z1 = std::abs((aParticle->GetPDGCharge())/eplus) ;
146 G4double m1 = (aParticle->GetPDGMass())*factorPDG2AMU ;
147
148 G4double nloss = StoppingPower(material, kineticEnergy, z1, m1)
149 * theZieglerFactor;
150
151 return nloss;
152}
153
154//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155
156G4double G4hNuclearStoppingModel::StoppingPower(
157 const G4Material* material,
158 G4double kineticEnergy,
159 G4double z1, G4double m1) const
160{
161 // Target nucleus
162 G4int NumberOfElements = material->GetNumberOfElements() ;
163 if(0 == NumberOfElements) return 0.0 ;
164
165 const G4ElementVector* theElementVector =
166 material->GetElementVector() ;
167 const G4double* theAtomicNumDensityVector =
168 material->GetAtomicNumDensityVector() ;
169
170 // loop for the elements in the material
171
172 G4double nloss = 0.0;
173
174 for (G4int iel=0; iel<NumberOfElements; iel++) {
175 const G4Element* element = (*theElementVector)[iel] ;
176 G4double z2 = element->GetZ();
177 G4double m2Local = element->GetA()*mole/g ;
178 nloss += (nStopingPowerTable->
179 NuclearStoppingPower(kineticEnergy, z1, z2, m1, m2Local))
180 * theAtomicNumDensityVector[iel] ;
181 }
182
183 return nloss;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
188
189
190
191
192
193
194
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetMass() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetA() const
Definition: G4Element.hh:138
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
G4double GetPDGCharge() const
G4hNuclearStoppingModel(const G4String &name)
G4double TheValue(const G4DynamicParticle *particle, const G4Material *material)