Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4hImpactIonisation.hh
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// G4hImpactIonisation
29//
30// $Id$
31//
32// Author: Maria Grazia Pia (MariaGrazia.Pia@ge.infn.it)
33//
34// 08 Sep 2008 - MGP - Created (initially based on G4hLowEnergyIonisation)
35// Added PIXE capabilities
36// Partial clean-up of the implementation (more needed)
37// Calculation of MicroscopicCrossSection delegated to specialised class
38//
39// ------------------------------------------------------------
40
41// Class Description:
42// Impact Ionisation process of charged hadrons and ions
43// Initially based on G4hLowEnergyIonisation, to be subject to redesign
44// and further evolution of physics capabilities
45//
46// The physics model of G4hLowEnergyIonisation is described in
47// CERN-OPEN-99-121 and CERN-OPEN-99-300.
48//
49// Documentation available in:
50// M.G. Pia et al., PIXE Simulation With Geant4,
51// IEEE Trans. Nucl. Sci., vol. 56, no. 6, pp. 3614-3649, Dec. 2009.
52
53// ------------------------------------------------------------
54
55#ifndef G4HIMPACTIONISATION
56#define G4HIMPACTIONISATION 1
57
58#include <map>
60
61#include "globals.hh"
62#include "G4hRDEnergyLoss.hh"
63#include "G4DataVector.hh"
66
70class G4PhysicsTable;
72class G4Track;
73class G4Step;
74
76{
77public: // With description
78
79 G4hImpactIonisation(const G4String& processName = "hImpactIoni");
80 // The ionisation process for hadrons/ions to be include in the
81 // UserPhysicsList
82
84 // Destructor
85
87 // True for all charged hadrons/ions
88
89 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ;
90 // Build physics table during initialisation
91
92 G4double GetMeanFreePath(const G4Track& track,
93 G4double previousStepSize,
95 // Return MeanFreePath until delta-electron production
96
97 void PrintInfoDefinition() const;
98 // Print out of the class parameters
99
100 void SetHighEnergyForProtonParametrisation(G4double energy) {protonHighEnergy = energy;} ;
101 // Definition of the boundary proton energy. For higher energies
102 // Bethe-Bloch formula is used, for lower energies a parametrisation
103 // of the energy losses is performed. Default is 2 MeV.
104
105 void SetLowEnergyForProtonParametrisation(G4double energy) {protonLowEnergy = energy;} ;
106 // Set of the boundary proton energy. For lower energies
107 // the Free Electron Gas model is used for the energy losses.
108 // Default is 1 keV.
109
110 void SetHighEnergyForAntiProtonParametrisation(G4double energy) {antiprotonHighEnergy = energy;} ;
111 // Set of the boundary antiproton energy. For higher energies
112 // Bethe-Bloch formula is used, for lower energies parametrisation
113 // of the energy losses is performed. Default is 2 MeV.
114
115 void SetLowEnergyForAntiProtonParametrisation(G4double energy) {antiprotonLowEnergy = energy;} ;
116 // Set of the boundary antiproton energy. For lower energies
117 // the Free Electron Gas model is used for the energy losses.
118 // Default is 1 keV.
119
121 G4double previousStepSize,
122 G4double currentMinimumStep,
123 G4double& currentSafety);
124 // Calculation of the step limit due to ionisation losses
125
127 const G4String& dedxTable);
128 // This method defines the electron ionisation parametrisation method
129 // via the name of the table. Default is "ICRU_49p".
130
132 {theNuclearTable = dedxTable; SetNuclearStoppingOn();};
133 // This method defines the nuclear ionisation parametrisation method
134 // via the name of the table. Default is "ICRU_49".
135
136 // ---- MGP ---- The following design of On/Off is nonsense; to be modified
137 // in a following design iteration
138
139 void SetNuclearStoppingOn() {nStopping = true;};
140 // This method switch on calculation of the nuclear stopping power.
141
142 void SetNuclearStoppingOff() {nStopping = false;};
143 // This method switch off calculation of the nuclear stopping power.
144
145 void SetBarkasOn() {theBarkas = true;};
146 // This method switch on calculation of the Barkas and Bloch effects.
147
148 void SetBarkasOff() {theBarkas = false;};
149 // This method switch off calculation of the Barkas and Bloch effects.
150
151 void SetPixe(const G4bool /* val */ ) {pixeIsActive = true;};
152 // This method switches atomic relaxation on/off; currently always on
153
154 G4VParticleChange* AlongStepDoIt(const G4Track& trackData ,
155 const G4Step& stepData ) ;
156 // Function to determine total energy deposition on the step
157
159 const G4Step& Step ) ;
160 // Simulation of delta-ray production.
161
163 const G4MaterialCutsCouple* couple,
164 G4double kineticEnergy);
165 // This method returns electronic dE/dx for protons or antiproton
166
168 // Set threshold energy for fluorescence
169
171 // Set threshold energy for Auger electron production
172
174 // Set Auger electron production flag on/off
175
176 // Accessors to configure PIXE
177 void SetPixeCrossSectionK(const G4String& name) { modelK = name; }
178 void SetPixeCrossSectionL(const G4String& name) { modelL = name; }
179 void SetPixeCrossSectionM(const G4String& name) { modelM = name; }
180 void SetPixeProjectileMinEnergy(G4double energy) { eMinPixe = energy; }
181 void SetPixeProjectileMaxEnergy(G4double energy) { eMaxPixe = energy; }
182
183protected:
184
185private:
186
187 void InitializeMe();
188 void InitializeParametrisation();
189 void BuildLossTable(const G4ParticleDefinition& aParticleType);
190 // void BuildDataForFluorescence(const G4ParticleDefinition& aParticleType);
191 void BuildLambdaTable(const G4ParticleDefinition& aParticleType);
192 void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable)
193 {protonTable = dedxTable ;};
194 // This method defines the ionisation parametrisation method via its name
195
196 void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable)
197 {antiprotonTable = dedxTable;};
198
199 G4double MicroscopicCrossSection(const G4ParticleDefinition& aParticleType,
200 G4double kineticEnergy,
201 G4double atomicNumber,
202 G4double deltaCutInEnergy) const;
203
204 G4double GetConstraints(const G4DynamicParticle* particle,
205 const G4MaterialCutsCouple* couple);
206 // Function to determine StepLimit
207
208 G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
209 G4double kineticEnergy) const;
210
211 G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
212 G4double kineticEnergy) const;
213
214 G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple,
215 G4double kineticEnergy,
216 G4double particleMass) const;
217 // This method returns average energy loss due to delta-rays emission with
218 // energy higher than the cut energy for given material.
219
220 G4double BarkasTerm(const G4Material* material,
221 G4double kineticEnergy) const;
222 // Function to compute the Barkas term for protons
223
224 G4double BlochTerm(const G4Material* material,
225 G4double kineticEnergy,
226 G4double cSquare) const;
227 // Function to compute the Bloch term for protons
228
229 G4double ElectronicLossFluctuation(const G4DynamicParticle* particle,
230 const G4MaterialCutsCouple* material,
231 G4double meanLoss,
232 G4double step) const;
233 // Function to sample electronic losses
234
235 // hide assignment operator
236 G4hImpactIonisation & operator=(const G4hImpactIonisation &right);
238
239private:
240 // private data members ...............................
241 G4VLowEnergyModel* betheBlochModel;
242 G4VLowEnergyModel* protonModel;
243 G4VLowEnergyModel* antiprotonModel;
244 G4VLowEnergyModel* theIonEffChargeModel;
245 G4VLowEnergyModel* theNuclearStoppingModel;
246 G4VLowEnergyModel* theIonChuFluctuationModel;
247 G4VLowEnergyModel* theIonYangFluctuationModel;
248
249 // std::map<G4int,G4double,std::less<G4int> > totalCrossSectionMap;
250
251 // name of parametrisation table of electron stopping power
252 G4String protonTable;
253 G4String antiprotonTable;
254 G4String theNuclearTable;
255
256 // interval of parametrisation of electron stopping power
257 G4double protonLowEnergy;
258 G4double protonHighEnergy;
259 G4double antiprotonLowEnergy;
260 G4double antiprotonHighEnergy;
261
262 // flag of parametrisation of nucleus stopping power
263 G4bool nStopping;
264 G4bool theBarkas;
265
266 G4DataVector cutForDelta;
267 G4DataVector cutForGamma;
268 G4double minGammaEnergy;
269 G4double minElectronEnergy;
270 G4PhysicsTable* theMeanFreePathTable;
271
272 const G4double paramStepLimit; // parameter limits the step at low energy
273
274 G4double fdEdx; // computed in GetContraints
275 G4double fRangeNow ; //
276 G4double charge; //
277 G4double chargeSquare; //
278 G4double initialMass; // mass to calculate Lambda tables
279 G4double fBarkas;
280
281 G4PixeCrossSectionHandler* pixeCrossSectionHandler;
282 G4AtomicDeexcitation atomicDeexcitation;
283 G4String modelK;
284 G4String modelL;
285 G4String modelM;
286 G4double eMinPixe;
287 G4double eMaxPixe;
288
289 G4bool pixeIsActive;
290
291};
292
293
295 G4double,
296 G4double currentMinimumStep,
297 G4double&)
298{
299 G4double step = GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
300
301 // ---- MGP ---- The following line, taken as is from G4hLowEnergyIonisation,
302 // is meaningless: currentMinimumStep is passed by value,
303 // therefore any local modification to it has no effect
304
305 if ((step > 0.) && (step < currentMinimumStep)) currentMinimumStep = step ;
306
307 return step ;
308}
309
310
312{
313 // ---- MGP ---- Better criterion for applicability to be defined;
314 // now hard-coded particle mass > 0.1 * proton_mass
315
316 return (particle.GetPDGCharge() != 0.0 && particle.GetPDGMass() > CLHEP::proton_mass_c2*0.1);
317}
318
319#endif
320
321
322
323
324
325
326
327
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4double GetPDGCharge() const
Definition: G4Step.hh:78
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
void SetPixeCrossSectionM(const G4String &name)
void SetPixeCrossSectionL(const G4String &name)
void SetLowEnergyForProtonParametrisation(G4double energy)
void SetHighEnergyForAntiProtonParametrisation(G4double energy)
void SetNuclearStoppingPowerModel(const G4String &dedxTable)
G4double ComputeDEDX(const G4ParticleDefinition *aParticle, const G4MaterialCutsCouple *couple, G4double kineticEnergy)
void SetLowEnergyForAntiProtonParametrisation(G4double energy)
G4bool IsApplicable(const G4ParticleDefinition &)
void SetCutForAugerElectrons(G4double cut)
G4VParticleChange * AlongStepDoIt(const G4Track &trackData, const G4Step &stepData)
void SetPixeProjectileMaxEnergy(G4double energy)
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
void SetPixeProjectileMinEnergy(G4double energy)
void SetElectronicStoppingPowerModel(const G4ParticleDefinition *aParticle, const G4String &dedxTable)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step)
void SetHighEnergyForProtonParametrisation(G4double energy)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, enum G4ForceCondition *condition)
void SetPixeCrossSectionK(const G4String &name)
void SetCutForSecondaryPhotons(G4double cut)
void ActivateAugerElectronProduction(G4bool val)
void SetPixe(const G4bool)