Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4BraggIonModel
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 13.10.2004
38//
39// Modifications:
40// 09-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
41// 11-05-05 Major optimisation of internal interfaces (V.Ivantchenko)
42// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
43// 25-04-06 Add stopping data from ASTAR (V.Ivanchenko)
44// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
45// CorrectionsAlongStep needed for ions(V.Ivanchenko)
46
47//
48// Class Description:
49//
50// Implementation of energy loss and delta-electron production
51// by heavy slow charged particles using ICRU'49 and NIST evaluated data
52// for He4 ions
53
54// -------------------------------------------------------------------
55//
56
57#ifndef G4BraggIonModel_h
58#define G4BraggIonModel_h 1
59
61
62#include "G4VEmModel.hh"
63#include "G4ASTARStopping.hh"
64
66class G4EmCorrections;
67
69{
70
71public:
72
74 const G4String& nam = "BraggIon");
75
76 virtual ~G4BraggIonModel();
77
78 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
79
82 G4double kineticEnergy,
83 G4double cutEnergy,
84 G4double maxEnergy);
85
88 G4double kineticEnergy,
89 G4double Z, G4double A,
90 G4double cutEnergy,
91 G4double maxEnergy);
92
95 G4double kineticEnergy,
96 G4double cutEnergy,
97 G4double maxEnergy);
98
101 G4double kineticEnergy,
102 G4double cutEnergy);
103
104 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
106 const G4DynamicParticle*,
107 G4double tmin,
108 G4double maxEnergy);
109
110 // Compute ion charge
112 const G4Material*,
113 G4double kineticEnergy);
114
116 const G4Material* mat,
117 G4double kineticEnergy);
118
119 // add correction to energy loss and ompute non-ionizing energy loss
120 virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
121 const G4DynamicParticle*,
122 G4double& eloss,
123 G4double& niel,
124 G4double length);
125
126protected:
127
129 G4double kinEnergy);
130
131private:
132
133 void SetParticle(const G4ParticleDefinition* p);
134
135 G4double HeEffChargeSquare(G4double z, G4double kinEnergyInMeV) const;
136
137 // hide assignment operator
138 G4BraggIonModel & operator=(const G4BraggIonModel &right);
140
141 G4bool HasMaterial(const G4Material* material);
142
143 G4double StoppingPower(const G4Material* material,
144 G4double kineticEnergy);
145
146 G4double ElectronicStoppingPower(G4double z,
147 G4double kineticEnergy) const;
148
149 G4double DEDX(const G4Material* material, G4double kineticEnergy);
150
151 G4EmCorrections* corr;
152
153 const G4ParticleDefinition* particle;
154 G4ParticleDefinition* theElectron;
155 G4ParticleChangeForLoss* fParticleChange;
156
157 G4ASTARStopping astar;
158
159 const G4Material* currentMaterial;
160
161 G4double mass;
162 G4double spin;
163 G4double chargeSquare;
164 G4double massRate;
165 G4double ratio;
166 G4double lowestKinEnergy;
167 G4double HeMass;
168 G4double massFactor;
169 G4double corrFactor;
170 G4double rateMassHe2p;
171 G4double theZieglerFactor;
172
173 G4int iMolecula; // index in the molecula's table
174 G4int iASTAR; // index in ASTAR
175 G4bool isIon;
176 G4bool isInitialised;
177};
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
181inline void G4BraggIonModel::SetParticle(const G4ParticleDefinition* p)
182{
183 particle = p;
184 mass = particle->GetPDGMass();
185 spin = particle->GetPDGSpin();
186 G4double q = particle->GetPDGCharge()/CLHEP::eplus;
187 chargeSquare = q*q;
188 massRate = mass/CLHEP::proton_mass_c2;
189 ratio = CLHEP::electron_mass_c2/mass;
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193
194#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
virtual G4double GetParticleCharge(const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual ~G4BraggIonModel()
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4double GetPDGCharge() const