Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BraggIonModel.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: G4BraggIonModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 13.10.2004
37//
38// Modifications:
39// 09-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
40// 11-05-05 Major optimisation of internal interfaces (V.Ivantchenko)
41// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
42// 25-04-06 Add stopping data from ASTAR (V.Ivanchenko)
43// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
44// CorrectionsAlongStep needed for ions(V.Ivanchenko)
45
46//
47// Class Description:
48//
49// Implementation of energy loss and delta-electron production
50// by heavy slow charged particles using ICRU'49 and NIST evaluated data
51// for He4 ions
52
53// -------------------------------------------------------------------
54//
55
56#ifndef G4BraggIonModel_h
57#define G4BraggIonModel_h 1
58
60
61#include "G4VEmModel.hh"
62#include "G4ASTARStopping.hh"
63
65class G4EmCorrections;
67
69{
70
71public:
72
73 explicit G4BraggIonModel(const G4ParticleDefinition* p = nullptr,
74 const G4String& nam = "BraggIon");
75
76 virtual ~G4BraggIonModel();
77
78 virtual void Initialise(const G4ParticleDefinition*,
79 const G4DataVector&) override;
80
82 const G4MaterialCutsCouple* couple) 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
123 // add correction to energy loss and ompute non-ionizing energy loss
124 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
125 const G4DynamicParticle*,
126 G4double& eloss,
127 G4double& niel,
128 G4double length) override;
129
130protected:
131
133 G4double kinEnergy) final;
134
135private:
136
137 void SetParticle(const G4ParticleDefinition* p);
138
139 G4double HeEffChargeSquare(G4double z, G4double kinEnergyInMeV) const;
140
141 // hide assignment operator
142 G4BraggIonModel & operator=(const G4BraggIonModel &right) = delete;
143 G4BraggIonModel(const G4BraggIonModel&) = delete;
144
145 void HasMaterial(const G4Material* material);
146
147 G4double StoppingPower(const G4Material* material,
148 G4double kineticEnergy);
149
150 G4double ElectronicStoppingPower(G4double z,
151 G4double kineticEnergy) const;
152
153 G4double DEDX(const G4Material* material, G4double kineticEnergy);
154
155 G4EmCorrections* corr;
156
157 const G4ParticleDefinition* particle;
158 G4ParticleDefinition* theElectron;
159 G4ParticleChangeForLoss* fParticleChange;
160
161 static G4ASTARStopping* fASTAR;
162 G4ICRU90StoppingData* fICRU90;
163
164 const G4Material* currentMaterial;
165 const G4Material* baseMaterial;
166
167 G4double mass;
168 G4double spin;
169 G4double chargeSquare;
170 G4double massRate;
171 G4double ratio;
172 G4double lowestKinEnergy;
173 G4double HeMass;
174 G4double massFactor;
175 G4double corrFactor;
176 G4double rateMassHe2p;
177 G4double theZieglerFactor;
178
179 G4int iMolecula; // index in the molecula's table
180 G4int iASTAR; // index in ASTAR
181 G4int iICRU90;
182 G4bool isIon;
183};
184
185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
186
187inline void G4BraggIonModel::SetParticle(const G4ParticleDefinition* p)
188{
189 particle = p;
190 mass = particle->GetPDGMass();
191 spin = particle->GetPDGSpin();
192 G4double q = particle->GetPDGCharge()/CLHEP::eplus;
193 chargeSquare = q*q;
194 massRate = mass/CLHEP::proton_mass_c2;
195 ratio = CLHEP::electron_mass_c2/mass;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199
200#endif
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) override
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, 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
virtual ~G4BraggIonModel()
G4double GetPDGCharge() const