Geant4-11
G4MuBetheBlochModel.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: G4MuBetheBlochModel
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 09.08.2002
37//
38// Modifications:
39//
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 24-01-03 Make models region aware (V.Ivanchenko)
42// 13-02-03 Add Nama (V.Ivanchenko)
43// 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
44// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
45// 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
46// 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
47//
48
49//
50// Class Description:
51//
52// Implementation of Bethe-Bloch model of energy loss and
53// delta-electron production by heavy charged particles
54
55// -------------------------------------------------------------------
56//
57
58#ifndef G4MuBetheBlochModel_h
59#define G4MuBetheBlochModel_h 1
60
62
63#include "G4VEmModel.hh"
64
66class G4EmCorrections;
67
69{
70
71public:
72
73 explicit G4MuBetheBlochModel(const G4ParticleDefinition* p = nullptr,
74 const G4String& nam = "MuBetheBloch");
75
77
78 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
79
81 const G4MaterialCutsCouple*) override;
82
85 G4double kineticEnergy,
86 G4double cutEnergy,
87 G4double maxEnergy);
88
91 G4double kineticEnergy,
93 G4double cutEnergy,
94 G4double maxEnergy) override;
95
98 G4double kineticEnergy,
99 G4double cutEnergy,
100 G4double maxEnergy) override;
101
104 G4double kineticEnergy,
105 G4double cutEnergy) override;
106
107 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4DynamicParticle*,
110 G4double tmin,
111 G4double maxEnergy) override;
112
113 // hide assignment operator
116
117protected:
118
120 G4double kinEnergy) override;
121
122private:
123
124 inline void SetParticle(const G4ParticleDefinition* p);
125
130
138 static G4double xgi[8],wgi[8];
139};
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142
144{
145 if(nullptr == particle) {
146 particle = p;
150 }
151}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
154
155#endif
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
const G4double A[17]
static G4double wgi[8]
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
static G4double xgi[8]
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
~G4MuBetheBlochModel()=default
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4MuBetheBlochModel(const G4MuBetheBlochModel &)=delete
G4ParticleDefinition * theElectron
void SetParticle(const G4ParticleDefinition *p)
G4EmCorrections * corr
G4MuBetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBetheBloch")
const G4ParticleDefinition * particle
G4MuBetheBlochModel & operator=(const G4MuBetheBlochModel &right)=delete
G4ParticleChangeForLoss * fParticleChange
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
static constexpr double electron_mass_c2