Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ionIonisation.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 file
30//
31//
32// File name: G4ionIonisation
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 07.05.2002
37//
38// Modifications:
39//
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
42// 13-02-03 SubCutoff regime is assigned to a region (V.Ivanchenko)
43// 18-04-03 Use IonFluctuations (V.Ivanchenko)
44// 03-08-03 Add effective charge (V.Ivanchenko)
45// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
46// 27-05-04 Set integral to be a default regime (V.Ivanchenko)
47// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
48// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
49// 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
50// 10-05-06 Add a possibility to download user data (V.Ivantchenko)
51// 13-05-06 Add data for light ion stopping in water (V.Ivantchenko)
52// 14-01-07 use SetEmModel() and SetFluctModel() from G4VEnergyLossProcess (mma)
53// 16-05-07 Add data for light ion stopping only for GenericIon (V.Ivantchenko)
54// 07-11-07 Fill non-ionizing energy loss (V.Ivantchenko)
55// 12-09-08 Removed InitialiseMassCharge and CorrectionsAlongStep (VI)
56//
57//
58// -------------------------------------------------------------------
59//
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63#include "G4ionIonisation.hh"
65#include "G4SystemOfUnits.hh"
66#include "G4Electron.hh"
67#include "G4GenericIon.hh"
68#include "G4BraggModel.hh"
69#include "G4BraggIonModel.hh"
70#include "G4BetheBlochModel.hh"
71#include "G4LossTableManager.hh"
72#include "G4WaterStopping.hh"
73#include "G4EmCorrections.hh"
74#include "G4EmParameters.hh"
75#include "G4EmStandUtil.hh"
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78
81 theParticle(nullptr),
82 isInitialised(false),
83 stopDataActive(false)
84{
89 eth = 2*CLHEP::MeV;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97
99{
100 return true;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104
106 const G4Material*,
107 G4double cut)
108{
109 return p->GetPDGMass()*(std::sqrt(1. + 0.5*cut/CLHEP::electron_mass_c2) - 1.0);
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
115 const G4ParticleDefinition* part,
116 const G4ParticleDefinition* bpart)
117{
119
120 if(!isInitialised) {
121 theParticle = part;
122
123 // define base particle
124 const G4ParticleDefinition* theBaseParticle = nullptr;
125 const G4int pdg = part->GetPDGEncoding();
126
127 if(part == bpart) {
128 theBaseParticle = nullptr;
129 } else if(nullptr != bpart) {
130 theBaseParticle = bpart;
131 } else if(part == ion || pdg == 1000020040) {
132 theBaseParticle = nullptr;
133 } else {
134 theBaseParticle = ion;
135 }
136 SetBaseParticle(theBaseParticle);
137
138 // model limit defined for protons
139 eth = 2*CLHEP::MeV*part->GetPDGMass()/CLHEP::proton_mass_c2;
140
142 G4double emin = param->MinKinEnergy();
143 G4double emax = param->MaxKinEnergy();
144
145 // define model of energy loss fluctuations
146 if (nullptr == FluctModel()) {
148 }
149
150 if (nullptr == EmModel(0)) { SetEmModel(new G4BraggIonModel()); }
151 // to compute ranges correctly we have to use low-energy
152 // model even if activation limit is high
153 EmModel(0)->SetLowEnergyLimit(emin);
154
155 // high energy limit may be eth or DBL_MAX
156 G4double emax1 = (EmModel(0)->HighEnergyLimit() < emax) ? eth : emax;
157 EmModel(0)->SetHighEnergyLimit(emax1);
158 AddEmModel(1, EmModel(0), FluctModel());
159
160 // second model is used if the first does not cover energy range
161 if(emax1 < emax) {
162 if (nullptr == EmModel(1)) { SetEmModel(new G4BetheBlochModel()); }
163 EmModel(1)->SetLowEnergyLimit(emax1);
164
165 // for extremely heavy particles upper limit of the model
166 // should be increased
167 emax = std::max(emax, eth*10);
168 EmModel(1)->SetHighEnergyLimit(emax);
169 AddEmModel(2, EmModel(1), FluctModel());
170
171 // Add ion stoping tables for Generic Ion if the default
172 // model is used (with eth ~= 2 MeV)
173 if(part == ion && (EmModel(1)->GetName() == "BetheBloch" ||
174 EmModel(1)->GetName() == "BetheBlochGasIon")) {
175 stopDataActive = true;
176 G4WaterStopping ws(corr, true);
177 corr->SetIonisationModels(EmModel(0), EmModel(1));
178 }
179 }
180 isInitialised = true;
181 }
182 // reinitialisation of corrections for the new run
183 if(part == ion) { corr->InitialiseForNewRun(); }
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
188void G4ionIonisation::StreamProcessInfo(std::ostream& out) const
189{
190 if (stopDataActive && G4GenericIon::GenericIon() == theParticle) {
191 out << " Stopping Power data for "
193 << " ion/material pairs" << G4endl;
194 }
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198
200 const G4String& mname,
201 G4PhysicsVector* dVector)
202{
203 corr->AddStoppingData(Z, A, mname, dVector);
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207
208void G4ionIonisation::ProcessDescription(std::ostream& out) const
209{
210 out << " Ion ionisation";
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fIonisation
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetIonisationModels(G4VEmModel *mod1=nullptr, G4VEmModel *mod2=nullptr)
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
G4int GetNumberOfStoppingVectors() const
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void ProcessDescription(std::ostream &outFile) const override
G4VEmModel * EmModel(std::size_t index=0) const
void SetEmModel(G4VEmModel *, G4int index=0)
void SetBaseParticle(const G4ParticleDefinition *p)
G4VEmFluctuationModel * FluctModel() const
void SetLinearLossLimit(G4double val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
~G4ionIonisation() override
void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) final
void StreamProcessInfo(std::ostream &outFile) const override
G4ionIonisation(const G4String &name="ionIoni")
G4bool IsApplicable(const G4ParticleDefinition &p) final
void AddStoppingData(G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
void ProcessDescription(std::ostream &) const override