Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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$
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 //
51 
52 //
53 // Class Description:
54 //
55 // Implementation of bremssrahlung by muons
56 
57 // -------------------------------------------------------------------
58 //
59 
60 #ifndef G4MuBremsstrahlungModel_h
61 #define G4MuBremsstrahlungModel_h 1
62 
64 
65 #include "G4VEmModel.hh"
66 #include "G4NistManager.hh"
67 
68 class G4Element;
70 
72 {
73 
74 public:
75 
77  const G4String& nam = "MuBrem");
78 
79  virtual ~G4MuBremsstrahlungModel();
80 
81  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82 
84  const G4MaterialCutsCouple*);
85 
87  const G4ParticleDefinition*,
88  G4double kineticEnergy,
89  G4double Z, G4double A,
90  G4double cutEnergy,
91  G4double maxEnergy);
92 
94  const G4ParticleDefinition*,
95  G4double kineticEnergy,
96  G4double cutEnergy);
97 
98  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
99  const G4MaterialCutsCouple*,
100  const G4DynamicParticle*,
101  G4double tmin,
102  G4double maxEnergy);
103 
104  inline void SetLowestKineticEnergy(G4double e);
105 
106 protected:
107 
109 
111  G4double Z,
112  G4double cut);
113 
115  G4double Z,
116  G4double gammaEnergy);
117 
118  inline void SetParticle(const G4ParticleDefinition*);
119 
120 private:
121 
122  G4DataVector* ComputePartialSumSigma(const G4Material* material,
123  G4double tkin, G4double cut);
124 
125  const G4Element* SelectRandomAtom(const G4MaterialCutsCouple* couple) const;
126 
127  // hide assignment operator
130 
131 protected:
132 
144 
145 private:
146 
147  G4ParticleDefinition* theGamma;
148  G4ParticleChangeForLoss* fParticleChange;
149 
150  G4double lowestKinEnergy;
151  G4double minThreshold;
152 
153  G4double fDN[93];
154 
155  std::vector<G4DataVector*> partialSumSigma;
156 };
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
161 {
162  lowestKinEnergy = e;
163 }
164 
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166 
167 inline
169 {
170  if(!particle) {
171  particle = p;
172  mass = particle->GetPDGMass();
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