Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4BetheBlochModel.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: G4BetheBlochModel
34//
35// Author: Vladimir Ivanchenko on base of Laszlo Urban code
36//
37// Creation date: 03.01.2002
38//
39// Modifications:
40//
41// 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
42// 24-01-03 Make models region aware (V.Ivanchenko)
43// 13-02-03 Add name (V.Ivanchenko)
44// 12-11-03 Fix for GenericIons (V.Ivanchenko)
45// 24-03-05 Add G4EmCorrections (V.Ivanchenko)
46// 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
47// 11-04-04 Move MaxSecondaryEnergy to models (V.Ivanchenko)
48// 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
49// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
50// CorrectionsAlongStep needed for ions(V.Ivanchenko)
51
52//
53// Class Description:
54//
55// Implementation of Bethe-Bloch model of energy loss and
56// delta-electron production by heavy charged particles
57
58// -------------------------------------------------------------------
59//
60
61#ifndef G4BetheBlochModel_h
62#define G4BetheBlochModel_h 1
63
65
66#include "G4VEmModel.hh"
67#include "G4NistManager.hh"
68
69class G4EmCorrections;
71
73{
74
75public:
76
78 const G4String& nam = "BetheBloch");
79
80 virtual ~G4BetheBlochModel();
81
82 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83
86 G4double kineticEnergy,
87 G4double cutEnergy,
88 G4double maxEnergy);
89
92 G4double kineticEnergy,
93 G4double Z, G4double A,
94 G4double cutEnergy,
95 G4double maxEnergy);
96
99 G4double kineticEnergy,
100 G4double cutEnergy,
101 G4double maxEnergy);
102
105 G4double kineticEnergy,
106 G4double cutEnergy);
107
109 const G4Material* mat,
110 G4double kineticEnergy);
111
113 const G4Material* mat,
114 G4double kineticEnergy);
115
116 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
117 const G4DynamicParticle* dp,
118 G4double& eloss,
119 G4double&,
120 G4double length);
121
122 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
124 const G4DynamicParticle*,
125 G4double tmin,
126 G4double maxEnergy);
127
128protected:
129
131 G4double kinEnergy);
132
133 inline G4double GetChargeSquareRatio() const;
134
135 inline void SetChargeSquareRatio(G4double val);
136
137private:
138
139 void SetupParameters();
140
141 inline void SetParticle(const G4ParticleDefinition* p);
142
143 inline void SetGenericIon(const G4ParticleDefinition* p);
144
145 // hide assignment operator
146 G4BetheBlochModel & operator=(const G4BetheBlochModel &right);
148
149 const G4ParticleDefinition* particle;
150 G4ParticleDefinition* theElectron;
151 G4EmCorrections* corr;
152 G4ParticleChangeForLoss* fParticleChange;
153 G4NistManager* nist;
154
155 G4double mass;
156 G4double tlimit;
157 G4double spin;
158 G4double magMoment2;
159 G4double chargeSquare;
160 G4double ratio;
161 G4double formfact;
162 G4double twoln10;
163 G4double bg2lim;
164 G4double taulim;
165 G4double corrFactor;
166 G4bool isIon;
167 G4bool isInitialised;
168};
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171
172inline void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p)
173{
174 if(particle != p) {
175 particle = p;
176 if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > 1.5*CLHEP::eplus)
177 { isIon = true; }
178 SetupParameters();
179 }
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
184inline void G4BetheBlochModel::SetGenericIon(const G4ParticleDefinition* p)
185{
186 if(p && particle != p) {
187 if(p->GetParticleName() == "GenericIon") { isIon = true; }
188 }
189}
190
191//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192
194{
195 return chargeSquare;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199
201{
202 chargeSquare = val;
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
207#endif
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual ~G4BetheBlochModel()
G4double GetChargeSquareRatio() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
void SetChargeSquareRatio(G4double val)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4double GetPDGCharge() const
const G4String & GetParticleName() const