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

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

import all except CVS

File size: 7.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: G4EmModelManager.hh,v 1.22 2007/11/09 11:35:54 vnivanch Exp $
27// GEANT4 tag $Name: $
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
77class G4RegionModels
78{
79
80friend class G4EmModelManager;
81
82public:
83
84private:
85
86 G4RegionModels(G4int nMod, std::vector<G4int>& list, G4DataVector& lowE);
87
88 ~G4RegionModels();
89
90 G4int SelectIndex(G4double e) const {
91 G4int idx = 0;
92 if (nModelsForRegion>1) {
93 idx = nModelsForRegion;
94 do {idx--;} while (idx && e <= lowKineticEnergy[idx]);
95 }
96 return theListOfModelIndexes[idx];
97 };
98
99 G4int ModelIndex(G4int n) const {
100 return theListOfModelIndexes[n];
101 };
102
103 G4int NumberOfModels() const {
104 return nModelsForRegion;
105 };
106
107 G4int nModelsForRegion;
108 G4int* theListOfModelIndexes;
109 G4double* lowKineticEnergy;
110
111};
112
113#include "G4VEmModel.hh"
114#include "G4VEmFluctuationModel.hh"
115#include "G4DynamicParticle.hh"
116
117class G4Region;
118class G4ParticleDefinition;
119class G4DataVector;
120class G4PhysicsVector;
121class G4MaterialCutsCouple;
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124
125class G4EmModelManager
126{
127public:
128
129 G4EmModelManager();
130
131 ~G4EmModelManager();
132
133 void Clear();
134
135 const G4DataVector* Initialise(const G4ParticleDefinition*,
136 const G4ParticleDefinition*,
137 G4double,
138 G4int);
139
140 const G4DataVector* Cuts() const;
141
142 const G4DataVector* SubCutoff() const;
143
144 void FillDEDXVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
145 G4EmTableType t = fRestricted);
146
147 void FillLambdaVector(G4PhysicsVector*, const G4MaterialCutsCouple*,
148 G4bool startFromNull = true, G4EmTableType t = fRestricted);
149
150 G4VEmModel* SelectModel(G4double& energy, size_t& index);
151
152 G4VEmModel* GetModel(G4int);
153
154 G4int NumberOfModels() const;
155
156 void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel*, const G4Region*);
157
158 void UpdateEmModel(const G4String&, G4double, G4double);
159
160private:
161
162 // hide assignment operator
163
164 G4EmModelManager(G4EmModelManager &);
165 G4EmModelManager & operator=(const G4EmModelManager &right);
166
167// =====================================================================
168
169private:
170
171 G4DataVector theCuts;
172 G4DataVector theSubCuts;
173
174 std::vector<G4VEmModel*> models;
175 std::vector<G4VEmFluctuationModel*> flucModels;
176 std::vector<const G4Region*> regions;
177 std::vector<G4int> orderOfModels;
178 G4DataVector upperEkin;
179
180 G4int nEmModels;
181 G4int nRegions;
182 G4int nCouples;
183
184 G4int* idxOfRegionModels;
185 G4RegionModels** setOfRegionModels;
186
187 G4double minSubRange;
188 G4double maxCutInRange;
189 G4double maxSubCutInRange;
190
191 const G4ParticleDefinition* particle;
192 const G4ParticleDefinition* secondaryParticle;
193 const G4ParticleDefinition* theGamma;
194 const G4ParticleDefinition* thePositron;
195
196 G4int verboseLevel;
197
198 // cash
199 G4int currentIdx;
200
201};
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
205
206inline G4VEmModel* G4EmModelManager::SelectModel(G4double& kinEnergy,
207 size_t& index)
208{
209 currentIdx =
210 (setOfRegionModels[idxOfRegionModels[index]])->SelectIndex(kinEnergy);
211 return models[currentIdx];
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215
216inline const G4DataVector* G4EmModelManager::Cuts() const
217{
218 return &theCuts;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
223inline const G4DataVector* G4EmModelManager::SubCutoff() const
224{
225 return &theSubCuts;
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
230inline G4int G4EmModelManager::NumberOfModels() const
231{
232 return nEmModels;
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236
237#endif
238
Note: See TracBrowser for help on using the repository browser.