Geant4 10.7.0
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 "G4BohrFluctuations.hh"
56#include "G4IonFluctuations.hh"
57#include "G4UnitsTable.hh"
58#include "G4Electron.hh"
59#include "G4EmParameters.hh"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
65 theParticle(nullptr),
66 //theBaseParticle(nullptr),
67 isInitialised(false)
68{
69 SetStepFunction(0.1, 0.1*mm);
73 mass = 0.0;
74 ratio = 0.0;
75 flucModel = nullptr;
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
81{}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86{
87 return (p.GetPDGCharge() != 0.0 && p.GetPDGMass() > 100.0*MeV &&
88 !p.IsShortLived());
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
94 const G4Material*,
95 G4double cut)
96{
97 G4double x = 0.5*cut/electron_mass_c2;
98 G4double y = electron_mass_c2/mass;
99 G4double gam = x*y + std::sqrt((1. + x)*(1. + x*y*y));
100 return mass*(gam - 1.0);
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104
106 const G4ParticleDefinition* part,
107 const G4ParticleDefinition* bpart)
108{
109 if(isInitialised) { return; }
110
111 theParticle = part;
112 if(bpart) {
113 G4cout << "G4hhIonisation::InitialiseEnergyLossProcess WARNING: no "
114 << "base particle should be defined for the process "
115 << GetProcessName() << G4endl;
116 }
118 mass = theParticle->GetPDGMass();
119 ratio = electron_mass_c2/mass;
120 G4double eth = 2*MeV*mass/proton_mass_c2;
121 flucModel = new G4IonFluctuations();
122
124 G4double emin = std::min(param->MinKinEnergy(), 0.1*eth);
125 G4double emax = std::max(param->MaxKinEnergy(), 100*eth);
126
127 SetMinKinEnergy(emin);
128 SetMaxKinEnergy(emax);
129 G4int bin = G4lrint(param->NumberOfBinsPerDecade()*std::log10(emax/emin));
130 SetDEDXBinning(bin);
131
132 G4VEmModel* em = nullptr;
133 if(part->GetPDGCharge() > 0.0) { em = new G4BraggNoDeltaModel(); }
134 else { em = new G4ICRU73NoDeltaModel(); }
135 em->SetLowEnergyLimit(emin);
136 em->SetHighEnergyLimit(eth);
137 AddEmModel(1, em, flucModel);
138
139 em = new G4BetheBlochNoDeltaModel();
140 em->SetLowEnergyLimit(eth);
141 em->SetHighEnergyLimit(emax);
142 SetEmModel(em);
143 AddEmModel(1, em, flucModel);
144
145 if(verboseLevel>1) {
146 G4cout << "G4hhIonisation is initialised" << G4endl;
147 }
148 isInitialised = true;
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152
154{
155 G4cout << " Delta-ray will not be produced; "
156 << G4endl;
157}
158
159//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
160
161void G4hhIonisation::ProcessDescription(std::ostream& out) const
162{
163 out << "No description available." << G4endl;
165}
166
167//....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:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
G4double GetPDGCharge() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetMaxKinEnergy(G4double e)
virtual void ProcessDescription(std::ostream &outFile) const override
void SetDEDXBinning(G4int nbins)
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetBaseParticle(const G4ParticleDefinition *p)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetMinKinEnergy(G4double e)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
virtual void ProcessDescription(std::ostream &) const override
G4hhIonisation(const G4String &name="hhIoni")
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override
virtual void PrintInfo() override
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *p, const G4Material *, G4double cut) override
virtual ~G4hhIonisation()
int G4lrint(double ad)
Definition: templates.hh:134