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
G4EmModelManager.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// File name: G4EmModelManager
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 07.05.2002
37//
38// Modifications:
39//
40// 03-12-02 V.Ivanchenko fix a bug in model selection
41// 20-01-03 Migrade to cut per region (V.Ivanchenko)
42// 27-01-03 Make models region aware (V.Ivanchenko)
43// 13-02-03 The set of models is defined for region (V.Ivanchenko)
44// 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
45// 13-04-03 Add startFromNull (V.Ivanchenko)
46// 13-05-03 Add calculation of precise range (V.Ivanchenko)
47// 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
48// 03-11-03 Substitute STL vector for G4RegionModels (V.Ivanchenko)
49// 11-04-05 Remove access to fluctuation models (V.Ivanchenko)
50// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
51// 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
52// 13-05-06 Add GetModel by index method (VI)
53// 15-03-07 Add maxCutInRange (V.Ivanchenko)
54// 08-04-08 Simplify Select method for only one G4RegionModel (VI)
55// 03-08-09 Removed unused members and simplify model search if only one
56// model is used (VI)
57// 14-07-11 Use pointer to the vector of cuts and not local copy (VI)
58//
59// Class Description:
60//
61// It is the unified energy loss process it calculates the continuous
62// energy loss for charged particles using a set of Energy Loss
63// models valid for different energy regions. There are a possibility
64// to create and access to dE/dx and range tables, or to calculate
65// that information on fly.
66
67// -------------------------------------------------------------------
68//
69
70
71#ifndef G4EmModelManager_h
72#define G4EmModelManager_h 1
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
76#include "globals.hh"
77#include "G4DataVector.hh"
78#include "G4EmTableType.hh"
79#include "G4EmProcessSubType.hh"
80#include "G4Region.hh"
81
83{
84
85friend class G4EmModelManager;
86
87public:
88
89private:
90
91 G4RegionModels(G4int nMod, std::vector<G4int>& indx,
92 G4DataVector& lowE, const G4Region* reg);
93
95
96 inline G4int SelectIndex(G4double e) const {
97 G4int idx = 0;
98 if (nModelsForRegion>1) {
99 idx = nModelsForRegion;
100 do {--idx;} while (idx > 0 && e <= lowKineticEnergy[idx]);
101 }
102 return theListOfModelIndexes[idx];
103 };
104
105 inline G4int ModelIndex(G4int n) const {
106 return theListOfModelIndexes[n];
107 };
108
109 inline G4int NumberOfModels() const {
110 return nModelsForRegion;
111 };
112
113 inline G4double LowEdgeEnergy(G4int n) const {
114 return lowKineticEnergy[n];
115 };
116
117 inline const G4Region* Region() const {
118 return theRegion;
119 };
120
122 G4RegionModels & operator=(const G4RegionModels &right);
123
124 const G4Region* theRegion;
125 G4int nModelsForRegion;
126 G4int* theListOfModelIndexes;
127 G4double* lowKineticEnergy;
128
129};
130
131#include "G4VEmModel.hh"
133#include "G4DynamicParticle.hh"
134
135class G4Region;
137class G4PhysicsVector;
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
143{
144public:
145
147
149
150 void Clear();
151
154 G4double,
155 G4int);
156
159
161 G4bool startFromNull = true, G4EmTableType t = fRestricted);
162
163 G4VEmModel* GetModel(G4int, G4bool ver = false);
164
166
168
169 void DumpModelList(G4int verb);
170
171 inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
172
173 inline const G4DataVector* Cuts() const;
174
175 inline const G4DataVector* SubCutoff() const;
176
177 inline void SetFluoFlag(G4bool val);
178
179 inline G4int NumberOfModels() const;
180
181private:
182
183 inline G4double ComputeDEDX(G4VEmModel* model,
185 G4double kinEnergy,
186 G4double cutEnergy,
187 G4double minEnergy);
188
189 // hide assignment operator
190
192 G4EmModelManager & operator=(const G4EmModelManager &right);
193
194// =====================================================================
195
196private:
197
198 const G4DataVector* theCuts;
199 G4DataVector* theSubCuts;
200
201 std::vector<G4VEmModel*> models;
202 std::vector<G4VEmFluctuationModel*> flucModels;
203 std::vector<const G4Region*> regions;
204 std::vector<G4int> orderOfModels;
205 std::vector<G4int> isUsed;
206
207 G4int nEmModels;
208 G4int nRegions;
209
210 std::vector<G4int> idxOfRegionModels;
211 std::vector<G4RegionModels*> setOfRegionModels;
212
213 G4double maxSubCutInRange;
214
215 const G4ParticleDefinition* particle;
216
217 G4int verboseLevel;
218 G4bool severalModels;
219 G4bool fluoFlag;
220
221 // may be changed in run time
222 G4RegionModels* currRegionModel;
223 G4VEmModel* currModel;
224};
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
230 size_t& index)
231{
232 if(severalModels) {
233 if(nRegions > 1) {
234 currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
235 }
236 currModel = models[currRegionModel->SelectIndex(kinEnergy)];
237 }
238 return currModel;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242
244{
245 return theCuts;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
251{
252 return theSubCuts;
253}
254
255//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256
258{
259 fluoFlag = val;
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
263
265{
266 return nEmModels;
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
270
271inline G4double
272G4EmModelManager::ComputeDEDX(G4VEmModel* model,
273 const G4MaterialCutsCouple* couple,
274 G4double e,
275 G4double cut,
276 G4double emin)
277{
278 G4double dedx = 0.0;
279 if(model && cut > emin) {
280 dedx = model->ComputeDEDX(couple,particle,e,cut);
281 if(emin > 0.0) {dedx -= model->ComputeDEDX(couple,particle,e,emin);}
282 }
283 return dedx;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
288#endif
289
G4EmTableType
@ fRestricted
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void UpdateEmModel(const G4String &, G4double, G4double)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4int NumberOfModels() const
void SetFluoFlag(G4bool val)
void DumpModelList(G4int verb)
G4VEmModel * GetModel(G4int, G4bool ver=false)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
G4VEmModel * SelectModel(G4double &energy, size_t &index)
const G4DataVector * SubCutoff() const
const G4DataVector * Cuts() const
const G4DataVector * Initialise(const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4int)
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:407