Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BetheBlochModel.hh
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 header file
30//
31//
32// File name: G4BetheBlochModel
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 03.01.2002
37//
38// Modifications:
39//
40// Modified by Michel Maire and Vladimir Ivanchenko
41//
42// Class Description:
43//
44// Implementation of Bethe-Bloch model of energy loss and
45// delta-electron production by heavy charged particles
46
47// -------------------------------------------------------------------
48//
49
50#ifndef G4BetheBlochModel_h
51#define G4BetheBlochModel_h 1
52
54
55#include "G4VEmModel.hh"
56#include "G4NistManager.hh"
57
58class G4EmCorrections;
61
63{
64
65public:
66
67 explicit G4BetheBlochModel(const G4ParticleDefinition* p = nullptr,
68 const G4String& nam = "BetheBloch");
69
70 virtual ~G4BetheBlochModel();
71
72 virtual void Initialise(const G4ParticleDefinition*,
73 const G4DataVector&) override;
74
76 const G4MaterialCutsCouple* couple) override;
77
80 G4double kineticEnergy,
81 G4double cutEnergy,
82 G4double maxEnergy);
83
86 G4double kineticEnergy,
88 G4double cutEnergy,
89 G4double maxEnergy) override;
90
93 G4double kineticEnergy,
94 G4double cutEnergy,
95 G4double maxEnergy) override;
96
99 G4double kineticEnergy,
100 G4double cutEnergy) override;
101
103 const G4Material* mat,
104 G4double kineticEnergy) override;
105
107 const G4Material* mat,
108 G4double kineticEnergy) override;
109
110 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
111 const G4DynamicParticle* dp,
112 G4double& eloss,
113 G4double&,
114 G4double length) override;
115
116 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
118 const G4DynamicParticle*,
119 G4double tmin,
120 G4double maxEnergy) override;
121
122protected:
123
125 G4double kinEnergy) override;
126
127 inline G4double GetChargeSquareRatio() const;
128
129 inline void SetChargeSquareRatio(G4double val);
130
131private:
132
133 void SetupParameters();
134
135 inline void SetParticle(const G4ParticleDefinition* p);
136
137 inline void SetGenericIon(const G4ParticleDefinition* p);
138
139 // hide assignment operator
140 G4BetheBlochModel & operator=(const G4BetheBlochModel &right) = delete;
141 G4BetheBlochModel(const G4BetheBlochModel&) = delete;
142
143 const G4ParticleDefinition* particle;
144 G4ParticleDefinition* theElectron;
145 G4EmCorrections* corr;
146 G4ParticleChangeForLoss* fParticleChange;
147 G4NistManager* nist;
148 G4ICRU90StoppingData* fICRU90;
149 const G4Material* currentMaterial;
150 const G4Material* baseMaterial;
151
152 G4double mass;
153 G4double tlimit;
154 G4double spin;
155 G4double magMoment2;
156 G4double chargeSquare;
157 G4double ratio;
158 G4double formfact;
159 G4double twoln10;
160 G4double corrFactor;
161 G4double fAlphaTlimit;
162 G4double fProtonTlimit;
163
164 G4int iICRU90;
165 G4bool isIon;
166};
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169
170inline void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p)
171{
172 if(particle != p) {
173 particle = p;
174 if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > CLHEP::eplus)
175 { isIon = true; }
176 SetupParameters();
177 }
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
182inline void G4BetheBlochModel::SetGenericIon(const G4ParticleDefinition* p)
183{
184 if(p && p->GetParticleName() == "GenericIon") { isIon = true; }
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
190{
191 return chargeSquare;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
195
197{
198 chargeSquare = val;
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
202
203#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
virtual ~G4BetheBlochModel()
G4double GetChargeSquareRatio() const
void SetChargeSquareRatio(G4double val)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double GetPDGCharge() const
const G4String & GetParticleName() const