Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4hIonisation.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 // $Id: G4hIonisation.cc 105121 2017-07-13 13:38:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4hIonisation
34 //
35 // Author: Laszlo Urban
36 //
37 // Creation date: 30.05.1997
38 //
39 // Modifications:
40 //
41 // corrected by L.Urban on 24/09/97
42 // several bugs corrected by L.Urban on 13/01/98
43 // 07-04-98 remove 'tracking cut' of the ionizing particle, mma
44 // 22-10-98 cleanup L.Urban
45 // 02-02-99 bugs fixed , L.Urban
46 // 29-07-99 correction in BuildLossTable for low energy, L.Urban
47 // 10-02-00 modifications , new e.m. structure, L.Urban
48 // 10-08-00 V.Ivanchenko change BuildLambdaTable, in order to
49 // simulate energy losses of ions; correction to
50 // cross section for particles with spin 1 is inserted as well
51 // 28-05-01 V.Ivanchenko minor changes to provide ANSI -wall compilation
52 // 10-08-01 new methods Store/Retrieve PhysicsTable (mma)
53 // 14-08-01 new function ComputeRestrictedMeandEdx() + 'cleanup' (mma)
54 // 29-08-01 PostStepDoIt: correction for spin 1/2 (instead of 1) (mma)
55 // 17-09-01 migration of Materials to pure STL (mma)
56 // 25-09-01 completion of RetrievePhysicsTable() (mma)
57 // 29-10-01 all static functions no more inlined
58 // 08-11-01 Charge renamed zparticle; added to the dedx
59 // 27-03-02 Bug fix in scaling of lambda table (V.Ivanchenko)
60 // 09-04-02 Update calculation of tables for GenericIons (V.Ivanchenko)
61 // 30-04-02 V.Ivanchenko update to new design
62 // 04-12-02 Add verbose level definition (VI)
63 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
64 // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
65 // 13-02-03 SubCutoff regime is assigned to a region (V.Ivanchenko)
66 // 23-05-03 Define default integral + BohrFluctuations (V.Ivanchenko)
67 // 03-06-03 Fix initialisation problem for STD ionisation (V.Ivanchenko)
68 // 04-08-03 Set integral=false to be default (V.Ivanchenko)
69 // 08-08-03 STD substitute standard (V.Ivanchenko)
70 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
71 // 27-05-04 Set integral to be a default regime (V.Ivanchenko)
72 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
73 // 24-03-05 Optimize internal interfaces (V.Ivantchenko)
74 // 12-08-05 SetStepLimits(0.2, 0.1*mm) (mma)
75 // 10-01-06 SetStepLimits -> SetStepFunction (V.Ivanchenko)
76 // 26-05-06 scale negative particles from pi- and pbar,
77 // positive from pi+ and p (VI)
78 // 14-01-07 use SetEmModel() and SetFluctModel() from G4VEnergyLossProcess (mma)
79 // 12-09-08 Removed CorrectionsAlongStep (VI)
80 // 27-05-10 Added G4ICRU73QOModel for anti-protons (VI)
81 //
82 // -------------------------------------------------------------------
83 //
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
87 #include "G4hIonisation.hh"
88 #include "G4PhysicalConstants.hh"
89 #include "G4SystemOfUnits.hh"
90 #include "G4Electron.hh"
91 #include "G4Proton.hh"
92 #include "G4AntiProton.hh"
93 #include "G4BraggModel.hh"
94 #include "G4BetheBlochModel.hh"
95 #include "G4IonFluctuations.hh"
97 #include "G4BohrFluctuations.hh"
98 #include "G4UnitsTable.hh"
99 #include "G4PionPlus.hh"
100 #include "G4PionMinus.hh"
101 #include "G4KaonPlus.hh"
102 #include "G4KaonMinus.hh"
103 #include "G4ICRU73QOModel.hh"
104 #include "G4EmParameters.hh"
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
108 using namespace std;
109 
111  : G4VEnergyLossProcess(name),
112  isInitialised(false)
113 {
116  mass = 0.0;
117  ratio = 0.0;
118  eth = 2*MeV;
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122 
124 {}
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
129 {
130  return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 10.0*MeV &&
131  !p.IsShortLived());
132 }
133 
134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135 
137  const G4Material*,
138  G4double cut)
139 {
140  G4double x = 0.5*cut/electron_mass_c2;
141  G4double gam = x*ratio + std::sqrt((1. + x)*(1. + x*ratio*ratio));
142  return mass*(gam - 1.0);
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146 
148  const G4ParticleDefinition* part,
149  const G4ParticleDefinition* bpart)
150 {
151  if(!isInitialised) {
152 
153  const G4ParticleDefinition* theBaseParticle = nullptr;
154  G4String pname = part->GetParticleName();
155  G4double q = part->GetPDGCharge();
156 
157  //G4cout << " G4hIonisation::InitialiseEnergyLossProcess " << pname
158  // << " " << bpart << G4endl;
159 
160  // standard base particles
161  if(part == bpart || pname == "proton" ||
162  pname == "anti_proton" ||
163  pname == "pi+" || pname == "pi-" ||
164  pname == "kaon+" || pname == "kaon-" || pname == "GenericIon"
165  || pname == "He3" || pname == "alpha")
166  {
167  theBaseParticle = nullptr;
168  }
169  // select base particle
170  else if(bpart == nullptr) {
171 
172  if(part->GetPDGSpin() == 0.0) {
173  if(q > 0.0) { theBaseParticle = G4KaonPlus::KaonPlus(); }
174  else { theBaseParticle = G4KaonMinus::KaonMinus(); }
175  } else {
176  if(q > 0.0) { theBaseParticle = G4Proton::Proton(); }
177  else { theBaseParticle = G4AntiProton::AntiProton(); }
178  }
179 
180  // base particle defined by interface
181  } else {
182  theBaseParticle = bpart;
183  }
184  SetBaseParticle(theBaseParticle);
185 
186  mass = part->GetPDGMass();
187  ratio = electron_mass_c2/mass;
188  eth = 2.0*MeV*mass/proton_mass_c2;
189 
191  G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
192  G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
193 
194  if(emin != param->MinKinEnergy() || emax != param->MaxKinEnergy()) {
195  SetMinKinEnergy(emin);
196  SetMaxKinEnergy(emax);
197  G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
198  SetDEDXBinning(bin);
199  }
200 
201  if (!EmModel(1)) {
202  if(q > 0.0) { SetEmModel(new G4BraggModel(),1); }
203  else { SetEmModel(new G4ICRU73QOModel(),1); }
204  }
205  EmModel(1)->SetLowEnergyLimit(param->MinKinEnergy());
206  EmModel(1)->SetHighEnergyLimit(eth);
207  AddEmModel(1, EmModel(1), new G4IonFluctuations());
208 
210 
211  if (!EmModel(2)) { SetEmModel(new G4BetheBlochModel(),2); }
212  EmModel(2)->SetLowEnergyLimit(eth);
214  AddEmModel(2, EmModel(2), FluctModel());
215 
216  isInitialised = true;
217  }
218 }
219 
220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
221 
223 {}
224 
225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4hIonisation(const G4String &name="hIoni")
G4int NumberOfBinsPerDecade() const
const XML_Char * name
Definition: expat.h:151
G4double MaxKinEnergy() const
tuple bin
Definition: plottest35.py:22
const char * p
Definition: xmltok.h:285
virtual ~G4hIonisation()
tuple x
Definition: test.py:50
void SetFluctModel(G4VEmFluctuationModel *)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4VEmFluctuationModel * FluctModel()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
virtual void PrintInfo() final
bool G4bool
Definition: G4Types.hh:79
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
float proton_mass_c2
Definition: hepunit.py:275
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) final
void SetMaxKinEnergy(G4double e)
float electron_mass_c2
Definition: hepunit.py:274
G4double MinKinEnergy() const
void SetSecondaryParticle(const G4ParticleDefinition *p)
string pname
Definition: eplot.py:33
static const G4double emax
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4EmParameters * Instance()
void SetEmModel(G4VEmModel *, G4int index=1)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetPDGSpin() const
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
static G4Electron * Electron()
Definition: G4Electron.cc:94
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double GetPDGCharge() const
void SetBaseParticle(const G4ParticleDefinition *p)
G4VEmModel * EmModel(G4int index=1) const
void SetDEDXBinning(G4int nbins)
void SetMinKinEnergy(G4double e)