Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4hhIonisation.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: G4hhIonisation
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 30.09.2005
37//
38// Modifications:
39// 10-01-06 SetStepLimits -> SetStepFunction (V.Ivantchenko)
40// 27-10-06 Add maxKinEnergy (V.Ivantchenko)
41//
42//
43// -------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
48#include "G4hhIonisation.hh"
50#include "G4SystemOfUnits.hh"
55#include "G4IonFluctuations.hh"
56#include "G4UnitsTable.hh"
57#include "G4Electron.hh"
58#include "G4EmParameters.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78{
79 return (p.GetPDGCharge() != 0.0);
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
85 const G4Material*,
86 G4double cut)
87{
88 G4double x = 0.5*cut/electron_mass_c2;
89 G4double y = electron_mass_c2/mass;
90 G4double gam = x*y + std::sqrt((1. + x)*(1. + x*y*y));
91 return mass*(gam - 1.0);
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95
97 const G4ParticleDefinition* part,
98 const G4ParticleDefinition* bpart)
99{
100 if(isInitialised) { return; }
101
102 theParticle = part;
103 if(nullptr != bpart) {
104 G4cout << "G4hhIonisation::InitialiseEnergyLossProcess WARNING: no "
105 << "base particle should be defined for the process "
106 << GetProcessName() << G4endl;
107 }
108 mass = theParticle->GetPDGMass();
109 ratio = electron_mass_c2/mass;
110 G4double eth = 2*CLHEP::MeV*mass/proton_mass_c2;
111 flucModel = new G4IonFluctuations();
112
114 G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
115 G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
116
117 SetMinKinEnergy(emin);
118 SetMaxKinEnergy(emax);
119 G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
120 SetDEDXBinning(bin);
121
122 G4VEmModel* em = EmModel(0);
123 if (nullptr == em) {
124 if(part->GetPDGCharge() > 0.0) { em = new G4BraggNoDeltaModel(); }
125 else { em = new G4ICRU73NoDeltaModel(); }
126 }
127 em->SetLowEnergyLimit(emin);
128 em->SetHighEnergyLimit(eth);
129 AddEmModel(1, em, flucModel);
130
131 em = EmModel(1);
132 if(nullptr == em) { em = new G4BetheBlochNoDeltaModel(); }
133 em->SetLowEnergyLimit(eth);
134 em->SetHighEnergyLimit(emax);
135 AddEmModel(1, em, flucModel);
136
137 if(verboseLevel>1) {
138 G4cout << "G4hhIonisation is initialised" << G4endl;
139 }
140 isInitialised = true;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144
145void G4hhIonisation::ProcessDescription(std::ostream& out) const
146{
147 out << "G4hhIonisation: no delta rays" << G4endl;
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fIonisation
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void ProcessDescription(std::ostream &outFile) const override
G4VEmModel * EmModel(std::size_t index=0) const
void SetDEDXBinning(G4int nbins)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetVerboseLevel(G4int value)
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const
void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
void ProcessDescription(std::ostream &) const override
G4hhIonisation(const G4String &name="hhIoni")
G4bool IsApplicable(const G4ParticleDefinition &p) override
~G4hhIonisation() override
G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) override
int G4lrint(double ad)
Definition templates.hh:134