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

Last change on this file since 1340 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 8.4 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.34 2009/08/11 10:29:30 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
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// 08-04-08 Simplify Select method for only one G4RegionModel (VI)
56// 03-08-09 Removed unused members and simplify model search if only one
57//          model is used (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
82class G4RegionModels
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
94  ~G4RegionModels();
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
121  const G4Region*    theRegion;
122  G4int              nModelsForRegion;
123  G4int*             theListOfModelIndexes;
124  G4double*          lowKineticEnergy;
125
126};
127
128#include "G4VEmModel.hh"
129#include "G4VEmFluctuationModel.hh"
130#include "G4DynamicParticle.hh"
131
132class G4Region;
133class G4ParticleDefinition;
134class G4PhysicsVector;
135class G4MaterialCutsCouple;
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139class G4EmModelManager
140{
141public:
142
143  G4EmModelManager();
144
145  ~G4EmModelManager();
146
147  void Clear();
148
149  const G4DataVector* Initialise(const G4ParticleDefinition*,
150                                 const G4ParticleDefinition*,
151                                       G4double,
152                                       G4int);
153
154  void FillDEDXVector(G4PhysicsVector*, const G4MaterialCutsCouple*, 
155                      G4EmTableType t = fRestricted);
156
157  void FillLambdaVector(G4PhysicsVector*, const G4MaterialCutsCouple*, 
158                        G4bool startFromNull = true, G4EmTableType t = fRestricted);
159
160  G4VEmModel* GetModel(G4int, G4bool ver = false);
161
162  void AddEmModel(G4int, G4VEmModel*, G4VEmFluctuationModel*, const G4Region*);
163 
164  void UpdateEmModel(const G4String&, G4double, G4double);
165
166  void DumpModelList(G4int verb);
167
168  inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
169
170  inline const G4DataVector* Cuts() const;
171
172  inline const G4DataVector* SubCutoff() const;
173
174  inline G4int NumberOfModels() const;
175
176private:
177
178  inline G4double ComputeDEDX(G4VEmModel* model,
179                              const G4MaterialCutsCouple*,
180                              G4double kinEnergy,
181                              G4double cutEnergy,
182                              G4double minEnergy);
183
184  // hide  assignment operator
185
186  G4EmModelManager(G4EmModelManager &);
187  G4EmModelManager & operator=(const G4EmModelManager &right);
188
189// =====================================================================
190
191private:
192
193  G4DataVector                theCuts;
194  G4DataVector                theSubCuts;
195
196  std::vector<G4VEmModel*>                models;
197  std::vector<G4VEmFluctuationModel*>     flucModels;
198  std::vector<const G4Region*>            regions;
199  std::vector<G4int>                      orderOfModels;
200  std::vector<G4int>                      isUsed;
201
202  G4int                       nEmModels;
203  G4int                       nRegions;
204
205  std::vector<G4int>            idxOfRegionModels;
206  std::vector<G4RegionModels*>  setOfRegionModels;
207
208  G4double                    maxSubCutInRange;
209
210  const G4ParticleDefinition* particle;
211
212  G4int                       verboseLevel;
213  G4bool                      severalModels;
214
215  // may be changed in run time
216  G4RegionModels*             currRegionModel;
217  G4VEmModel*                 currModel;
218};
219
220//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
222
223inline G4VEmModel* G4EmModelManager::SelectModel(G4double& kinEnergy, 
224                                                 size_t& index)
225{
226  if(severalModels) {
227    if(nRegions > 1) {
228      currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
229    }
230    currModel = models[currRegionModel->SelectIndex(kinEnergy)];
231  }
232  return currModel;
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236
237inline const G4DataVector* G4EmModelManager::Cuts() const
238{
239  return &theCuts;
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
243
244inline const G4DataVector* G4EmModelManager::SubCutoff() const
245{
246  return &theSubCuts;
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
250
251inline G4int G4EmModelManager::NumberOfModels() const
252{
253  return nEmModels;
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
257
258inline G4double
259G4EmModelManager::ComputeDEDX(G4VEmModel* model,
260                              const G4MaterialCutsCouple* couple,
261                              G4double e,
262                              G4double cut,
263                              G4double emin)
264{
265  G4double dedx = 0.0;
266  if(model && cut > emin) {
267    dedx = model->ComputeDEDX(couple,particle,e,cut); 
268    if(emin > 0.0) {dedx -= model->ComputeDEDX(couple,particle,e,emin);} 
269  }
270  return dedx;
271}
272
273//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
274
275#endif
276
Note: See TracBrowser for help on using the repository browser.