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

Last change on this file since 1066 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

File size: 7.8 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//
[1055]26// $Id: G4EmModelManager.hh,v 1.28 2009/04/09 15:53:17 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
[819]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)
[1055]55// 08-04-08 Simplify Select method for only one G4RegionModel (VI)
[819]56//
57// Class Description:
58//
59// It is the unified energy loss process it calculates the continuous
60// energy loss for charged particles using a set of Energy Loss
61// models valid for different energy regions. There are a possibility
62// to create and access to dE/dx and range tables, or to calculate
63// that information on fly.
64
65// -------------------------------------------------------------------
66//
67
68
69#ifndef G4EmModelManager_h
70#define G4EmModelManager_h 1
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
74#include "globals.hh"
75#include "G4DataVector.hh"
76#include "G4EmTableType.hh"
[961]77#include "G4EmProcessSubType.hh"
78#include "G4Region.hh"
[819]79
80class G4RegionModels
81{
82
83friend class G4EmModelManager;
84
85public:
86
87private:
88
[961]89 G4RegionModels(G4int nMod, std::vector<G4int>& indx,
90 G4DataVector& lowE, const G4Region* reg);
[819]91
92 ~G4RegionModels();
93
[1055]94 inline G4int SelectIndex(G4double e) const {
[819]95 G4int idx = 0;
96 if (nModelsForRegion>1) {
97 idx = nModelsForRegion;
98 do {idx--;} while (idx && e <= lowKineticEnergy[idx]);
99 }
100 return theListOfModelIndexes[idx];
101 };
102
[1055]103 inline G4int ModelIndex(G4int n) const {
[819]104 return theListOfModelIndexes[n];
105 };
106
[1055]107 inline G4int NumberOfModels() const {
[819]108 return nModelsForRegion;
109 };
110
[1055]111 inline G4double LowEdgeEnergy(G4int n) const {
[961]112 return lowKineticEnergy[n];
113 };
114
[1055]115 inline const G4Region* Region() const {
[961]116 return theRegion;
117 };
118
119 const G4Region* theRegion;
[819]120 G4int nModelsForRegion;
121 G4int* theListOfModelIndexes;
122 G4double* lowKineticEnergy;
123
124};
125
126#include "G4VEmModel.hh"
127#include "G4VEmFluctuationModel.hh"
128#include "G4DynamicParticle.hh"
129
130class G4Region;
131class G4ParticleDefinition;
132class G4DataVector;
133class G4PhysicsVector;
134class G4MaterialCutsCouple;
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138class G4EmModelManager
139{
140public:
141
142 G4EmModelManager();
143
144 ~G4EmModelManager();
145
146 void Clear();
147
148 const G4DataVector* Initialise(const G4ParticleDefinition*,
149 const G4ParticleDefinition*,
150 G4double,
151 G4int);
152
153 void FillDEDXVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
154 G4EmTableType t = fRestricted);
155
156 void FillLambdaVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
157 G4bool startFromNull = true, G4EmTableType t = fRestricted);
158
[961]159 G4VEmModel* GetModel(G4int, G4bool ver = false);
[819]160
161 void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel*, const G4Region*);
162
163 void UpdateEmModel(const G4String&, G4double, G4double);
164
[961]165 void DumpModelList(G4int verb);
166
[1055]167 inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
168
169 inline const G4DataVector* Cuts() const;
170
171 inline const G4DataVector* SubCutoff() const;
172
173 inline G4int NumberOfModels() const;
174
[819]175private:
176
177 // hide assignment operator
178
179 G4EmModelManager(G4EmModelManager &);
180 G4EmModelManager & operator=(const G4EmModelManager &right);
181
182// =====================================================================
183
184private:
185
186 G4DataVector theCuts;
187 G4DataVector theSubCuts;
188
189 std::vector<G4VEmModel*> models;
190 std::vector<G4VEmFluctuationModel*> flucModels;
191 std::vector<const G4Region*> regions;
192 std::vector<G4int> orderOfModels;
[1055]193 std::vector<G4int> isUsed;
[819]194
195 G4int nEmModels;
196 G4int nRegions;
197 G4int nCouples;
198
199 G4int* idxOfRegionModels;
200 G4RegionModels** setOfRegionModels;
201
202 G4double minSubRange;
203 G4double maxCutInRange;
204 G4double maxSubCutInRange;
205
206 const G4ParticleDefinition* particle;
207 const G4ParticleDefinition* secondaryParticle;
208 const G4ParticleDefinition* theGamma;
209 const G4ParticleDefinition* thePositron;
210
211 G4int verboseLevel;
212
[1055]213 // may be changed in run time
214 G4RegionModels* currRegionModel;
[819]215};
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
219
220inline G4VEmModel* G4EmModelManager::SelectModel(G4double& kinEnergy,
221 size_t& index)
222{
[1055]223 if(nRegions > 1) currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
224 return models[currRegionModel->SelectIndex(kinEnergy)];
[819]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.