source: trunk/source/processes/electromagnetic/lowenergy/src/G4hNuclearStoppingModel.cc@ 1197

Last change on this file since 1197 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.5 KB
Line 
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 (Vladimir.Ivanchenko@cern.ch)
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
54#include "G4hNuclearStoppingModel.hh"
55#include "G4UnitsTable.hh"
56#include "globals.hh"
57#include "G4hZiegler1977Nuclear.hh"
58#include "G4hZiegler1985Nuclear.hh"
59#include "G4hICRU49Nuclear.hh"
60#include "G4DynamicParticle.hh"
61#include "G4ParticleDefinition.hh"
62#include "G4ElementVector.hh"
63#include "G4Material.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67G4hNuclearStoppingModel::G4hNuclearStoppingModel(const G4String& name)
68 :G4VLowEnergyModel(name), modelName(name)
69{
70 InitializeMe() ;
71}
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75void G4hNuclearStoppingModel::InitializeMe()
76{
77 // Constants
78 highEnergyLimit = 100.*MeV ;
79 lowEnergyLimit = 1.*eV ;
80 factorPDG2AMU = 1.007276/proton_mass_c2 ;
81 theZieglerFactor= eV*cm2*1.0e-15 ;
82
83 // Registration of parametrisation models of nuclear energy losses
84 G4String blank = G4String(" ") ;
85 G4String zi77 = G4String("Ziegler1977") ;
86 G4String ir49 = G4String("ICRU_R49") ;
87 G4String zi85 = G4String("Ziegler1985") ;
88 if(zi77 == modelName) {
89 nStopingPowerTable = new G4hZiegler1977Nuclear();
90
91 } else if(ir49 == modelName || blank == modelName) {
92 nStopingPowerTable = new G4hICRU49Nuclear();
93
94 } else if(zi85 == modelName) {
95 nStopingPowerTable = new G4hZiegler1985Nuclear();
96
97 } else {
98 G4cout <<
99 "G4hLowEnergyIonisation warning: There is no table with the modelName <"
100 << modelName << ">"
101 << " for nuclear stopping, <ICRU_R49> is applied "
102 << G4endl;
103 nStopingPowerTable = new G4hICRU49Nuclear();
104 }
105
106 // Default is nuclear stopping fluctuations On
107 // nStopingPowerTable->SetNuclearStoppingFluctuationsOn();
108 nStopingPowerTable->SetNuclearStoppingFluctuationsOff();
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112
113G4hNuclearStoppingModel::~G4hNuclearStoppingModel()
114{
115 delete nStopingPowerTable;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120G4double G4hNuclearStoppingModel::TheValue(
121 const G4DynamicParticle* particle,
122 const G4Material* material)
123{
124 // Projectile nucleus
125 G4double energy = particle->GetKineticEnergy() ;
126 G4double z1 = std::abs((particle->GetCharge())/eplus) ;
127 G4double m1 = (particle->GetMass())*factorPDG2AMU ;
128
129 G4double nloss = StoppingPower(material, energy, z1, m1) * theZieglerFactor;
130
131 return nloss;
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135
136G4double G4hNuclearStoppingModel::TheValue(
137 const G4ParticleDefinition* aParticle,
138 const G4Material* material,
139 G4double kineticEnergy)
140{
141 // Projectile nucleus
142 G4double z1 = std::abs((aParticle->GetPDGCharge())/eplus) ;
143 G4double m1 = (aParticle->GetPDGMass())*factorPDG2AMU ;
144
145 G4double nloss = StoppingPower(material, kineticEnergy, z1, m1)
146 * theZieglerFactor;
147
148 return nloss;
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152
153G4double G4hNuclearStoppingModel::StoppingPower(
154 const G4Material* material,
155 G4double kineticEnergy,
156 G4double z1, G4double m1) const
157{
158 // Target nucleus
159 G4int NumberOfElements = material->GetNumberOfElements() ;
160 if(0 == NumberOfElements) return 0.0 ;
161
162 const G4ElementVector* theElementVector =
163 material->GetElementVector() ;
164 const G4double* theAtomicNumDensityVector =
165 material->GetAtomicNumDensityVector() ;
166
167 // loop for the elements in the material
168
169 G4double nloss = 0.0;
170
171 for (G4int iel=0; iel<NumberOfElements; iel++) {
172 const G4Element* element = (*theElementVector)[iel] ;
173 G4double z2 = element->GetZ();
174 G4double m2 = element->GetA()*mole/g ;
175 nloss += (nStopingPowerTable->
176 NuclearStoppingPower(kineticEnergy, z1, z2, m1, m2))
177 * theAtomicNumDensityVector[iel] ;
178 }
179
180 return nloss;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184
185
186
187
188
189
190
191
Note: See TracBrowser for help on using the repository browser.