Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IonDEDXScalingICRU73.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 source file
30//
31// Class: G4IonDEDXScalingICRU73
32//
33// Base class: G4VIonDEDXScalingAlgorithm
34//
35// Author: Anton Lechner ([email protected])
36//
37// First implementation: 10. 05. 2009
38//
39// Modifications: 06. 08. 2009 - Minor bug fix (initialization of cache) AL
40// 12. 11. 2009 - Moved all decision logic concerning ICRU 73
41// scaling for heavy ions into this class.
42// Adapting ScalingFactorEnergy class according
43// to changes in base class (AL).
44//
45// Class description:
46// dE/dx scaling algorithm applied on top of ICRU 73 data (for ions not
47// covered by the ICRU 73 report)
48//
49// Comments:
50//
51// ===========================================================================
52
55#include "G4ParticleTable.hh"
56#include "G4Material.hh"
57
58
59// ###########################################################################
60
62 G4int minAtomicNumberIon,
63 G4int maxAtomicNumberIon) :
64 minAtomicNumber( minAtomicNumberIon ),
65 maxAtomicNumber( maxAtomicNumberIon ),
66 referenceFe( 0 ),
67 atomicNumberRefFe( 26 ),
68 massNumberRefFe( 56 ),
69 atomicNumberRefPow23Fe( 0 ),
70 chargeRefFe( 0 ),
71 massRefFe( 0 ),
72 referenceAr( 0 ),
73 atomicNumberRefAr( 18 ),
74 massNumberRefAr( 40 ),
75 atomicNumberRefPow23Ar( 0 ),
76 chargeRefAr( 0 ),
77 massRefAr( 0 ),
78 useFe( true ),
79 cacheParticle( 0 ),
80 cacheMassNumber( 0 ),
81 cacheAtomicNumber( 0 ),
82 cacheAtomicNumberPow23( 0 ),
83 cacheCharge( 0 ),
84 cacheMass( 0 ),
85 cacheMaterial( 0 ) {
86
87}
88
89// ###########################################################################
90
92
93}
94
95// ###########################################################################
96
97void G4IonDEDXScalingICRU73::CreateReferenceParticles() {
98
100
101 G4double excitationEnergy = 0.0;
102
103 referenceFe = particleTable -> GetIon(atomicNumberRefFe, massNumberRefFe,
104 excitationEnergy);
105 referenceAr = particleTable -> GetIon(atomicNumberRefAr, massNumberRefAr,
106 excitationEnergy);
107
108 massRefFe = referenceFe -> GetPDGMass();
109 massRefAr = referenceAr -> GetPDGMass();
110
111 chargeRefFe = referenceFe -> GetPDGCharge();
112 chargeRefAr = referenceAr -> GetPDGCharge();
113
114 atomicNumberRefPow23Fe = std::pow(G4double(atomicNumberRefFe), 2./3.);
115 atomicNumberRefPow23Ar = std::pow(G4double(atomicNumberRefAr), 2./3.);
116}
117
118
119// ###########################################################################
120
121
123 const G4ParticleDefinition* particle, // Projectile (ion)
124 const G4Material* material) { // Target material
125
126 G4double factor = 1.0;
127
128 UpdateCacheParticle(particle);
129 UpdateCacheMaterial(material);
130
131 if(cacheAtomicNumber >= minAtomicNumber &&
132 cacheAtomicNumber <= maxAtomicNumber &&
133 cacheAtomicNumber != atomicNumberRefFe &&
134 cacheAtomicNumber != atomicNumberRefAr) {
135
136 if(referenceFe == 0 || referenceAr == 0) CreateReferenceParticles();
137
138 if( useFe )
139 factor = cacheMassNumber * (massRefFe / cacheMass) / massNumberRefFe;
140 else
141 factor = cacheMassNumber * (massRefAr / cacheMass) / massNumberRefAr;
142 }
143
144 return factor;
145}
146
147// ###########################################################################
148
150 const G4ParticleDefinition* particle, // Projectile (ion)
151 const G4Material* material, // Target material
152 G4double kineticEnergy) { // Kinetic energy
153
154 G4double factor = 1.0;
155
156 UpdateCacheParticle(particle);
157 UpdateCacheMaterial(material);
158
159 if(cacheAtomicNumber >= minAtomicNumber &&
160 cacheAtomicNumber <= maxAtomicNumber &&
161 cacheAtomicNumber != atomicNumberRefFe &&
162 cacheAtomicNumber != atomicNumberRefAr) {
163
164 if(referenceFe == 0 || referenceAr == 0) CreateReferenceParticles();
165
166 if( useFe ) {
167
168 G4double equilibriumCharge = EquilibriumCharge(cacheMass,
169 cacheCharge,
170 cacheAtomicNumberPow23,
171 kineticEnergy);
172
173 G4double scaledKineticEnergy = kineticEnergy * (massRefFe / cacheMass);
174
175 G4double equilibriumChargeRefFe = EquilibriumCharge(massRefFe,
176 chargeRefFe,
177 atomicNumberRefPow23Fe,
178 scaledKineticEnergy);
179
180 factor = equilibriumCharge * equilibriumCharge/
181 ( equilibriumChargeRefFe * equilibriumChargeRefFe );
182
183 }
184 else {
185
186 G4double equilibriumCharge = EquilibriumCharge(cacheMass,
187 cacheCharge,
188 cacheAtomicNumberPow23,
189 kineticEnergy);
190
191 G4double scaledKineticEnergy = kineticEnergy * (massRefAr / cacheMass);
192
193 G4double equilibriumChargeRefAr = EquilibriumCharge(massRefAr,
194 chargeRefAr,
195 atomicNumberRefPow23Ar,
196 scaledKineticEnergy);
197
198 factor = equilibriumCharge * equilibriumCharge/
199 ( equilibriumChargeRefAr * equilibriumChargeRefAr );
200
201 }
202 }
203
204 return factor;
205}
206
207// ###########################################################################
208
210 G4int atomicNumberIon, // Atomic number of ion
211 const G4Material* material) { // Target material
212
213 UpdateCacheMaterial(material);
214
215 G4int atomicNumber = atomicNumberIon;
216
217 if(atomicNumberIon >= minAtomicNumber &&
218 atomicNumberIon <= maxAtomicNumber &&
219 atomicNumberIon != atomicNumberRefFe &&
220 atomicNumberIon != atomicNumberRefAr) {
221
222 if(referenceFe == 0 || referenceAr == 0) CreateReferenceParticles();
223
224 if( useFe ) atomicNumber = atomicNumberRefFe;
225 else atomicNumber = atomicNumberRefAr;
226 }
227
228 return atomicNumber;
229}
230
231// ###########################################################################
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double ScalingFactorDEDX(const G4ParticleDefinition *particle, const G4Material *, G4double kineticEnergy)
G4int AtomicNumberBaseIon(G4int atomicNumberIon, const G4Material *)
G4IonDEDXScalingICRU73(G4int minAtomicNumberIon=19, G4int maxAtomicNumberIon=102)
G4double ScalingFactorEnergy(const G4ParticleDefinition *particle, const G4Material *material)
static G4ParticleTable * GetParticleTable()