source: trunk/source/processes/electromagnetic/standard/include/G4BraggIonModel.hh@ 966

Last change on this file since 966 was 961, checked in by garnier, 17 years ago

update processes

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