source: trunk/source/processes/electromagnetic/utils/include/G4EmCalculator.hh@ 1197

Last change on this file since 1197 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 10.9 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: G4EmCalculator.hh,v 1.19 2009/11/11 23:59:48 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29//
30// -------------------------------------------------------------------
31//
32// GEANT4 Class header file
33//
34//
35// File name: G4EmCalculator
36//
37// Author: Vladimir Ivanchenko
38//
39// Creation date: 27.06.2004
40//
41// Modifications:
42// 17.11.2004 Change signature of methods, add new methods (V.Ivanchenko)
43// 11.01.2006 Add GetCSDARange (V.Ivanchenko)
44// 26.01.2006 Rename GetRange -> GetRangeFromRestricteDEDX (V.Ivanchenko)
45// 22.03.2006 Add ComputeElectronicDEDX and ComputeTotalDEDX (V.Ivanchenko)
46// 29.09.2006 Add member loweModel (V.Ivanchenko)
47// 15.03.2007 Add ComputeEnergyCutFromRangeCut methods (V.Ivanchenko)
48//
49// Class Description:
50//
51// Provide access to dE/dx and cross sections
52
53// -------------------------------------------------------------------
54//
55
56#ifndef G4EmCalculator_h
57#define G4EmCalculator_h 1
58
59#include <vector>
60#include "globals.hh"
61#include "G4DataVector.hh"
62#include "G4DynamicParticle.hh"
63
64class G4LossTableManager;
65class G4Material;
66class G4MaterialCutsCouple;
67class G4ParticleDefinition;
68class G4PhysicsTable;
69class G4VEmModel;
70class G4VEnergyLossProcess;
71class G4ionEffectiveCharge;
72class G4Region;
73class G4Element;
74class G4EmCorrections;
75
76class G4EmCalculator
77{
78
79public:
80
81 G4EmCalculator();
82
83 ~G4EmCalculator();
84
85 //==================================================================================
86 // Methods to access precalculated dE/dx and cross sections
87 // Materials should exist in the list of the G4MaterialCutsCouple
88 //==================================================================================
89
90 G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*, const G4Material*,
91 const G4Region* r = 0);
92 G4double GetDEDX(G4double kinEnergy, const G4String& part, const G4String& mat,
93 const G4String& s = "world");
94
95 G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4ParticleDefinition*,
96 const G4Material*,
97 const G4Region* r = 0);
98 G4double GetRangeFromRestricteDEDX(G4double kinEnergy, const G4String& part,
99 const G4String& mat,
100 const G4String& s = "world");
101
102 G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition*,
103 const G4Material*,
104 const G4Region* r = 0);
105 G4double GetCSDARange(G4double kinEnergy, const G4String& part,
106 const G4String& mat,
107 const G4String& s = "world");
108
109 G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*,
110 const G4Material*,
111 const G4Region* r = 0);
112 G4double GetRange(G4double kinEnergy, const G4String& part,
113 const G4String& mat,
114 const G4String& s = "world");
115
116 G4double GetKinEnergy(G4double range, const G4ParticleDefinition*, const G4Material*,
117 const G4Region* r = 0);
118 G4double GetKinEnergy(G4double range, const G4String& part, const G4String& mat,
119 const G4String& s = "world");
120
121 G4double GetCrossSectionPerVolume(
122 G4double kinEnergy, const G4ParticleDefinition*,
123 const G4String& processName, const G4Material*,
124 const G4Region* r = 0);
125 G4double GetCrossSectionPerVolume(
126 G4double kinEnergy, const G4String& part, const G4String& proc,
127 const G4String& mat, const G4String& s = "world");
128
129 G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition*,
130 const G4String& processName, const G4Material*,
131 const G4Region* r = 0);
132 G4double GetMeanFreePath(G4double kinEnergy, const G4String& part, const G4String& proc,
133 const G4String& mat, const G4String& s = "world");
134
135 void PrintDEDXTable(const G4ParticleDefinition*);
136
137 void PrintRangeTable(const G4ParticleDefinition*);
138
139 void PrintInverseRangeTable(const G4ParticleDefinition*);
140
141 //==================================================================================
142 // Methods to calculate dE/dx and cross sections "on fly"
143 // Existing tables and G4MaterialCutsCouples are not used
144 //==================================================================================
145
146 G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition*,
147 const G4String& processName, const G4Material*,
148 G4double cut = DBL_MAX);
149 G4double ComputeDEDX(G4double kinEnergy, const G4String& part, const G4String& proc,
150 const G4String& mat, G4double cut = DBL_MAX);
151
152 G4double ComputeElectronicDEDX(G4double kinEnergy, const G4ParticleDefinition*,
153 const G4Material* mat, G4double cut = DBL_MAX);
154 G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
155 const G4String& mat, G4double cut = DBL_MAX);
156
157 G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition*,
158 const G4Material*);
159 G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
160 const G4String& mat);
161
162 G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition*,
163 const G4Material*, G4double cut = DBL_MAX);
164 G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
165 const G4String& mat, G4double cut = DBL_MAX);
166
167 G4double ComputeCrossSectionPerVolume(
168 G4double kinEnergy, const G4ParticleDefinition*,
169 const G4String& processName, const G4Material*,
170 G4double cut = 0.0);
171 G4double ComputeCrossSectionPerVolume(
172 G4double kinEnergy, const G4String& part, const G4String& proc,
173 const G4String& mat, G4double cut = 0.0);
174
175 G4double ComputeCrossSectionPerAtom(
176 G4double kinEnergy, const G4ParticleDefinition*,
177 const G4String& processName, G4double Z, G4double A,
178 G4double cut = 0.0);
179 G4double ComputeCrossSectionPerAtom(
180 G4double kinEnergy, const G4String& part,
181 const G4String& processName, const G4Element*,
182 G4double cut = 0.0);
183
184 G4double ComputeMeanFreePath(
185 G4double kinEnergy, const G4ParticleDefinition*,
186 const G4String& processName, const G4Material*,
187 G4double cut = 0.0);
188 G4double ComputeMeanFreePath(
189 G4double kinEnergy, const G4String&, const G4String&,
190 const G4String& processName, G4double cut = 0.0);
191
192 G4double ComputeEnergyCutFromRangeCut(
193 G4double range, const G4ParticleDefinition*,
194 const G4Material*);
195 G4double ComputeEnergyCutFromRangeCut(
196 G4double range, const G4String&,
197 const G4String&);
198
199 //==================================================================================
200 // Methods to access particles, materials, regions
201 //==================================================================================
202
203 const G4ParticleDefinition* FindParticle(const G4String&);
204
205 const G4ParticleDefinition* FindIon(G4int Z, G4int A);
206
207 const G4Material* FindMaterial(const G4String&);
208
209 const G4Region* FindRegion(const G4String&);
210
211 const G4MaterialCutsCouple* FindCouple(const G4Material*, const G4Region* r = 0);
212
213 void SetVerbose(G4int val);
214
215 //==================================================================================
216 // Private methods
217 //==================================================================================
218
219private:
220
221 G4bool UpdateParticle(const G4ParticleDefinition*, G4double kinEnergy);
222
223 G4bool UpdateCouple(const G4Material*, G4double cut);
224
225 void FindLambdaTable(const G4ParticleDefinition*, const G4String& processName);
226
227 G4bool FindEmModel(const G4ParticleDefinition*,
228 const G4String& processName,
229 G4double kinEnergy);
230
231 G4VEnergyLossProcess* FindEnergyLossProcess(const G4ParticleDefinition*);
232
233 G4EmCalculator & operator=(const G4EmCalculator &right);
234 G4EmCalculator(const G4EmCalculator&);
235
236 std::vector<const G4Material*> localMaterials;
237 std::vector<const G4MaterialCutsCouple*> localCouples;
238
239 G4LossTableManager* manager;
240 G4EmCorrections* corr;
241 G4DataVector localCuts;
242 G4int nLocalMaterials;
243
244 G4int verbose;
245
246 // cash
247 G4int currentCoupleIndex;
248 const G4MaterialCutsCouple* currentCouple;
249 const G4Material* currentMaterial;
250 const G4ParticleDefinition* currentParticle;
251 const G4ParticleDefinition* baseParticle;
252 const G4PhysicsTable* currentLambda;
253 G4VEmModel* currentModel;
254 G4VEmModel* loweModel;
255 G4VEnergyLossProcess* currentProcess;
256
257 const G4ParticleDefinition* theGenericIon;
258 G4ionEffectiveCharge* ionEffCharge;
259 G4DynamicParticle dynParticle;
260
261 G4String currentName;
262 G4double currentCut;
263 G4double chargeSquare;
264 G4double massRatio;
265 G4double mass;
266 G4bool isIon;
267 G4bool isApplicable;
268
269 G4String currentParticleName;
270 G4String currentMaterialName;
271 G4String currentProcessName;
272};
273
274//....oooOO0OOooo.......oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
275
276#endif
Note: See TracBrowser for help on using the repository browser.