Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4MuBetheBlochModel.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: G4MuBetheBlochModel.hh 97392 2016-06-02 10:10:32Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4MuBetheBlochModel
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 09.08.2002
38 //
39 // Modifications:
40 //
41 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
42 // 24-01-03 Make models region aware (V.Ivanchenko)
43 // 13-02-03 Add Nama (V.Ivanchenko)
44 // 10-02-04 Calculation of radiative corrections using R.Kokoulin model (V.I)
45 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
46 // 12-04-05 Add usage of G4EmCorrections (V.Ivanchenko)
47 // 13-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
48 //
49 
50 //
51 // Class Description:
52 //
53 // Implementation of Bethe-Bloch model of energy loss and
54 // delta-electron production by heavy charged particles
55 
56 // -------------------------------------------------------------------
57 //
58 
59 #ifndef G4MuBetheBlochModel_h
60 #define G4MuBetheBlochModel_h 1
61 
63 
64 #include "G4VEmModel.hh"
65 
67 class G4EmCorrections;
68 
70 {
71 
72 public:
73 
74  explicit G4MuBetheBlochModel(const G4ParticleDefinition* p = nullptr,
75  const G4String& nam = "MuBetheBloch");
76 
77  virtual ~G4MuBetheBlochModel();
78 
79  virtual void Initialise(const G4ParticleDefinition*,
80  const G4DataVector&) override;
81 
83  const G4MaterialCutsCouple*) override;
84 
86  const G4ParticleDefinition*,
87  G4double kineticEnergy,
88  G4double cutEnergy,
89  G4double maxEnergy);
90 
92  const G4ParticleDefinition*,
93  G4double kineticEnergy,
95  G4double cutEnergy,
96  G4double maxEnergy) override;
97 
99  const G4ParticleDefinition*,
100  G4double kineticEnergy,
101  G4double cutEnergy,
102  G4double maxEnergy) override;
103 
105  const G4ParticleDefinition*,
106  G4double kineticEnergy,
107  G4double cutEnergy) override;
108 
109  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
110  const G4MaterialCutsCouple*,
111  const G4DynamicParticle*,
112  G4double tmin,
113  G4double maxEnergy) override;
114 
115 protected:
116 
118  G4double kinEnergy) override;
119 
120 private:
121 
122  inline void SetParticle(const G4ParticleDefinition* p);
123 
124  // hide assignment operator
125  G4MuBetheBlochModel & operator=(const G4MuBetheBlochModel &right) = delete;
126  G4MuBetheBlochModel(const G4MuBetheBlochModel&) = delete;
127 
128  const G4ParticleDefinition* particle;
129  G4ParticleDefinition* theElectron;
130  G4ParticleChangeForLoss* fParticleChange;
131  G4EmCorrections* corr;
132 
133  G4double limitKinEnergy;
134  G4double logLimitKinEnergy;
135  G4double mass;
136  G4double massSquare;
137  G4double ratio;
138  G4double twoln10;
139  //G4double bg2lim;
140  //G4double taulim;
141  G4double alphaprime;
142  static G4double xgi[8],wgi[8];
143 };
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146 
147 inline void G4MuBetheBlochModel::SetParticle(const G4ParticleDefinition* p)
148 {
149  if(!particle) {
150  particle = p;
151  mass = particle->GetPDGMass();
152  massSquare = mass*mass;
153  ratio = CLHEP::electron_mass_c2/mass;
154  }
155 }
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
158 
159 #endif
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
G4MuBetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MuBetheBloch")
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
const char * p
Definition: xmltok.h:285
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
static constexpr double electron_mass_c2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
double A(double temperature)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double GetPDGMass() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
double G4double
Definition: G4Types.hh:76