source: trunk/source/processes/electromagnetic/muons/include/G4MuBremsstrahlungModel.hh@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 6.2 KB
RevLine 
[819]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//
[1007]26// $Id: G4MuBremsstrahlungModel.hh,v 1.21 2008/07/22 16:11:34 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
[819]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)
[961]47// 13-02-06 Add ComputeCrossSectionPerAtom (mma)
[819]48// 11-10-07 Add ignoreCut flag (V.Ivanchenko)
[961]49// 28-02-08 Reorganized protected methods and members (V.Ivanchenko)
50// 06-03-08 Remove obsolete methods and members (V.Ivanchenko)
[819]51//
52
53//
54// Class Description:
55//
[961]56// Implementation of bremssrahlung by muons
[819]57
58// -------------------------------------------------------------------
59//
60
61#ifndef G4MuBremsstrahlungModel_h
62#define G4MuBremsstrahlungModel_h 1
63
64#include "G4VEmModel.hh"
[961]65#include "G4NistManager.hh"
[819]66
67class G4Element;
68class G4ParticleChangeForLoss;
69
70class G4MuBremsstrahlungModel : public G4VEmModel
71{
72
73public:
74
75 G4MuBremsstrahlungModel(const G4ParticleDefinition* p = 0,
76 const G4String& nam = "MuBrem");
77
78 virtual ~G4MuBremsstrahlungModel();
79
[1007]80 void SetParticle(const G4ParticleDefinition*);
[819]81
[1007]82 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83
84 G4double MinEnergyCut(const G4ParticleDefinition*,
85 const G4MaterialCutsCouple*);
[819]86
87 virtual G4double ComputeCrossSectionPerAtom(
88 const G4ParticleDefinition*,
89 G4double kineticEnergy,
90 G4double Z, G4double A,
91 G4double cutEnergy,
92 G4double maxEnergy);
[961]93
[819]94 virtual G4double ComputeDEDXPerVolume(const G4Material*,
95 const G4ParticleDefinition*,
96 G4double kineticEnergy,
97 G4double cutEnergy);
98
[1007]99 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
100 const G4MaterialCutsCouple*,
101 const G4DynamicParticle*,
102 G4double tmin,
103 G4double maxEnergy);
[819]104
[961]105 inline void SetLowestKineticEnergy(G4double e);
106
[819]107protected:
108
[961]109 G4double ComputMuBremLoss(G4double Z, G4double tkin, G4double cut);
110
111 G4double ComputeMicroscopicCrossSection(G4double tkin,
112 G4double Z,
113 G4double cut);
[819]114
[961]115 virtual G4double ComputeDMicroscopicCrossSection(G4double tkin,
116 G4double Z,
117 G4double gammaEnergy);
[819]118
[1007]119 G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
120 G4double kineticEnergy);
[819]121
122private:
123
[961]124 G4DataVector* ComputePartialSumSigma(const G4Material* material,
125 G4double tkin, G4double cut);
[819]126
[961]127 const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple) const;
[819]128
129 // hide assignment operator
130 G4MuBremsstrahlungModel & operator=(const G4MuBremsstrahlungModel &right);
131 G4MuBremsstrahlungModel(const G4MuBremsstrahlungModel&);
132
[961]133protected:
134
135 const G4ParticleDefinition* particle;
136 G4NistManager* nist;
137 G4double mass;
138 G4double rmass;
139 G4double cc;
140 G4double coeff;
141 G4double sqrte;
142 G4double bh;
143 G4double bh1;
144 G4double btf;
145 G4double btf1;
146
147private:
148
[819]149 G4ParticleDefinition* theGamma;
150 G4ParticleChangeForLoss* fParticleChange;
151
152 G4double highKinEnergy;
153 G4double lowKinEnergy;
154 G4double lowestKinEnergy;
155 G4double minThreshold;
156
157 std::vector<G4DataVector*> partialSumSigma;
158};
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
[1007]162inline
163G4double G4MuBremsstrahlungModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
164 G4double kineticEnergy)
165{
166 return kineticEnergy;
167}
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
[961]171inline void G4MuBremsstrahlungModel::SetLowestKineticEnergy(G4double e)
[819]172{
[961]173 lowestKinEnergy = e;
[819]174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
[961]178inline
[1007]179G4double G4MuBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*,
180 const G4MaterialCutsCouple*)
181{
182 return minThreshold;
183}
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186
187inline
[961]188void G4MuBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p)
[819]189{
[961]190 if(!particle) {
191 particle = p;
192 mass = particle->GetPDGMass();
193 rmass=mass/electron_mass_c2 ;
194 cc=classic_electr_radius/rmass ;
195 coeff= 16.*fine_structure_const*cc*cc/3. ;
196 }
[819]197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200
201#endif
Note: See TracBrowser for help on using the repository browser.