source: trunk/source/processes/electromagnetic/utils/include/G4EmModelManager.hh@ 1046

Last change on this file since 1046 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 7.5 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: G4EmModelManager.hh,v 1.25 2008/10/13 14:56:56 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33// File name: G4EmModelManager
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 07.05.2002
38//
39// Modifications:
40//
41// 03-12-02 V.Ivanchenko fix a bug in model selection
42// 20-01-03 Migrade to cut per region (V.Ivanchenko)
43// 27-01-03 Make models region aware (V.Ivanchenko)
44// 13-02-03 The set of models is defined for region (V.Ivanchenko)
45// 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
46// 13-04-03 Add startFromNull (V.Ivanchenko)
47// 13-05-03 Add calculation of precise range (V.Ivanchenko)
48// 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
49// 03-11-03 Substitute STL vector for G4RegionModels (V.Ivanchenko)
50// 11-04-05 Remove access to fluctuation models (V.Ivanchenko)
51// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
52// 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
53// 13-05-06 Add GetModel by index method (VI)
54// 15-03-07 Add maxCutInRange (V.Ivanchenko)
55//
56// Class Description:
57//
58// It is the unified energy loss process it calculates the continuous
59// energy loss for charged particles using a set of Energy Loss
60// models valid for different energy regions. There are a possibility
61// to create and access to dE/dx and range tables, or to calculate
62// that information on fly.
63
64// -------------------------------------------------------------------
65//
66
67
68#ifndef G4EmModelManager_h
69#define G4EmModelManager_h 1
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73#include "globals.hh"
74#include "G4DataVector.hh"
75#include "G4EmTableType.hh"
76#include "G4EmProcessSubType.hh"
77#include "G4Region.hh"
78
79class G4RegionModels
80{
81
82friend class G4EmModelManager;
83
84public:
85
86private:
87
88 G4RegionModels(G4int nMod, std::vector<G4int>& indx,
89 G4DataVector& lowE, const G4Region* reg);
90
91 ~G4RegionModels();
92
93 G4int SelectIndex(G4double e) const {
94 G4int idx = 0;
95 if (nModelsForRegion>1) {
96 idx = nModelsForRegion;
97 do {idx--;} while (idx && e <= lowKineticEnergy[idx]);
98 }
99 return theListOfModelIndexes[idx];
100 };
101
102 G4int ModelIndex(G4int n) const {
103 return theListOfModelIndexes[n];
104 };
105
106 G4int NumberOfModels() const {
107 return nModelsForRegion;
108 };
109
110 G4double LowEdgeEnergy(G4int n) const {
111 return lowKineticEnergy[n];
112 };
113
114 const G4Region* Region() const {
115 return theRegion;
116 };
117
118 const G4Region* theRegion;
119 G4int nModelsForRegion;
120 G4int* theListOfModelIndexes;
121 G4double* lowKineticEnergy;
122
123};
124
125#include "G4VEmModel.hh"
126#include "G4VEmFluctuationModel.hh"
127#include "G4DynamicParticle.hh"
128
129class G4Region;
130class G4ParticleDefinition;
131class G4DataVector;
132class G4PhysicsVector;
133class G4MaterialCutsCouple;
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137class G4EmModelManager
138{
139public:
140
141 G4EmModelManager();
142
143 ~G4EmModelManager();
144
145 void Clear();
146
147 const G4DataVector* Initialise(const G4ParticleDefinition*,
148 const G4ParticleDefinition*,
149 G4double,
150 G4int);
151
152 const G4DataVector* Cuts() const;
153
154 const G4DataVector* SubCutoff() const;
155
156 void FillDEDXVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
157 G4EmTableType t = fRestricted);
158
159 void FillLambdaVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
160 G4bool startFromNull = true, G4EmTableType t = fRestricted);
161
162 G4VEmModel* SelectModel(G4double& energy, size_t& index);
163
164 G4VEmModel* GetModel(G4int, G4bool ver = false);
165
166 G4int NumberOfModels() const;
167
168 void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel*, const G4Region*);
169
170 void UpdateEmModel(const G4String&, G4double, G4double);
171
172 void DumpModelList(G4int verb);
173
174private:
175
176 // hide assignment operator
177
178 G4EmModelManager(G4EmModelManager &);
179 G4EmModelManager & operator=(const G4EmModelManager &right);
180
181// =====================================================================
182
183private:
184
185 G4DataVector theCuts;
186 G4DataVector theSubCuts;
187
188 std::vector<G4VEmModel*> models;
189 std::vector<G4VEmFluctuationModel*> flucModels;
190 std::vector<const G4Region*> regions;
191 std::vector<G4int> orderOfModels;
192
193 G4int nEmModels;
194 G4int nRegions;
195 G4int nCouples;
196
197 G4int* idxOfRegionModels;
198 G4RegionModels** setOfRegionModels;
199
200 G4double minSubRange;
201 G4double maxCutInRange;
202 G4double maxSubCutInRange;
203
204 const G4ParticleDefinition* particle;
205 const G4ParticleDefinition* secondaryParticle;
206 const G4ParticleDefinition* theGamma;
207 const G4ParticleDefinition* thePositron;
208
209 G4int verboseLevel;
210
211 // cash
212 G4int currentIdx;
213
214};
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218
219inline G4VEmModel* G4EmModelManager::SelectModel(G4double& kinEnergy,
220 size_t& index)
221{
222 currentIdx =
223 (setOfRegionModels[idxOfRegionModels[index]])->SelectIndex(kinEnergy);
224 return models[currentIdx];
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
228
229inline const G4DataVector* G4EmModelManager::Cuts() const
230{
231 return &theCuts;
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
235
236inline const G4DataVector* G4EmModelManager::SubCutoff() const
237{
238 return &theSubCuts;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242
243inline G4int G4EmModelManager::NumberOfModels() const
244{
245 return nEmModels;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
250#endif
251
Note: See TracBrowser for help on using the repository browser.