Geant4  10.00.p02
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: G4BetheBlochModel.hh 72939 2013-08-14 13:25:36Z gcosmo $
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 
64 #include <CLHEP/Units/SystemOfUnits.h>
65 
66 #include "G4VEmModel.hh"
67 #include "G4NistManager.hh"
68 
69 class G4EmCorrections;
71 
73 {
74 
75 public:
76 
78  const G4String& nam = "BetheBloch");
79 
80  virtual ~G4BetheBlochModel();
81 
82  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
83 
85  const G4MaterialCutsCouple* couple);
86 
88  const G4ParticleDefinition*,
89  G4double kineticEnergy,
90  G4double cutEnergy,
91  G4double maxEnergy);
92 
94  const G4ParticleDefinition*,
95  G4double kineticEnergy,
96  G4double Z, G4double A,
97  G4double cutEnergy,
98  G4double maxEnergy);
99 
101  const G4ParticleDefinition*,
102  G4double kineticEnergy,
103  G4double cutEnergy,
104  G4double maxEnergy);
105 
107  const G4ParticleDefinition*,
108  G4double kineticEnergy,
109  G4double cutEnergy);
110 
112  const G4Material* mat,
113  G4double kineticEnergy);
114 
116  const G4Material* mat,
117  G4double kineticEnergy);
118 
119  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
120  const G4DynamicParticle* dp,
121  G4double& eloss,
122  G4double&,
123  G4double length);
124 
125  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
126  const G4MaterialCutsCouple*,
127  const G4DynamicParticle*,
128  G4double tmin,
129  G4double maxEnergy);
130 
131 protected:
132 
134  G4double kinEnergy);
135 
136  inline G4double GetChargeSquareRatio() const;
137 
138  inline void SetChargeSquareRatio(G4double val);
139 
140 private:
141 
142  void SetupParameters();
143 
144  inline void SetParticle(const G4ParticleDefinition* p);
145 
146  inline void SetGenericIon(const G4ParticleDefinition* p);
147 
148  // hide assignment operator
151 
157 
171 };
172 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
176 {
177  if(particle != p) {
178  particle = p;
179  if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > 1.5*CLHEP::eplus)
180  { isIon = true; }
181  SetupParameters();
182  }
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186 
188 {
189  if(p && particle != p) {
190  if(p->GetParticleName() == "GenericIon") { isIon = true; }
191  }
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
195 
197 {
198  return chargeSquare;
199 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
202 
204 {
205  chargeSquare = val;
206 }
207 
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
209 
210 #endif
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
void SetChargeSquareRatio(G4double val)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
G4BetheBlochModel & operator=(const G4BetheBlochModel &right)
const G4ParticleDefinition * particle
const G4String & GetParticleName() const
void SetGenericIon(const G4ParticleDefinition *p)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple)
G4NistManager * nist
G4EmCorrections * corr
G4BetheBlochModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheBloch")
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const G4double A[nN]
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)
void SetParticle(const G4ParticleDefinition *p)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4ParticleChangeForLoss * fParticleChange
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
G4double GetChargeSquareRatio() const
static const double eplus
Definition: G4SIunits.hh:178
G4double GetPDGCharge() const