Geant4 11.2.2
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
62#include "G4VEmModel.hh"
63#include "G4PSTARStopping.hh"
64
66class G4EmCorrections;
68class G4PSTARStopping;
69
71{
72
73public:
74
75 explicit G4BraggModel(const G4ParticleDefinition* p = nullptr,
76 const G4String& nam = "Bragg");
77
78 ~G4BraggModel() override;
79
80 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
81
83 const G4MaterialCutsCouple* couple) override;
84
87 G4double kineticEnergy,
88 G4double cutEnergy,
89 G4double maxEnergy);
90
93 G4double kineticEnergy,
95 G4double cutEnergy,
96 G4double maxEnergy) override;
97
100 G4double kineticEnergy,
101 G4double cutEnergy,
102 G4double maxEnergy) override;
103
106 G4double kineticEnergy,
107 G4double cutEnergy) override;
108
109 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
111 const G4DynamicParticle*,
112 G4double tmin,
113 G4double maxEnergy) override;
114
115 // Compute ion charge
117 const G4Material*,
118 G4double kineticEnergy) override;
119
121 const G4Material* mat,
122 G4double kineticEnergy) override;
123
124 // hide assignment operator
125 G4BraggModel & operator=(const G4BraggModel &right) = delete;
126 G4BraggModel(const G4BraggModel&) = delete;
127
128protected:
129
130 void SetParticle(const G4ParticleDefinition* p);
131
133 G4double kinEnergy) final;
134
135 inline void SetChargeSquareRatio(G4double val);
136
137private:
138
139 void HasMaterial(const G4Material* material);
140
141 G4double StoppingPower(const G4Material* material,
142 G4double kineticEnergy);
143
144 G4double ElectronicStoppingPower(G4double z,
145 G4double kineticEnergy) const;
146
147 G4double DEDX(const G4Material* material, G4double kineticEnergy);
148
149 G4bool MolecIsInZiegler1988(const G4Material* material);
150
151 G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
152
153protected:
154
158
159 const G4Material* currentMaterial = nullptr;
160 const G4Material* baseMaterial = nullptr;
161
163
166
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 NIST PSTAR
180
181private:
182
183 G4bool isIon = false;
184 G4bool isFirst = false;
185};
186
188{
189 chargeSquare = val;
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193
194#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4BraggModel(const G4BraggModel &)=delete
G4ParticleChangeForLoss * fParticleChange
G4double protonMassAMU
~G4BraggModel() override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
const G4Material * baseMaterial
G4double massRate
G4double lowestKinEnergy
G4BraggModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Bragg")
G4double theZieglerFactor
void SetParticle(const G4ParticleDefinition *p)
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4EmCorrections * corr
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)
G4ParticleDefinition * theElectron
G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
G4double chargeSquare
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
const G4ParticleDefinition * particle
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
const G4Material * currentMaterial
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double expStopPower125
static G4PSTARStopping * fPSTAR
G4BraggModel & operator=(const G4BraggModel &right)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static G4ICRU90StoppingData * fICRU90