Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 96934 2016-05-18 09:10:41Z 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 
65 
66 #include "G4VEmModel.hh"
67 #include "G4NistManager.hh"
68 
69 class G4EmCorrections;
71 
73 {
74 
75 public:
76 
77  explicit G4BetheBlochModel(const G4ParticleDefinition* p = nullptr,
78  const G4String& nam = "BetheBloch");
79 
80  virtual ~G4BetheBlochModel();
81 
82  virtual void Initialise(const G4ParticleDefinition*,
83  const G4DataVector&) override;
84 
86  const G4MaterialCutsCouple* couple) override;
87 
89  const G4ParticleDefinition*,
90  G4double kineticEnergy,
91  G4double cutEnergy,
92  G4double maxEnergy);
93 
95  const G4ParticleDefinition*,
96  G4double kineticEnergy,
98  G4double cutEnergy,
99  G4double maxEnergy) override;
100 
102  const G4ParticleDefinition*,
103  G4double kineticEnergy,
104  G4double cutEnergy,
105  G4double maxEnergy) override;
106 
108  const G4ParticleDefinition*,
109  G4double kineticEnergy,
110  G4double cutEnergy) override;
111 
113  const G4Material* mat,
114  G4double kineticEnergy) override;
115 
117  const G4Material* mat,
118  G4double kineticEnergy) override;
119 
120  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
121  const G4DynamicParticle* dp,
122  G4double& eloss,
123  G4double&,
124  G4double length) override;
125 
126  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
127  const G4MaterialCutsCouple*,
128  const G4DynamicParticle*,
129  G4double tmin,
130  G4double maxEnergy) override;
131 
132 protected:
133 
135  G4double kinEnergy) override;
136 
137  inline G4double GetChargeSquareRatio() const;
138 
139  inline void SetChargeSquareRatio(G4double val);
140 
141 private:
142 
143  void SetupParameters();
144 
145  inline void SetParticle(const G4ParticleDefinition* p);
146 
147  inline void SetGenericIon(const G4ParticleDefinition* p);
148 
149  // hide assignment operator
150  G4BetheBlochModel & operator=(const G4BetheBlochModel &right) = delete;
151  G4BetheBlochModel(const G4BetheBlochModel&) = delete;
152 
153  const G4ParticleDefinition* particle;
154  G4ParticleDefinition* theElectron;
155  G4EmCorrections* corr;
156  G4ParticleChangeForLoss* fParticleChange;
157  G4NistManager* nist;
158 
159  G4double mass;
160  G4double tlimit;
161  G4double spin;
162  G4double magMoment2;
163  G4double chargeSquare;
164  G4double ratio;
165  G4double formfact;
166  G4double twoln10;
167  G4double corrFactor;
168  G4bool isIon;
169 };
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
173 inline void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p)
174 {
175  if(particle != p) {
176  particle = p;
177  if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > CLHEP::eplus)
178  { isIon = true; }
179  SetupParameters();
180  }
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184 
185 inline void G4BetheBlochModel::SetGenericIon(const G4ParticleDefinition* p)
186 {
187  if(p && particle != p) {
188  if(p->GetParticleName() == "GenericIon") { isIon = true; }
189  }
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
195 {
196  return chargeSquare;
197 }
198 
199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
200 
202 {
203  chargeSquare = val;
204 }
205 
206 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
207 
208 #endif
void SetChargeSquareRatio(G4double val)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
const char * p
Definition: xmltok.h:285
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &eloss, G4double &, G4double length) override
const G4String & GetParticleName() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
double A(double temperature)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
static constexpr double eplus
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual ~G4BetheBlochModel()
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4BetheBlochModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="BetheBloch")
double G4double
Definition: G4Types.hh:76
G4double GetChargeSquareRatio() const
G4double GetPDGCharge() const
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override