Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BraggModel.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4BraggModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 03.01.2002
37//
38// Modifications:
39// 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
40// 24-01-03 Make models region aware (V.Ivanchenko)
41// 13-02-03 Add name (V.Ivanchenko)
42// 12-11-03 Fix for GenericIons (V.Ivanchenko)
43// 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
44// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
45// 25-04-06 Added stopping data from PSTAR (V.Ivanchenko)
46// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
47// CorrectionsAlongStep needed for ions(V.Ivanchenko)
48
49//
50// Class Description:
51//
52// Implementation of energy loss and delta-electron production
53// by heavy slow charged particles using ICRU'49 and NIST evaluated data
54// for protons
55
56// -------------------------------------------------------------------
57//
58
59#ifndef G4BraggModel_h
60#define G4BraggModel_h 1
61
63
64#include "G4VEmModel.hh"
65#include "G4PSTARStopping.hh"
66
68class G4EmCorrections;
70
72{
73
74public:
75
76 explicit G4BraggModel(const G4ParticleDefinition* p = nullptr,
77 const G4String& nam = "Bragg");
78
79 virtual ~G4BraggModel();
80
81 virtual void Initialise(const G4ParticleDefinition*,
82 const G4DataVector&) override;
83
86 G4double kineticEnergy,
87 G4double cutEnergy,
88 G4double maxEnergy);
89
92 G4double kineticEnergy,
94 G4double cutEnergy,
95 G4double maxEnergy) override;
96
99 G4double kineticEnergy,
100 G4double cutEnergy,
101 G4double maxEnergy) override;
102
105 G4double kineticEnergy,
106 G4double cutEnergy) override;
107
108 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
110 const G4DynamicParticle*,
111 G4double tmin,
112 G4double maxEnergy) override;
113
114 // Compute ion charge
116 const G4Material*,
117 G4double kineticEnergy) override;
118
120 const G4Material* mat,
121 G4double kineticEnergy) override;
122
123protected:
124
126 G4double kinEnergy) final;
127
128 inline G4double GetChargeSquareRatio() const;
129
130 inline void SetChargeSquareRatio(G4double val);
131
132private:
133
134 inline void SetParticle(const G4ParticleDefinition* p);
135
136 void HasMaterial(const G4Material* material);
137
138 G4double StoppingPower(const G4Material* material,
139 G4double kineticEnergy);
140
141 G4double ElectronicStoppingPower(G4double z,
142 G4double kineticEnergy) const;
143
144 G4double DEDX(const G4Material* material, G4double kineticEnergy);
145
146 G4bool MolecIsInZiegler1988(const G4Material* material);
147
148 G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
149
150 // hide assignment operator
151 G4BraggModel & operator=(const G4BraggModel &right) = delete;
152 G4BraggModel(const G4BraggModel&) = delete;
153
154 G4EmCorrections* corr;
155
156 const G4ParticleDefinition* particle;
157 G4ParticleDefinition* theElectron;
158 G4ParticleChangeForLoss* fParticleChange;
159 static G4PSTARStopping* fPSTAR;
160 G4ICRU90StoppingData* fICRU90;
161
162 const G4Material* currentMaterial;
163 const G4Material* baseMaterial;
164
165 G4double mass;
166 G4double spin;
167 G4double chargeSquare;
168 G4double massRate;
169 G4double ratio;
170 G4double lowestKinEnergy;
171 G4double protonMassAMU;
172 G4double theZieglerFactor;
173 G4double expStopPower125; // Experimental Stopping power at 125keV
174
175 G4int iMolecula; // index in the molecula's table
176 G4int iPSTAR; // index in PSTAR
177 G4int iICRU90;
178 G4bool isIon;
179};
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183inline void G4BraggModel::SetParticle(const G4ParticleDefinition* p)
184{
185 particle = p;
186 mass = particle->GetPDGMass();
187 spin = particle->GetPDGSpin();
188 G4double q = particle->GetPDGCharge()/CLHEP::eplus;
189 chargeSquare = q*q;
190 massRate = mass/CLHEP::proton_mass_c2;
191 ratio = CLHEP::electron_mass_c2/mass;
192}
193
195{
196 return chargeSquare;
197}
198
200{
201 chargeSquare = val;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
206#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void SetChargeSquareRatio(G4double val)
G4double GetChargeSquareRatio() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual ~G4BraggModel()
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double GetPDGCharge() const