Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmUtility.cc
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// Geant4 class G4EmUtility
28//
29// Author V.Ivanchenko 14.03.2022
30//
31
32#include "G4EmUtility.hh"
33#include "G4RegionStore.hh"
35#include "G4VEmProcess.hh"
36#include "G4EmParameters.hh"
37#include "G4PhysicsVector.hh"
38#include "Randomize.hh"
39#include "G4Log.hh"
40#include "G4Exp.hh"
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44static const G4double g4log10 = G4Log(10.);
45
46const G4Region*
47G4EmUtility::FindRegion(const G4String& regionName, const G4int verbose)
48{
49 const G4Region* reg = nullptr;
51 G4String r = regionName;
52 if(r == "") { r = "DefaultRegionForTheWorld"; }
53 reg = regStore->GetRegion(r, true);
54 if(nullptr == reg && verbose > 0) {
55 G4cout << "### G4EmUtility WARNING: fails to find a region <"
56 << r << G4endl;
57 } else if(verbose > 1) {
58 G4cout << "### G4EmUtility finds out G4Region <" << r << ">"
59 << G4endl;
60 }
61 return reg;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
65
67{
68 const G4Element* elm = mat->GetElement(0);
69 std::size_t nElements = mat->GetNumberOfElements();
70 if(1 < nElements) {
72 const G4double* y = mat->GetVecNbOfAtomsPerVolume();
73 for(std::size_t i=0; i<nElements; ++i) {
74 elm = mat->GetElement((G4int)i);
75 x -= y[i]*elm->GetZ();
76 if(x <= 0.0) { break; }
77 }
78 }
79 return elm;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
83
85{
86 const std::size_t ni = elm->GetNumberOfIsotopes();
87 const G4Isotope* iso = elm->GetIsotope(0);
88 if(ni > 1) {
89 const G4double* ab = elm->GetRelativeAbundanceVector();
91 for(std::size_t idx=0; idx<ni; ++idx) {
92 x -= ab[idx];
93 if (x <= 0.0) {
94 iso = elm->GetIsotope((G4int)idx);
95 break;
96 }
97 }
98 }
99 return iso;
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
103
105{
106 std::vector<G4double>* ptr = nullptr;
107 if(nullptr == p) { return ptr; }
108
109 const std::size_t n = p->length();
110 ptr = new std::vector<G4double>;
111 ptr->resize(n, DBL_MAX);
112
113 G4bool isPeak = false;
114 G4double e, ss, ee, xs;
115
116 // first loop on existing vectors
117 for (std::size_t i=0; i<n; ++i) {
118 const G4PhysicsVector* pv = (*p)[i];
119 xs = ee = 0.0;
120 if(nullptr != pv) {
121 G4int nb = (G4int)pv->GetVectorLength();
122 for (G4int j=0; j<nb; ++j) {
123 e = pv->Energy(j);
124 ss = (*pv)(j);
125 if(ss >= xs) {
126 xs = ss;
127 ee = e;
128 continue;
129 } else {
130 isPeak = true;
131 (*ptr)[i] = ee;
132 break;
133 }
134 }
135 }
136 }
137
138 // there is no peak for any material
139 if(!isPeak) {
140 delete ptr;
141 ptr = nullptr;
142 }
143 return ptr;
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
147
148std::vector<G4double>*
150 const G4ParticleDefinition* part)
151{
152 std::vector<G4double>* ptr = nullptr;
153 if (nullptr == p || nullptr == part) { return ptr; }
154
155 G4EmParameters* theParameters = G4EmParameters::Instance();
156 const G4double tmin = theParameters->MinKinEnergy();
157 const G4double tmax = theParameters->MaxKinEnergy();
158 const G4double ee = G4Log(tmax/tmin);
159 const G4double scale = theParameters->NumberOfBinsPerDecade()/g4log10;
160 G4int nbin = static_cast<G4int>(ee*scale);
161 nbin = std::max(nbin, 4);
162 G4double x = G4Exp(ee/(G4double)nbin);
163
164 const G4ProductionCutsTable* theCoupleTable=
166 std::size_t n = theCoupleTable->GetTableSize();
167 ptr = new std::vector<G4double>;
168 ptr->resize(n, DBL_MAX);
169
170 G4bool isPeak = false;
171
172 // first loop on existing vectors
173 const G4int nn = static_cast<G4int>(n);
174 for (G4int i=0; i<nn; ++i) {
175 G4double sm = 0.0;
176 G4double em = 0.0;
177 G4double e = tmin;
178 for (G4int j=0; j<=nbin; ++j) {
179 G4double sig = p->GetCrossSection(e, theCoupleTable->GetMaterialCutsCouple(i));
180 if (sig >= sm) {
181 em = e;
182 sm = sig;
183 e = (j+1 < nbin) ? e*x : tmax;
184 } else {
185 isPeak = true;
186 (*ptr)[i] = em;
187 break;
188 }
189 }
190 }
191 // there is no peak for any couple
192 if (!isPeak) {
193 delete ptr;
194 ptr = nullptr;
195 }
196 return ptr;
197}
198
199//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
200
201std::vector<G4TwoPeaksXS*>*
203{
204 std::vector<G4TwoPeaksXS*>* ptr = nullptr;
205 if(nullptr == p) { return ptr; }
206
207 const G4int n = (G4int)p->length();
208 ptr = new std::vector<G4TwoPeaksXS*>;
209 ptr->resize(n, nullptr);
210
211 G4double e, ss, xs, ee;
212 G4double e1peak, e1deep, e2peak, e2deep, e3peak;
213 G4bool isDeep = false;
214
215 // first loop on existing vectors
216 for (G4int i=0; i<n; ++i) {
217 const G4PhysicsVector* pv = (*p)[i];
218 ee = xs = 0.0;
219 e1peak = e1deep = e2peak = e2deep = e3peak = DBL_MAX;
220 if(nullptr != pv) {
221 G4int nb = (G4int)pv->GetVectorLength();
222 for (G4int j=0; j<nb; ++j) {
223 e = pv->Energy(j);
224 ss = (*pv)(j);
225 // find out 1st peak
226 if(e1peak == DBL_MAX) {
227 if(ss >= xs) {
228 xs = ss;
229 ee = e;
230 continue;
231 } else {
232 e1peak = ee;
233 }
234 }
235 // find out the deep
236 if(e1deep == DBL_MAX) {
237 if(ss <= xs) {
238 xs = ss;
239 ee = e;
240 continue;
241 } else {
242 e1deep = ee;
243 isDeep = true;
244 }
245 }
246 // find out 2nd peak
247 if(e2peak == DBL_MAX) {
248 if(ss >= xs) {
249 xs = ss;
250 ee = e;
251 continue;
252 } else {
253 e2peak = ee;
254 }
255 }
256 if(e2deep == DBL_MAX) {
257 if(ss <= xs) {
258 xs = ss;
259 ee = e;
260 continue;
261 } else {
262 e2deep = ee;
263 break;
264 }
265 }
266 // find out 3d peak
267 if(e3peak == DBL_MAX) {
268 if(ss >= xs) {
269 xs = ss;
270 ee = e;
271 continue;
272 } else {
273 e3peak = ee;
274 }
275 }
276 }
277 }
278 G4TwoPeaksXS* x = (*ptr)[i];
279 if(nullptr == x) {
280 x = new G4TwoPeaksXS();
281 (*ptr)[i] = x;
282 }
283 x->e1peak = e1peak;
284 x->e1deep = e1deep;
285 x->e2peak = e2peak;
286 x->e2deep = e2deep;
287 x->e3peak = e3peak;
288 }
289 // case of no 1st peak in all vectors
290 if(!isDeep) {
291 delete ptr;
292 ptr = nullptr;
293 return ptr;
294 }
295 // check base particles
296 if(!bld->GetBaseMaterialFlag()) { return ptr; }
297
298 auto theDensityIdx = bld->GetCoupleIndexes();
299 // second loop using base materials
300 for (G4int i=0; i<n; ++i) {
301 const G4PhysicsVector* pv = (*p)[i];
302 if (nullptr == pv) {
303 G4int j = (*theDensityIdx)[i];
304 if(j == i) { continue; }
305 G4TwoPeaksXS* x = (*ptr)[i];
306 G4TwoPeaksXS* y = (*ptr)[j];
307 if(nullptr == x) {
308 x = new G4TwoPeaksXS();
309 (*ptr)[i] = x;
310 }
311 x->e1peak = y->e1peak;
312 x->e1deep = y->e1deep;
313 x->e2peak = y->e2peak;
314 x->e2deep = y->e2deep;
315 x->e3peak = y->e3peak;
316 }
317 }
318 return ptr;
319}
320
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
322
324 const G4ParticleDefinition* part,
325 const G4DataVector& cuts,
326 const G4double elow,
327 const G4double ehigh)
328{
329 // using spline for element selectors should be investigated in details
330 // because small number of points may provide biased results
331 // large number of points requires significant increase of memory
332 G4bool spline = false;
333
335
336 G4ProductionCutsTable* theCoupleTable=
338 std::size_t numOfCouples = theCoupleTable->GetTableSize();
339
340 // prepare vector
341 auto elmSelectors = mod->GetElementSelectors();
342 if(nullptr == elmSelectors) {
343 elmSelectors = new std::vector<G4EmElementSelector*>;
344 }
345 std::size_t nSelectors = elmSelectors->size();
346 if(numOfCouples > nSelectors) {
347 for(std::size_t i=nSelectors; i<numOfCouples; ++i) {
348 elmSelectors->push_back(nullptr);
349 }
350 nSelectors = numOfCouples;
351 }
352
353 // initialise vector
354 for(std::size_t i=0; i<numOfCouples; ++i) {
355
356 // no need in element selectors for infinite cuts
357 if(cuts[i] == DBL_MAX) { continue; }
358
359 auto couple = theCoupleTable->GetMaterialCutsCouple((G4int)i);
360 auto mat = couple->GetMaterial();
361 mod->SetCurrentCouple(couple);
362
363 // selector already exist then delete
364 delete (*elmSelectors)[i];
365
366 G4double emin = std::max(elow, mod->MinPrimaryEnergy(mat, part, cuts[i]));
367 G4double emax = std::max(ehigh, 10*emin);
368 static const G4double invlog106 = 1.0/(6*G4Log(10.));
369 G4int nbins = G4lrint(nbinsPerDec*G4Log(emax/emin)*invlog106);
370 nbins = std::max(nbins, 3);
371
372 (*elmSelectors)[i] = new G4EmElementSelector(mod,mat,nbins,
373 emin,emax,spline);
374 ((*elmSelectors)[i])->Initialise(part, cuts[i]);
375 /*
376 G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
377 << " " << part->GetParticleName()
378 << " for " << mod->GetName() << " cut= " << cuts[i]
379 << " " << (*elmSelectors)[i] << G4endl;
380 ((*elmSelectors)[i])->Dump(part);
381 */
382 }
383 mod->SetElementSelectors(elmSelectors);
384}
385
386//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
G4double GetZ() const
Definition G4Element.hh:119
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
static const G4Element * SampleRandomElement(const G4Material *)
static std::vector< G4TwoPeaksXS * > * FillPeaksStructure(G4PhysicsTable *, G4LossTableBuilder *)
static const G4Isotope * SampleRandomIsotope(const G4Element *)
static void InitialiseElementSelectors(G4VEmModel *, const G4ParticleDefinition *, const G4DataVector &cuts, const G4double emin, const G4double emax)
static const G4Region * FindRegion(const G4String &regionName, const G4int verbose=0)
static std::vector< G4double > * FindCrossSectionMax(G4PhysicsTable *)
const std::vector< G4int > * GetCoupleIndexes() const
const G4Material * GetMaterial() const
const G4Element * GetElement(G4int iel) const
G4double GetTotNbOfElectPerVolume() const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
std::size_t length() const
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
virtual G4double GetCrossSection(const G4double, const G4MaterialCutsCouple *)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()
void SetCurrentCouple(const G4MaterialCutsCouple *)
G4double e1peak
G4double e3peak
G4double e2deep
G4double e1deep
G4double e2peak
int G4lrint(double ad)
Definition templates.hh:134
#define DBL_MAX
Definition templates.hh:62