Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // $Id: G4MuBremsstrahlungModel.hh 97392 2016-06-02 10:10:32Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4MuBremsstrahlungModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 18.05.2002
38 //
39 // Modifications:
40 //
41 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42 // 27-01-03 Make models region aware (V.Ivanchenko)
43 // 13-02-03 Add name (V.Ivanchenko)
44 // 10-02-04 Add lowestKinEnergy (V.Ivanchenko)
45 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
46 // 13-02-06 Add ComputeCrossSectionPerAtom (mma)
47 // 11-10-07 Add ignoreCut flag (V.Ivanchenko)
48 // 28-02-08 Reorganized protected methods and members (V.Ivanchenko)
49 // 06-03-08 Remove obsolete methods and members (V.Ivanchenko)
50 // 31-05-13 Use element selectors instead of local data (V.Ivanchenko)
51 //
52 
53 //
54 // Class Description:
55 //
56 // Implementation of bremssrahlung by muons
57 
58 // -------------------------------------------------------------------
59 //
60 
61 #ifndef G4MuBremsstrahlungModel_h
62 #define G4MuBremsstrahlungModel_h 1
63 
65 
66 #include "G4VEmModel.hh"
67 #include "G4NistManager.hh"
68 
69 class G4Element;
71 
73 {
74 
75 public:
76 
77  explicit G4MuBremsstrahlungModel(const G4ParticleDefinition* p = nullptr,
78  const G4String& nam = "MuBrem");
79 
80  virtual ~G4MuBremsstrahlungModel();
81 
82  virtual void Initialise(const G4ParticleDefinition*,
83  const G4DataVector&) override;
84 
85  virtual void InitialiseLocal(const G4ParticleDefinition*,
86  G4VEmModel* masterModel) override;
87 
89  const G4MaterialCutsCouple*) override;
90 
92  const G4ParticleDefinition*,
93  G4double kineticEnergy,
95  G4double cutEnergy,
96  G4double maxEnergy) override;
97 
99  const G4ParticleDefinition*,
100  G4double kineticEnergy,
101  G4double cutEnergy) override;
102 
103  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
104  const G4MaterialCutsCouple*,
105  const G4DynamicParticle*,
106  G4double tmin,
107  G4double maxEnergy) override;
108 
109  inline void SetLowestKineticEnergy(G4double e);
110 
111  virtual G4double MinPrimaryEnergy(const G4Material*,
112  const G4ParticleDefinition*,
113  G4double) override;
114 
115 protected:
116 
118 
120  G4double Z,
121  G4double cut);
122 
124  G4double Z,
125  G4double gammaEnergy);
126 
127  inline void SetParticle(const G4ParticleDefinition*);
128 
129 private:
130 
131  // hide assignment operator
133  operator=(const G4MuBremsstrahlungModel &right) = delete;
135 
136 protected:
137 
149 
150 private:
151 
152  G4ParticleDefinition* theGamma;
153  G4ParticleChangeForLoss* fParticleChange;
154 
155  G4double lowestKinEnergy;
156  G4double minThreshold;
157 
158  static const G4double xgi[6];
159  static const G4double wgi[6];
160 
161  static G4double fDN[93];
162 };
163 
164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165 
167 {
168  lowestKinEnergy = e;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172 
173 inline
175 {
176  if(!particle) {
177  particle = p;
178  mass = particle->GetPDGMass();
182  }
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186 
187 #endif
G4double ComputeMicroscopicCrossSection(G4double tkin, G4double Z, G4double cut)
const char * p
Definition: xmltok.h:285
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void SetParticle(const G4ParticleDefinition *)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, G4double Z, G4double gammaEnergy)
static constexpr double electron_mass_c2
double A(double temperature)
const G4ParticleDefinition * particle
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
static constexpr double classic_electr_radius
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut)
G4double GetPDGMass() const
G4MuBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBrem")
double G4double
Definition: G4Types.hh:76
static constexpr double fine_structure_const