Geant4  10.02.p02
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: G4EmModelManager.hh 91745 2015-08-04 11:51:12Z gcosmo $
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 
85 friend class G4EmModelManager;
86 
87 public:
88 
89 private:
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  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
101  do {--idx;} while (idx > 0 && e <= lowKineticEnergy[idx]);
102  }
103  return theListOfModelIndexes[idx];
104  };
105 
106  inline G4int ModelIndex(G4int n) const {
107  return theListOfModelIndexes[n];
108  };
109 
110  inline G4int NumberOfModels() const {
111  return nModelsForRegion;
112  };
113 
114  inline G4double LowEdgeEnergy(G4int n) const {
115  return lowKineticEnergy[n];
116  };
117 
118  inline const G4Region* Region() const {
119  return theRegion;
120  };
121 
124 
129 
130 };
131 
132 #include "G4VEmModel.hh"
133 #include "G4VEmFluctuationModel.hh"
134 #include "G4DynamicParticle.hh"
135 
136 class G4Region;
138 class G4PhysicsVector;
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142 
144 {
145 public:
146 
148 
150 
151  void Clear();
152 
153  const G4DataVector* Initialise(const G4ParticleDefinition* part,
154  const G4ParticleDefinition* secPart,
155  G4double minSubRange,
156  G4int verb);
157 
160 
162  G4bool startFromNull = true,
164 
165  G4VEmModel* GetModel(G4int, G4bool ver = false);
166 
168 
169  void UpdateEmModel(const G4String&, G4double, G4double);
170 
171  void DumpModelList(G4int verb);
172 
173  inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
174 
175  inline const G4DataVector* Cuts() const;
176 
177  inline const G4DataVector* SubCutoff() const;
178 
179  inline void SetFluoFlag(G4bool val);
180 
181  inline G4int NumberOfModels() const;
182 
183 private:
184 
185  inline G4double ComputeDEDX(G4VEmModel* model,
186  const G4MaterialCutsCouple*,
187  G4double kinEnergy,
188  G4double cutEnergy,
189  G4double minEnergy);
190 
191  // hide assignment operator
192 
195 
196 // =====================================================================
197 
198 private:
199 
203 
204  std::vector<G4VEmModel*> models;
205  std::vector<G4VEmFluctuationModel*> flucModels;
206  std::vector<const G4Region*> regions;
207  std::vector<G4int> orderOfModels;
208  std::vector<G4int> isUsed;
209 
212 
213  std::vector<G4int> idxOfRegionModels;
214  std::vector<G4RegionModels*> setOfRegionModels;
215 
217 
219 
223 
224  // may be changed in run time
227 };
228 
229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
231 
233  size_t& index)
234 {
235  if(severalModels) {
236  if(nRegions > 1) {
238  }
240  }
241  return currModel;
242 }
243 
244 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
245 
247 {
248  return theCuts;
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
254 {
255  return theSubCuts;
256 }
257 
258 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259 
261 {
262  fluoFlag = val;
263 }
264 
265 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266 
268 {
269  return nEmModels;
270 }
271 
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
273 
274 inline G4double
276  const G4MaterialCutsCouple* couple,
277  G4double e,
278  G4double cut,
279  G4double emin)
280 {
281  G4double dedx = 0.0;
282  if(model && cut > emin) {
283  dedx = model->ComputeDEDX(couple,particle,e,cut);
284  if(emin > 0.0) {dedx -= model->ComputeDEDX(couple,particle,e,emin);}
285  }
286  return dedx;
287 }
288 
289 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
290 
291 #endif
292 
G4double * lowKineticEnergy
std::vector< G4int > idxOfRegionModels
void UpdateEmModel(const G4String &, G4double, G4double)
const G4DataVector * Cuts() const
G4VEmModel * GetModel(G4int, G4bool ver=false)
G4double LowEdgeEnergy(G4int n) const
const G4DataVector * theCuts
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
const G4DataVector * SubCutoff() const
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
int G4int
Definition: G4Types.hh:78
G4EmModelManager & operator=(const G4EmModelManager &right)
std::vector< G4int > orderOfModels
G4RegionModels & operator=(const G4RegionModels &right)
static const G4double reg
G4int ModelIndex(G4int n) const
G4int NumberOfModels() const
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
std::vector< G4int > isUsed
bool G4bool
Definition: G4Types.hh:79
G4RegionModels(G4int nMod, std::vector< G4int > &indx, G4DataVector &lowE, const G4Region *reg)
const G4ParticleDefinition * particle
const G4int n
void DumpModelList(G4int verb)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
G4DataVector * theSubCuts
std::vector< G4VEmModel * > models
G4int NumberOfModels() const
G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.hh:490
G4double energy(const ThreeVector &p, const G4double m)
G4double ComputeDEDX(G4VEmModel *model, const G4MaterialCutsCouple *, G4double kinEnergy, G4double cutEnergy, G4double minEnergy)
std::vector< const G4Region * > regions
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4EmTableType
void SetFluoFlag(G4bool val)
G4RegionModels * currRegionModel
G4VEmModel * currModel
double G4double
Definition: G4Types.hh:76
const G4Region * theRegion
const G4Region * Region() const
G4int * theListOfModelIndexes
std::vector< G4VEmFluctuationModel * > flucModels
G4int SelectIndex(G4double e) const
G4DataVector * theCutsNew
std::vector< G4RegionModels * > setOfRegionModels