Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
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#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"
77#include "G4EmProcessSubType.hh"
78#include "G4Region.hh"
79
80#include "G4VEmModel.hh"
82#include "G4DynamicParticle.hh"
83#include <iostream>
84
86{
87
88friend class G4EmModelManager;
89
90private:
91
92 G4RegionModels(G4int nMod, std::vector<G4int>& indx,
93 G4DataVector& lowE, const G4Region* reg);
94
96
97 inline G4int SelectIndex(G4double e) const {
98 G4int idx = 0;
99 if (nModelsForRegion>1) {
100 idx = nModelsForRegion;
101 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
102 do {--idx;} while (idx > 0 && e <= lowKineticEnergy[idx]);
103 }
104 return theListOfModelIndexes[idx];
105 };
106
107 inline G4int ModelIndex(G4int n) const {
108 return theListOfModelIndexes[n];
109 };
110
111 inline G4int NumberOfModels() const {
112 return nModelsForRegion;
113 };
114
115 inline G4double LowEdgeEnergy(G4int n) const {
116 return lowKineticEnergy[n];
117 };
118
119 inline const G4Region* Region() const {
120 return theRegion;
121 };
122
123 G4RegionModels(G4RegionModels &) = delete;
124 G4RegionModels & operator=(const G4RegionModels &right) = delete;
125
126 const G4Region* theRegion;
127 G4int nModelsForRegion;
128 G4int* theListOfModelIndexes;
129 G4double* lowKineticEnergy;
130
131};
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134
135class G4Region;
137class G4PhysicsVector;
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
143{
144public:
145
147
149
150 void Clear();
151
153 const G4ParticleDefinition* secPart,
154 G4int verb);
155
158
160 G4bool startFromNull = true,
162
164 const G4Region* r);
165
166 // Get model pointer from the model list
167 G4VEmModel* GetModel(G4int idx, G4bool ver = false) const;
168
169 // Get model pointer from the model list for a given material cuts couple
170 // no check on material cuts couple index
171 G4VEmModel* GetRegionModel(G4int idx, std::size_t index_couple);
172
173 // total number of models for material cut couples
174 // no check on material cuts couple index
175 G4int NumberOfRegionModels(std::size_t index_couple) const;
176
177 // Automatic documentation
178 void DumpModelList(std::ostream& out, G4int verb);
179
180 // Select model for given material cuts couple index
181 inline G4VEmModel* SelectModel(G4double energy, std::size_t index);
182
183 // Access to cuts
184 inline const G4DataVector* Cuts() const;
185
186 // Set flag of fluorescence
187 inline void SetFluoFlag(G4bool val);
188
189 // total number of models
190 inline G4int NumberOfModels() const;
191
192 // hide assignment operator
195
196private:
197
198 const G4ParticleDefinition* particle = nullptr;
199 const G4DataVector* theCuts = nullptr;
200 G4DataVector* theCutsNew = nullptr;
201
202 // may be changed in run time
203 G4RegionModels* currRegionModel = nullptr;
204 G4VEmModel* currModel = nullptr;
205
206 G4int nEmModels = 0;
207 G4int nRegions = 0;
208
209 G4int verboseLevel = 0;
210 G4bool severalModels = true;
211 G4bool fluoFlag = false;
212
213 std::vector<G4VEmModel*> models;
214 std::vector<G4VEmFluctuationModel*> flucModels;
215 std::vector<const G4Region*> regions;
216 std::vector<G4int> orderOfModels;
217 std::vector<G4int> isUsed;
218
219 std::vector<G4int> idxOfRegionModels;
220 std::vector<G4RegionModels*> setOfRegionModels;
221};
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
225
226inline
228{
229 if(severalModels) {
230 if(nRegions > 1) {
231 currRegionModel = setOfRegionModels[idxOfRegionModels[index]];
232 }
233 currModel = models[currRegionModel->SelectIndex(kinEnergy)];
234 }
235 return currModel;
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
241{
242 return theCuts;
243}
244
245//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
246
248{
249 fluoFlag = val;
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253
255{
256 return nEmModels;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260
261#endif
262
G4EmTableType
@ fRestricted
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4int verb)
G4VEmModel * SelectModel(G4double energy, std::size_t index)
G4EmModelManager & operator=(const G4EmModelManager &right)=delete
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
void SetFluoFlag(G4bool val)
G4VEmModel * GetRegionModel(G4int idx, std::size_t index_couple)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
const G4DataVector * Cuts() const
G4EmModelManager(G4EmModelManager &)=delete
G4int NumberOfRegionModels(std::size_t index_couple) const