ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmModelManager.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4EmModelManager.hh
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 //
27 // -------------------------------------------------------------------
28 //
29 // GEANT4 Class header file
30 //
31 // File name: G4EmModelManager
32 //
33 // Author: Vladimir Ivanchenko
34 //
35 // Creation date: 07.05.2002
36 //
37 // Modifications:
38 //
39 // 03-12-02 V.Ivanchenko fix a bug in model selection
40 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
41 // 27-01-03 Make models region aware (V.Ivanchenko)
42 // 13-02-03 The set of models is defined for region (V.Ivanchenko)
43 // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
44 // 13-04-03 Add startFromNull (V.Ivanchenko)
45 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
46 // 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
47 // 03-11-03 Substitute STL vector for G4RegionModels (V.Ivanchenko)
48 // 11-04-05 Remove access to fluctuation models (V.Ivanchenko)
49 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
50 // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
51 // 13-05-06 Add GetModel by index method (VI)
52 // 15-03-07 Add maxCutInRange (V.Ivanchenko)
53 // 08-04-08 Simplify Select method for only one G4RegionModel (VI)
54 // 03-08-09 Removed unused members and simplify model search if only one
55 // model is used (VI)
56 // 14-07-11 Use pointer to the vector of cuts and not local copy (VI)
57 //
58 // Class Description:
59 //
60 // It is the unified energy loss process it calculates the continuous
61 // energy loss for charged particles using a set of Energy Loss
62 // models valid for different energy regions. There are a possibility
63 // to create and access to dE/dx and range tables, or to calculate
64 // that information on fly.
65 
66 // -------------------------------------------------------------------
67 //
68 
69 
70 #ifndef G4EmModelManager_h
71 #define G4EmModelManager_h 1
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
75 #include "globals.hh"
76 #include "G4DataVector.hh"
77 #include "G4EmTableType.hh"
78 #include "G4EmProcessSubType.hh"
79 #include "G4Region.hh"
80 
82 {
83 
84 friend class G4EmModelManager;
85 
86 private:
87 
88  G4RegionModels(G4int nMod, std::vector<G4int>& indx,
89  G4DataVector& lowE, const G4Region* reg);
90 
92 
93  inline G4int SelectIndex(G4double e) const {
94  G4int idx = 0;
95  if (nModelsForRegion>1) {
96  idx = nModelsForRegion;
97  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
98  do {--idx;} while (idx > 0 && e <= lowKineticEnergy[idx]);
99  }
100  return theListOfModelIndexes[idx];
101  };
102 
103  inline G4int ModelIndex(G4int n) const {
104  return theListOfModelIndexes[n];
105  };
106 
107  inline G4int NumberOfModels() const {
108  return nModelsForRegion;
109  };
110 
111  inline G4double LowEdgeEnergy(G4int n) const {
112  return lowKineticEnergy[n];
113  };
114 
115  inline const G4Region* Region() const {
116  return theRegion;
117  };
118 
119  G4RegionModels(G4RegionModels &) = delete;
120  G4RegionModels & operator=(const G4RegionModels &right) = delete;
121 
126 
127 };
128 
129 #include "G4VEmModel.hh"
130 #include "G4VEmFluctuationModel.hh"
131 #include "G4DynamicParticle.hh"
132 #include <iostream>
133 
134 class G4Region;
136 class G4PhysicsVector;
138 
139 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140 
142 {
143 public:
144 
146 
148 
149  void Clear();
150 
152  const G4ParticleDefinition* secPart,
153  G4double minSubRange,
154  G4int verb);
155 
158 
160  G4bool startFromNull = true,
162 
164 
165  void UpdateEmModel(const G4String& model_name, G4double emin, G4double emax);
166 
167  // Get model pointer from the model list
168  G4VEmModel* GetModel(G4int idx, G4bool ver = false);
169 
170  // Get model pointer from the model list for a given material cuts couple
171  // no check on material cuts couple index
172  G4VEmModel* GetRegionModel(G4int idx, size_t index_couple);
173 
174  // total number of models for material cut couples
175  // no check on material cuts couple index
176  G4int NumberOfRegionModels(size_t index_couple) const;
177 
178  // Automatic documentation
179  void DumpModelList(std::ostream& out, G4int verb);
180 
181  // Select model for given material cuts couple index
182  inline G4VEmModel* SelectModel(G4double& energy, size_t& index);
183 
184  // Access to cuts
185  inline const G4DataVector* Cuts() const;
186  inline const G4DataVector* SubCutoff() const;
187 
188  // Set flag of fluorescence
189  inline void SetFluoFlag(G4bool val);
190 
191  // total number of models
192  inline G4int NumberOfModels() const;
193 
194 private:
195 
197  const G4MaterialCutsCouple*,
198  G4double kinEnergy,
199  G4double cutEnergy,
200  G4double minEnergy);
201 
202  // hide assignment operator
203 
206 
207 // =====================================================================
208 
209 private:
210 
214 
215  std::vector<G4VEmModel*> models;
216  std::vector<G4VEmFluctuationModel*> flucModels;
217  std::vector<const G4Region*> regions;
218  std::vector<G4int> orderOfModels;
219  std::vector<G4int> isUsed;
220 
223 
224  std::vector<G4int> idxOfRegionModels;
225  std::vector<G4RegionModels*> setOfRegionModels;
226 
228 
230 
234 
235  // may be changed in run time
238 };
239 
240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242 
244  size_t& index)
245 {
246  if(severalModels) {
247  if(nRegions > 1) {
249  }
251  }
252  return currModel;
253 }
254 
255 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
256 
258 {
259  return theCuts;
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
263 
265 {
266  return theSubCuts;
267 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
270 
272 {
273  fluoFlag = val;
274 }
275 
276 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
277 
279 {
280  return nEmModels;
281 }
282 
283 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
284 
285 inline G4double
287  const G4MaterialCutsCouple* couple,
288  G4double e,
289  G4double cut,
290  G4double emin)
291 {
292  G4double dedx = 0.0;
293  if(model && cut > emin) {
294  dedx = model->ComputeDEDX(couple,particle,e,cut);
295  if(emin > 0.0) {dedx -= model->ComputeDEDX(couple,particle,e,emin);}
296  }
297  return dedx;
298 }
299 
300 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
301 
302 #endif
303