source: trunk/source/processes/electromagnetic/standard/include/G4BraggModel.hh@ 900

Last change on this file since 900 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.3 KB
RevLine 
[819]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: G4BraggModel.hh,v 1.10 2007/05/22 17:34:36 vnivanch Exp $
27// GEANT4 tag $Name: $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4BraggModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41// 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
42// 24-01-03 Make models region aware (V.Ivanchenko)
43// 13-02-03 Add name (V.Ivanchenko)
44// 12-11-03 Fix for GenericIons (V.Ivanchenko)
45// 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
46// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
47// 25-04-06 Add stopping data from PSTAR (V.Ivanchenko)
48
49//
50// Class Description:
51//
52// Implementation of energy loss and delta-electron production
53// by heavy slow charged particles using eveluated data
54
55// -------------------------------------------------------------------
56//
57
58#ifndef G4BraggModel_h
59#define G4BraggModel_h 1
60
61#include "G4VEmModel.hh"
62#include "G4PSTARStopping.hh"
63
64class G4ParticleChangeForLoss;
65
66class G4BraggModel : public G4VEmModel
67{
68
69public:
70
71 G4BraggModel(const G4ParticleDefinition* p = 0,
72 const G4String& nam = "Bragg");
73
74 virtual ~G4BraggModel();
75
76 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
77
78 G4double MinEnergyCut(const G4ParticleDefinition*,
79 const G4MaterialCutsCouple*);
80
81 virtual G4double ComputeCrossSectionPerElectron(
82 const G4ParticleDefinition*,
83 G4double kineticEnergy,
84 G4double cutEnergy,
85 G4double maxEnergy);
86
87 virtual G4double ComputeCrossSectionPerAtom(
88 const G4ParticleDefinition*,
89 G4double kineticEnergy,
90 G4double Z, G4double A,
91 G4double cutEnergy,
92 G4double maxEnergy);
93
94 virtual G4double CrossSectionPerVolume(const G4Material*,
95 const G4ParticleDefinition*,
96 G4double kineticEnergy,
97 G4double cutEnergy,
98 G4double maxEnergy);
99
100 virtual G4double ComputeDEDXPerVolume(const G4Material*,
101 const G4ParticleDefinition*,
102 G4double kineticEnergy,
103 G4double cutEnergy);
104
105 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
106 const G4MaterialCutsCouple*,
107 const G4DynamicParticle*,
108 G4double tmin,
109 G4double maxEnergy);
110
111protected:
112
113 G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
114 G4double kinEnergy);
115
116private:
117
118 void SetParticle(const G4ParticleDefinition* p);
119
120 G4bool HasMaterial(const G4Material* material);
121
122 G4double StoppingPower(const G4Material* material,
123 G4double kineticEnergy);
124
125 G4double ElectronicStoppingPower(G4double z,
126 G4double kineticEnergy) const;
127
128 void SetMoleculaNumber(G4int number) {iMolecula = number;};
129
130 G4double DEDX(const G4Material* material, G4double kineticEnergy);
131
132 G4bool MolecIsInZiegler1988(const G4Material* material);
133
134 G4double ChemicalFactor(G4double kineticEnergy, G4double eloss125) const;
135
136 void SetExpStopPower125(G4double value) {expStopPower125 = value;};
137
138 // hide assignment operator
139 G4BraggModel & operator=(const G4BraggModel &right);
140 G4BraggModel(const G4BraggModel&);
141
142 const G4ParticleDefinition* particle;
143 G4ParticleDefinition* theElectron;
144 G4ParticleChangeForLoss* fParticleChange;
145 G4PSTARStopping pstar;
146
147 G4double mass;
148 G4double spin;
149 G4double chargeSquare;
150 G4double massRate;
151 G4double ratio;
152 G4double highKinEnergy;
153 G4double lowKinEnergy;
154 G4double lowestKinEnergy;
155 G4double protonMassAMU;
156 G4double theZieglerFactor;
157 G4double expStopPower125; // Experimental Stopping power at 125keV
158
159 G4int iMolecula; // index in the molecula's table
160 G4bool isIon;
161};
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
165inline G4double G4BraggModel::MaxSecondaryEnergy(
166 const G4ParticleDefinition* pd,
167 G4double kinEnergy)
168{
169 if(pd != particle) SetParticle(pd);
170 G4double tau = kinEnergy/mass;
171 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
172 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
173 return tmax;
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
178inline void G4BraggModel::SetParticle(const G4ParticleDefinition* p)
179{
180 particle = p;
181 mass = particle->GetPDGMass();
182 spin = particle->GetPDGSpin();
183 G4double q = particle->GetPDGCharge()/eplus;
184 chargeSquare = q*q;
185 massRate = mass/proton_mass_c2;
186 ratio = electron_mass_c2/mass;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
191#endif
Note: See TracBrowser for help on using the repository browser.