source: trunk/source/processes/electromagnetic/lowenergy/include/G4hLowEnergyIonisation.hh @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 13.9 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//      GEANT 4 class header file
29//
30//      History: based on object model of
31//      2nd December 1995, G.Cosmo
32//      ---------- G4hLowEnergyIonisation physics process -----
33//                by Vladimir Ivanchenko, 14 July 1999
34//                was made on the base of G4hIonisation class
35//                developed by Laszlo Urban
36// ************************************************************
37
38// ************************************************************
39// 28 July   1999 V.Ivanchenko cleen up
40// 17 August 1999 G.Mancinelli implemented ICRU parametrization (protons) 
41// 20 August 1999 G.Mancinelli implemented ICRU parametrization (alpha) 
42// 31 August 1999 V.Ivanchenko update and cleen up
43// 23 May    2000    MG Pia  Clean up for QAO model
44// 25 July   2000 V.Ivanchenko New design iteration
45// 09 August 2000 V.Ivanchenko Add GetContinuousStepLimit
46// 17 August 2000 V.Ivanchenko Add IonFluctuationModel
47// 23 Oct    2000 V.Ivanchenko Renew comments
48// 30 Oct    2001 V.Ivanchenko Add minGammaEnergy and minElectronEnergy
49// 07 Dec    2001 V.Ivanchenko Add SetFluorescence method
50// 26 Feb    2002 V.Ivanchenko Add initialMass for GenericIons
51// 21 Jan    2003 V.Ivanchenko Cut per region
52// 03 Oct    2009 ALF added SelectShellIonisationCS
53// ------------------------------------------------------------
54 
55// Class Description:
56// Ionisation process of charged hadrons and ions, including low energy
57// extensions
58// The physics model is described in CERN-OPEN-99-121 and CERN-OPEN-99-300.
59// The user may select parametrisation tables for electronic
60// stopping powers and nuclear stopping powers
61// The list of available tables:
62// Electronic stopping powers: "ICRU_49p" (default), "ICRU_49He",
63//                             "Ziegler1977p", "Ziegler1985p",
64//                             "Ziegler1977He"
65// Nuclear stopping powers:    "ICRU_49" (default), "Ziegler1977",
66//                             "Ziegler1985"
67// Further documentation available from http://www.ge.infn.it/geant4/lowE
68// and in the Physics Reference Manual
69
70// ------------------------------------------------------------
71
72#ifndef G4hLowEnergyIonisation_h
73#define G4hLowEnergyIonisation_h 1
74 
75#include "globals.hh"
76#include "G4hLowEnergyLoss.hh"
77#include "G4VLowEnergyModel.hh"
78#include "G4Track.hh"
79#include "G4Step.hh"
80#include "G4Electron.hh"
81#include "G4PhysicsLogVector.hh"
82#include "G4PhysicsLinearVector.hh"
83#include "G4hNuclearStoppingModel.hh"
84#include "G4hBetheBlochModel.hh"
85#include "G4hParametrisedLossModel.hh"
86#include "G4QAOLowEnergyLoss.hh"
87#include "G4hIonEffChargeSquare.hh"
88#include "G4IonChuFluctuationModel.hh"
89#include "G4IonYangFluctuationModel.hh"
90#include "G4AtomicDeexcitation.hh"
91#include "G4MaterialCutsCouple.hh"
92#include <map>
93
94class G4VEMDataSet;
95class G4ShellVacancy;
96class G4VhShellCrossSection;
97
98class G4hLowEnergyIonisation : public G4hLowEnergyLoss
99{
100public: // With description
101 
102  G4hLowEnergyIonisation(const G4String& processName = "hLowEIoni"); 
103  // The ionisation process for hadrons/ions to be include in the
104  // UserPhysicsList
105
106  ~G4hLowEnergyIonisation();
107  // Destructor
108 
109  G4bool IsApplicable(const G4ParticleDefinition&);
110  // True for all charged hadrons/ions
111   
112  void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ;
113  // Build physics table during initialisation
114
115  G4double GetMeanFreePath(const G4Track& track,
116                                 G4double previousStepSize,
117                            enum G4ForceCondition* condition );
118  // Return MeanFreePath until delta-electron production
119 
120  void PrintInfoDefinition() const;
121  // Print out of the class parameters
122
123  void SetHighEnergyForProtonParametrisation(G4double energy) 
124                             {protonHighEnergy = energy;} ;
125  // Definition of the boundary proton energy. For higher energies
126  // Bethe-Bloch formula is used, for lower energies a parametrisation
127  // of the energy losses is performed. Default is 2 MeV.
128
129  void SetLowEnergyForProtonParametrisation(G4double energy) 
130                             {protonLowEnergy = energy;} ;
131  // Set of the boundary proton energy. For lower energies
132  // the Free Electron Gas model is used for the energy losses.
133  // Default is 1 keV.
134
135  void SetHighEnergyForAntiProtonParametrisation(G4double energy) 
136                             {antiProtonHighEnergy = energy;} ;
137  // Set of the boundary antiproton energy. For higher energies
138  // Bethe-Bloch formula is used, for lower energies parametrisation
139  // of the energy losses is performed. Default is 2 MeV.
140
141  void SetLowEnergyForAntiProtonParametrisation(G4double energy) 
142                              {antiProtonLowEnergy = energy;} ;
143  // Set of the boundary antiproton energy. For lower energies
144  // the Free Electron Gas model is used for the energy losses.
145  // Default is 1 keV.
146
147  G4double GetContinuousStepLimit(const G4Track& track,
148                                        G4double previousStepSize,
149                                        G4double currentMinimumStep,
150                                        G4double& currentSafety); 
151  // Calculation of the step limit due to ionisation losses
152
153  void SetElectronicStoppingPowerModel(const G4ParticleDefinition* aParticle,
154                                       const G4String& dedxTable);
155  // This method defines the electron ionisation parametrisation method
156  // via the name of the table. Default is "ICRU_49p".
157
158  void SetNuclearStoppingPowerModel(const G4String& dedxTable)
159                 {theNuclearTable = dedxTable; SetNuclearStoppingOn();};
160  // This method defines the nuclear ionisation parametrisation method
161  // via the name of the table. Default is "ICRU_49".
162
163  void SetNuclearStoppingOn() {nStopping = true;};
164  // This method switch on calculation of the nuclear stopping power.
165 
166  void SetNuclearStoppingOff() {nStopping = false;};
167  // This method switch off calculation of the nuclear stopping power.
168
169  void SetBarkasOn() {theBarkas = true;};
170  // This method switch on calculation of the Barkas and Bloch effects.
171
172  void SetBarkasOff() {theBarkas = false;};
173  // This method switch off calculation of the Barkas and Bloch effects.
174
175  void SetFluorescence(const G4bool val) {theFluo = val;};
176  // This method switch on/off simulation of the fluorescence of the media.
177
178  void SelectShellIonisationCS(G4String);
179
180
181  G4VParticleChange* AlongStepDoIt(const G4Track& trackData ,
182                                   const G4Step& stepData ) ;
183  // Function to determine total energy deposition on the step
184
185  G4VParticleChange* PostStepDoIt(const G4Track& track,
186                                  const G4Step& Step  ) ;
187  // Simulation of delta rays production.
188
189  G4double ComputeDEDX(const G4ParticleDefinition* aParticle,
190                       const G4MaterialCutsCouple* couple,
191                             G4double kineticEnergy);
192  // This method returns electronic dE/dx for protons or antiproton.
193
194  void SetCutForSecondaryPhotons(G4double cut);
195  // Set threshold energy for fluorescence
196
197  void SetCutForAugerElectrons(G4double cut);
198  // Set threshold energy for Auger electron production
199
200  void ActivateAugerElectronProduction(G4bool val);
201  // Set Auger electron production flag on/off
202
203
204protected:
205
206private:
207
208  void InitializeMe();
209
210  void InitializeParametrisation();
211
212  void BuildLossTable(const G4ParticleDefinition& aParticleType);
213
214  void BuildDataForFluorescence(const G4ParticleDefinition& aParticleType);
215
216  void BuildLambdaTable(const G4ParticleDefinition& aParticleType);
217
218  void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable)
219                              {theProtonTable = dedxTable ;};
220  // This method defines the ionisation parametrisation method via its name
221
222  void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable)
223                              {theAntiProtonTable = dedxTable ;};
224
225  G4double ComputeMicroscopicCrossSection(
226                  const G4ParticleDefinition& aParticleType,
227                        G4double kineticEnergy,
228                        G4double atomicNumber,
229                        G4double deltaCutInEnergy) const;
230
231  G4double GetConstraints(const G4DynamicParticle* particle,
232                          const G4MaterialCutsCouple* couple);
233  // Function to determine StepLimit
234
235  G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
236                                        G4double kineticEnergy) const;
237
238  G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
239                                            G4double kineticEnergy) const;
240
241  G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple,
242                                 G4double kineticEnergy,
243                                 G4double particleMass) const;
244  // This method returns average energy loss due to delta-rays emission with
245  // energy higher than the cut energy for given material.
246
247  G4double BarkasTerm(const G4Material* material,
248                            G4double kineticEnergy) const;
249  // Function to compute the Barkas term for protons
250
251  G4double BlochTerm(const G4Material* material,
252                           G4double kineticEnergy,
253                           G4double cSquare) const;
254  // Function to compute the Bloch term for protons
255
256  G4double ElectronicLossFluctuation(const G4DynamicParticle* particle,
257                                     const G4MaterialCutsCouple* material,
258                                           G4double meanLoss,
259                                           G4double step) const;
260  // Function to sample electronic losses
261
262  std::vector<G4DynamicParticle*>* DeexciteAtom(const G4MaterialCutsCouple* couple,
263                                                  G4double incidentEnergy,
264                                                  G4double hMass,
265                                                  G4double eLoss);
266
267
268  G4int SelectRandomAtom(const G4MaterialCutsCouple* couple,
269                               G4double kineticEnergy) const;
270
271  // hide assignment operator
272  G4hLowEnergyIonisation & operator=(const G4hLowEnergyIonisation &right);
273  G4hLowEnergyIonisation(const G4hLowEnergyIonisation&);
274
275private:
276  //  private data members ...............................
277  G4VLowEnergyModel* theBetheBlochModel;
278  G4VLowEnergyModel* theProtonModel;
279  G4VLowEnergyModel* theAntiProtonModel;
280  G4VLowEnergyModel* theIonEffChargeModel;
281  G4VLowEnergyModel* theNuclearStoppingModel;
282  G4VLowEnergyModel* theIonChuFluctuationModel;
283  G4VLowEnergyModel* theIonYangFluctuationModel;
284  std::map<G4int,G4double,std::less<G4int> > totalCrossSectionMap;
285
286  // name of parametrisation table of electron stopping power
287  G4String theProtonTable;
288  G4String theAntiProtonTable;
289  G4String theNuclearTable;
290
291  // interval of parametrisation of electron stopping power
292  G4double protonLowEnergy;
293  G4double protonHighEnergy;
294  G4double antiProtonLowEnergy;
295  G4double antiProtonHighEnergy;
296
297  // flag of parametrisation of nucleus stopping power
298  G4bool nStopping;
299  G4bool theBarkas;
300
301  G4DataVector cutForDelta;
302  G4DataVector cutForGamma;
303  G4double minGammaEnergy;
304  G4double minElectronEnergy;
305  G4PhysicsTable* theMeanFreePathTable;
306
307  const G4double paramStepLimit; // parameter limits the step at low energy
308
309  G4double fdEdx;        // computed in GetContraints
310  G4double fRangeNow ;   //
311  G4double charge;       //
312  G4double chargeSquare; //
313  G4double initialMass;  // mass to calculate Lambda tables
314  G4double fBarkas;
315
316  G4AtomicDeexcitation deexcitationManager;
317  G4ShellVacancy* shellVacancy;
318  G4VhShellCrossSection* shellCS;
319  std::vector<G4VEMDataSet*> zFluoDataVector;
320  G4bool theFluo;
321};
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
324
325inline G4double G4hLowEnergyIonisation::GetContinuousStepLimit(
326                                        const G4Track& track,
327                                              G4double,
328                                              G4double currentMinimumStep,
329                                              G4double&)
330{
331  G4double Step =
332    GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
333
334  if((Step>0.0)&&(Step<currentMinimumStep))
335     currentMinimumStep = Step ;
336
337  return Step ;
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
341
342inline G4bool G4hLowEnergyIonisation::IsApplicable(
343                                const G4ParticleDefinition& particle)
344{
345   return(particle.GetPDGCharge() != 0.0
346       && particle.GetPDGMass() > proton_mass_c2*0.1);
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350
351#endif
352
353
354
355
356
357
358
359
Note: See TracBrowser for help on using the repository browser.