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,v 1.22 2009/02/20 14:48:16 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-cand-01 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class header file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4MuBremsstrahlungModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
37 | // |
---|
38 | // Creation date: 18.05.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) |
---|
43 | // 27-01-03 Make models region aware (V.Ivanchenko) |
---|
44 | // 13-02-03 Add name (V.Ivanchenko) |
---|
45 | // 10-02-04 Add lowestKinEnergy (V.Ivanchenko) |
---|
46 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
47 | // 13-02-06 Add ComputeCrossSectionPerAtom (mma) |
---|
48 | // 11-10-07 Add ignoreCut flag (V.Ivanchenko) |
---|
49 | // 28-02-08 Reorganized protected methods and members (V.Ivanchenko) |
---|
50 | // 06-03-08 Remove obsolete methods and members (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 | |
---|
64 | #include "G4VEmModel.hh" |
---|
65 | #include "G4NistManager.hh" |
---|
66 | |
---|
67 | class G4Element; |
---|
68 | class G4ParticleChangeForLoss; |
---|
69 | |
---|
70 | class G4MuBremsstrahlungModel : public G4VEmModel |
---|
71 | { |
---|
72 | |
---|
73 | public: |
---|
74 | |
---|
75 | G4MuBremsstrahlungModel(const G4ParticleDefinition* p = 0, |
---|
76 | const G4String& nam = "MuBrem"); |
---|
77 | |
---|
78 | virtual ~G4MuBremsstrahlungModel(); |
---|
79 | |
---|
80 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
81 | |
---|
82 | virtual G4double MinEnergyCut(const G4ParticleDefinition*, |
---|
83 | const G4MaterialCutsCouple*); |
---|
84 | |
---|
85 | virtual G4double ComputeCrossSectionPerAtom( |
---|
86 | const G4ParticleDefinition*, |
---|
87 | G4double kineticEnergy, |
---|
88 | G4double Z, G4double A, |
---|
89 | G4double cutEnergy, |
---|
90 | G4double maxEnergy); |
---|
91 | |
---|
92 | virtual G4double ComputeDEDXPerVolume(const G4Material*, |
---|
93 | const G4ParticleDefinition*, |
---|
94 | G4double kineticEnergy, |
---|
95 | G4double cutEnergy); |
---|
96 | |
---|
97 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
98 | const G4MaterialCutsCouple*, |
---|
99 | const G4DynamicParticle*, |
---|
100 | G4double tmin, |
---|
101 | G4double maxEnergy); |
---|
102 | |
---|
103 | inline void SetLowestKineticEnergy(G4double e); |
---|
104 | |
---|
105 | protected: |
---|
106 | |
---|
107 | G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut); |
---|
108 | |
---|
109 | G4double ComputeMicroscopicCrossSection(G4double tkin, |
---|
110 | G4double Z, |
---|
111 | G4double cut); |
---|
112 | |
---|
113 | virtual G4double ComputeDMicroscopicCrossSection(G4double tkin, |
---|
114 | G4double Z, |
---|
115 | G4double gammaEnergy); |
---|
116 | |
---|
117 | inline void SetParticle(const G4ParticleDefinition*); |
---|
118 | |
---|
119 | private: |
---|
120 | |
---|
121 | G4DataVector* ComputePartialSumSigma(const G4Material* material, |
---|
122 | G4double tkin, G4double cut); |
---|
123 | |
---|
124 | const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple) const; |
---|
125 | |
---|
126 | // hide assignment operator |
---|
127 | G4MuBremsstrahlungModel & operator=(const G4MuBremsstrahlungModel &right); |
---|
128 | G4MuBremsstrahlungModel(const G4MuBremsstrahlungModel&); |
---|
129 | |
---|
130 | protected: |
---|
131 | |
---|
132 | const G4ParticleDefinition* particle; |
---|
133 | G4NistManager* nist; |
---|
134 | G4double mass; |
---|
135 | G4double rmass; |
---|
136 | G4double cc; |
---|
137 | G4double coeff; |
---|
138 | G4double sqrte; |
---|
139 | G4double bh; |
---|
140 | G4double bh1; |
---|
141 | G4double btf; |
---|
142 | G4double btf1; |
---|
143 | |
---|
144 | private: |
---|
145 | |
---|
146 | G4ParticleDefinition* theGamma; |
---|
147 | G4ParticleChangeForLoss* fParticleChange; |
---|
148 | |
---|
149 | G4double highKinEnergy; |
---|
150 | G4double lowKinEnergy; |
---|
151 | G4double lowestKinEnergy; |
---|
152 | G4double minThreshold; |
---|
153 | |
---|
154 | std::vector<G4DataVector*> partialSumSigma; |
---|
155 | }; |
---|
156 | |
---|
157 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
158 | |
---|
159 | inline void G4MuBremsstrahlungModel::SetLowestKineticEnergy(G4double e) |
---|
160 | { |
---|
161 | lowestKinEnergy = e; |
---|
162 | } |
---|
163 | |
---|
164 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
165 | |
---|
166 | inline |
---|
167 | void G4MuBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p) |
---|
168 | { |
---|
169 | if(!particle) { |
---|
170 | particle = p; |
---|
171 | mass = particle->GetPDGMass(); |
---|
172 | rmass=mass/electron_mass_c2 ; |
---|
173 | cc=classic_electr_radius/rmass ; |
---|
174 | coeff= 16.*fine_structure_const*cc*cc/3. ; |
---|
175 | } |
---|
176 | } |
---|
177 | |
---|
178 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
179 | |
---|
180 | #endif |
---|