Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuBremsstrahlungModel.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: G4MuBremsstrahlungModel
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 18.05.2002
37//
38// Modifications:
39//
40// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
41// 27-01-03 Make models region aware (V.Ivanchenko)
42// 13-02-03 Add name (V.Ivanchenko)
43// 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
44// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
45// 13-02-06 Add ComputeCrossSectionPerAtom (mma)
46// 11-10-07 Add ignoreCut flag (V.Ivanchenko)
47// 28-02-08 Reorganized protected methods and members (V.Ivanchenko)
48// 06-03-08 Remove obsolete methods and members (V.Ivanchenko)
49// 31-05-13 Use element selectors instead of local data (V.Ivanchenko)
50//
51// Class Description:
52//
53// Implementation of bremssrahlung by muons
54// A.G. Bogdanov et al., IEEE Trans. Nuc. Sci., Vol.53, No.2, 2006
55//
56// -------------------------------------------------------------------
57//
58
59#ifndef G4MuBremsstrahlungModel_h
60#define G4MuBremsstrahlungModel_h 1
61
63
64#include "G4VEmModel.hh"
65#include "G4NistManager.hh"
66
67class G4Element;
69
71{
72
73public:
74
75 explicit G4MuBremsstrahlungModel(const G4ParticleDefinition* p = nullptr,
76 const G4String& nam = "MuBrem");
77
79
80 virtual void Initialise(const G4ParticleDefinition*,
81 const G4DataVector&) override;
82
83 virtual void InitialiseLocal(const G4ParticleDefinition*,
84 G4VEmModel* masterModel) override;
85
87 const G4MaterialCutsCouple*) override;
88
91 G4double kineticEnergy,
92 G4double Z, G4double A,
93 G4double cutEnergy,
94 G4double maxEnergy) override;
95
98 G4double kineticEnergy,
99 G4double cutEnergy) override;
100
101 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
103 const G4DynamicParticle*,
104 G4double tmin,
105 G4double maxEnergy) override;
106
107 inline void SetLowestKineticEnergy(G4double e);
108
109 virtual G4double MinPrimaryEnergy(const G4Material*,
110 const G4ParticleDefinition*,
111 G4double) override;
112
113 // hide assignment operator
115 operator=(const G4MuBremsstrahlungModel &right) = delete;
117
118protected:
119
121
123 G4double Z,
124 G4double cut);
125
127 G4double Z,
128 G4double gammaEnergy);
129
130 inline void SetParticle(const G4ParticleDefinition*);
131
132protected:
133
145
148
151
152 static const G4double xgi[6];
153 static const G4double wgi[6];
154
155 static G4double fDN[93];
156};
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159
161{
162 lowestKinEnergy = e;
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166
167inline
169{
170 if(!particle) {
171 particle = p;
173 rmass=mass/CLHEP::electron_mass_c2 ;
174 cc=CLHEP::classic_electr_radius/rmass ;
175 coeff= 16.*CLHEP::fine_structure_const*cc*cc/3. ;
176 }
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
181#endif
double G4double
Definition: G4Types.hh:83
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
~G4MuBremsstrahlungModel()=default
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void SetParticle(const G4ParticleDefinition *)
static const G4double xgi[6]
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static const G4double wgi[6]
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4MuBremsstrahlungModel & operator=(const G4MuBremsstrahlungModel &right)=delete
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4ParticleDefinition * theGamma
G4ParticleChangeForLoss * fParticleChange
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
const G4ParticleDefinition * particle
G4MuBremsstrahlungModel(const G4MuBremsstrahlungModel &)=delete