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