Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 ~G4BraggModel() override;
80
81 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
82
84 const G4MaterialCutsCouple* couple) override;
85
88 G4double kineticEnergy,
89 G4double cutEnergy,
90 G4double maxEnergy);
91
94 G4double kineticEnergy,
96 G4double cutEnergy,
97 G4double maxEnergy) override;
98
101 G4double kineticEnergy,
102 G4double cutEnergy,
103 G4double maxEnergy) override;
104
107 G4double kineticEnergy,
108 G4double cutEnergy) override;
109
110 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
112 const G4DynamicParticle*,
113 G4double tmin,
114 G4double maxEnergy) override;
115
116 // Compute ion charge
118 const G4Material*,
119 G4double kineticEnergy) override;
120
122 const G4Material* mat,
123 G4double kineticEnergy) override;
124
125 // hide assignment operator
126 G4BraggModel & operator=(const G4BraggModel &right) = delete;
127 G4BraggModel(const G4BraggModel&) = delete;
128
129protected:
130
132 G4double kinEnergy) final;
133
134 inline G4double GetChargeSquareRatio() const;
135
136 inline void SetChargeSquareRatio(G4double val);
137
138private:
139
140 inline void SetParticle(const G4ParticleDefinition* p);
141
142 void HasMaterial(const G4Material* material);
143
144 G4double StoppingPower(const G4Material* material,
145 G4double kineticEnergy);
146
147 G4double ElectronicStoppingPower(G4double z,
148 G4double kineticEnergy) const;
149
150 G4double DEDX(const G4Material* material, G4double kineticEnergy);
151
152 G4bool MolecIsInZiegler1988(const G4Material* material);
153
154 G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
155
156 G4EmCorrections* corr = nullptr;
157 const G4ParticleDefinition* particle = nullptr;
158 G4ParticleDefinition* theElectron = nullptr;
159 G4ParticleChangeForLoss* fParticleChange = nullptr;
160 G4ICRU90StoppingData* fICRU90 = nullptr;
161
162 const G4Material* currentMaterial = nullptr;
163 const G4Material* baseMaterial = nullptr;
164
165 static G4PSTARStopping* fPSTAR;
166
167 G4double mass;
168 G4double spin;
169 G4double chargeSquare;
170 G4double massRate;
171 G4double ratio;
172 G4double lowestKinEnergy;
173 G4double protonMassAMU;
174 G4double theZieglerFactor;
175 G4double expStopPower125; // Experimental Stopping power at 125keV
176
177 G4int iMolecula = -1; // index in the molecula's table
178 G4int iPSTAR = -1; // index in PSTAR
179 G4int iICRU90 = -1;
180 G4bool isIon = false;
181};
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184
185inline void G4BraggModel::SetParticle(const G4ParticleDefinition* p)
186{
187 particle = p;
188 mass = particle->GetPDGMass();
189 spin = particle->GetPDGSpin();
190 G4double q = particle->GetPDGCharge()/CLHEP::eplus;
191 chargeSquare = q*q;
192 massRate = mass/CLHEP::proton_mass_c2;
193 ratio = CLHEP::electron_mass_c2/mass;
194}
195
197{
198 return chargeSquare;
199}
200
202{
203 chargeSquare = val;
204}
205
206//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
207
208#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4BraggModel(const G4BraggModel &)=delete
~G4BraggModel() override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
void SetChargeSquareRatio(G4double val)
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
G4double GetChargeSquareRatio() const
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4BraggModel & operator=(const G4BraggModel &right)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double GetPDGCharge() const