Geant4-11
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
53#include "G4VEmModel.hh"
54
55class G4EmCorrections;
58class G4NistManager;
59
61{
62
63public:
64
65 explicit G4BetheBlochModel(const G4ParticleDefinition* p = nullptr,
66 const G4String& nam = "BetheBloch");
67
68 ~G4BetheBlochModel() override;
69
70 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
71
73 const G4MaterialCutsCouple* couple) override;
74
77 G4double kineticEnergy,
78 G4double cutEnergy,
79 G4double maxEnergy);
80
83 G4double kineticEnergy,
85 G4double cutEnergy,
86 G4double maxEnergy) override;
87
90 G4double kineticEnergy,
91 G4double cutEnergy,
92 G4double maxEnergy) override;
93
96 G4double kineticEnergy,
97 G4double cutEnergy) override;
98
100 const G4Material* mat,
101 G4double kineticEnergy) override;
102
104 const G4Material* mat,
105 G4double kineticEnergy) override;
106
108 const G4DynamicParticle* dp,
109 const G4double& length,
110 G4double& eloss) override;
111
112 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
114 const G4DynamicParticle*,
115 G4double tmin,
116 G4double maxEnergy) override;
117
118 // hide assignment operator
121
122protected:
123
125 G4double kinEnergy) override;
126
127 inline G4double GetChargeSquareRatio() const;
128
129 inline void SetChargeSquareRatio(G4double val);
130
131private:
132
134
141 const G4Material* currentMaterial = nullptr;
142 const G4Material* baseMaterial = nullptr;
143
154
156 G4bool isIon = false;
158};
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
163{
164 return chargeSquare;
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
170{
171 chargeSquare = val;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175
176#endif
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]
void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, const G4double &length, G4double &eloss) override
G4EmCorrections * corr
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
const G4ParticleDefinition * theElectron
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double GetChargeSquareRatio() const
G4NistManager * nist
void SetChargeSquareRatio(G4double val)
~G4BetheBlochModel() override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
const G4ParticleDefinition * particle
G4ParticleChangeForLoss * fParticleChange
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4BetheBlochModel(const G4BetheBlochModel &)=delete
const G4Material * currentMaterial
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4BetheBlochModel & operator=(const G4BetheBlochModel &right)=delete
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4BetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheBloch")
G4ICRU90StoppingData * fICRU90
void SetupParameters(const G4ParticleDefinition *p)
const G4Material * baseMaterial
#define DBL_MAX
Definition: templates.hh:62