Geant4-11
G4VEmModel.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//
28// GEANT4 Class file
29//
30//
31// File name: G4VEmModel
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 25.07.2005
36//
37// Modifications:
38// 25.10.2005 Set default highLimit=100.TeV (V.Ivanchenko)
39// 06.02.2006 add method ComputeMeanFreePath() (mma)
40// 16.02.2009 Move implementations of virtual methods to source (VI)
41//
42//
43// Class Description:
44//
45// Abstract interface to energy loss models
46
47// -------------------------------------------------------------------
48//
49
50#include "G4VEmModel.hh"
51#include "G4ElementData.hh"
52#include "G4LossTableManager.hh"
53#include "G4LossTableBuilder.hh"
57#include "G4EmParameters.hh"
58#include "G4SystemOfUnits.hh"
59#include "G4Log.hh"
60#include "Randomize.hh"
61#include <iostream>
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
67 inveplus(1.0/CLHEP::eplus),
68 lowLimit(0.1*CLHEP::keV),
69 highLimit(100.0*CLHEP::TeV),
70 polarAngleLimit(CLHEP::pi),
71 name(nam)
72{
73 xsec.resize(nsec);
75 fEmManager->Register(this);
76
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85{
87 for(G4int i=0; i<nSelectors; ++i) {
88 delete (*elmSelectors)[i];
89 }
90 delete elmSelectors;
91 }
92 delete anglModel;
93
94 if(localTable && xSectionTable != nullptr) {
96 delete xSectionTable;
97 xSectionTable = nullptr;
98 }
99 if(isMaster && fElementData != nullptr) {
100 delete fElementData;
101 fElementData = nullptr;
102 }
103 fEmManager->DeRegister(this);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109{
110 G4ParticleChangeForLoss* p = nullptr;
111 if (pParticleChange != nullptr) {
112 p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
113 } else {
114 p = new G4ParticleChangeForLoss();
115 pParticleChange = p;
116 }
117 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
118 return p;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
124{
125 G4ParticleChangeForGamma* p = nullptr;
126 if (pParticleChange != nullptr) {
127 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
128 } else {
129 p = new G4ParticleChangeForGamma();
130 pParticleChange = p;
131 }
132 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
133 return p;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139 const G4DataVector& cuts)
140{
141 // using spline for element selectors should be investigated in details
142 // because small number of points may provide biased results
143 // large number of points requires significant increase of memory
144 G4bool spline = false;
145
146 //G4cout << "IES: for " << GetName() << " Emin(MeV)= " << lowLimit/MeV
147 // << " Emax(MeV)= " << highLimit/MeV << G4endl;
148
149 // two times less bins because probability functon is normalized
150 // so correspondingly is more smooth
151 if(highLimit <= lowLimit) { return; }
152
154
155 G4ProductionCutsTable* theCoupleTable=
157 G4int numOfCouples = theCoupleTable->GetTableSize();
158
159 // prepare vector
160 if(!elmSelectors) {
161 elmSelectors = new std::vector<G4EmElementSelector*>;
162 }
163 if(numOfCouples > nSelectors) {
164 for(G4int i=nSelectors; i<numOfCouples; ++i) {
165 elmSelectors->push_back(nullptr);
166 }
167 nSelectors = numOfCouples;
168 }
169
170 // initialise vector
171 for(G4int i=0; i<numOfCouples; ++i) {
172
173 // no need in element selectors for infinite cuts
174 if(cuts[i] == DBL_MAX) { continue; }
175
176 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
177 auto material = couple->GetMaterial();
178 SetCurrentCouple(couple);
179
180 // selector already exist then delete
181 delete (*elmSelectors)[i];
182
183 G4double emin = std::max(lowLimit, MinPrimaryEnergy(material, part, cuts[i]));
184 G4double emax = std::max(highLimit, 10*emin);
185 static const G4double invlog106 = 1.0/(6*G4Log(10.));
186 G4int nbins = (G4int)(nbinsPerDec*G4Log(emax/emin)*invlog106);
187 nbins = std::max(nbins, 3);
188
189 (*elmSelectors)[i] = new G4EmElementSelector(this,material,nbins,
190 emin,emax,spline);
191 ((*elmSelectors)[i])->Initialise(part, cuts[i]);
192 /*
193 G4cout << "G4VEmModel::InitialiseElmSelectors i= " << i
194 << " " << part->GetParticleName()
195 << " for " << GetName() << " cut= " << cuts[i]
196 << " " << (*elmSelectors)[i] << G4endl;
197 ((*elmSelectors)[i])->Dump(part);
198 */
199 }
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203
205{}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
208
210 const G4Material* material)
211{
212 if(material != nullptr) {
213 size_t n = material->GetNumberOfElements();
214 for(size_t i=0; i<n; ++i) {
215 G4int Z = material->GetElement(i)->GetZasInt();
216 InitialiseForElement(part, Z);
217 }
218 }
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222
224{}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227
231{
232 return 0.0;
233}
234
235//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
236
238 const G4ParticleDefinition* p,
239 G4double ekin,
240 G4double emin,
242{
243 SetupForMaterial(p, mat, ekin);
244 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
245 G4int nelm = mat->GetNumberOfElements();
246 if(nelm > nsec) {
247 xsec.resize(nelm);
248 nsec = nelm;
249 }
250 G4double cross = 0.0;
251 for (G4int i=0; i<nelm; ++i) {
252 cross += theAtomNumDensityVector[i]*
253 ComputeCrossSectionPerAtom(p,mat->GetElement(i),ekin,emin,emax);
254 xsec[i] = cross;
255 }
256 return cross;
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260
263 G4double)
264{
265 return 0.0;
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
269
271{}
272
273//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
274
276 const G4ParticleDefinition* pd,
277 G4double kinEnergy,
278 G4double tcut,
279 G4double tmax)
280{
281 size_t n = mat->GetNumberOfElements();
282 fCurrentElement = mat->GetElement(0);
283 if (n > 1) {
284 const G4double x = G4UniformRand()*
285 G4VEmModel::CrossSectionPerVolume(mat,pd,kinEnergy,tcut,tmax);
286 for(size_t i=0; i<n; ++i) {
287 if (x <= xsec[i]) {
288 fCurrentElement = mat->GetElement(i);
289 break;
290 }
291 }
292 }
293 return fCurrentElement;
294}
295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
296
298{
299 // this algorith assumes that cross section is proportional to
300 // number electrons multiplied by number of atoms
301 const size_t nn = mat->GetNumberOfElements();
302 fCurrentElement = mat->GetElement(0);
303 if(1 < nn) {
304 const G4double* at = mat->GetVecNbOfAtomsPerVolume();
306 for(size_t i=0; i<nn; ++i) {
307 tot -= at[i];
308 if(tot <= 0.0) {
309 fCurrentElement = mat->GetElement(i);
310 break;
311 }
312 }
313 }
314 return fCurrentElement->GetZasInt();
315}
316
317//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
318
320{
322 const size_t ni = elm->GetNumberOfIsotopes();
323 fCurrentIsotope = elm->GetIsotope(0);
324 size_t idx = 0;
325 if(ni > 1) {
326 const G4double* ab = elm->GetRelativeAbundanceVector();
328 for(; idx<ni; ++idx) {
329 x -= ab[idx];
330 if (x <= 0.0) {
331 fCurrentIsotope = elm->GetIsotope(idx);
332 break;
333 }
334 }
335 }
336 return fCurrentIsotope->GetN();
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
340
344{
345 return 0.0;
346}
347
348//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
349
352 G4int, G4int,
354{
355 return 0.0;
356}
357
358//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
359
361{}
362
363//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
364
366 G4int& numberOfRecoil)
367{
368 numberOfTriplets = 0;
369 numberOfRecoil = 0;
370}
371
372//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
373
375{
377 track.GetMaterial(), track.GetKineticEnergy());
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
381
383 const G4Material*, G4double)
384{
385 const G4double q = p->GetPDGCharge()*inveplus;
386 return q*q;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
390
392 const G4Material*, G4double)
393{
394 return p->GetPDGCharge();
395}
396
397//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
398
400 const G4DynamicParticle*,
401 const G4double&,G4double&)
402{}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
405
407 const G4ParticleDefinition* p, G4double e)
408{
409 SetCurrentCouple(couple);
411}
412
413//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
414
417 G4double)
418{
419 return 0.0;
420}
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
423
426{
427 return 0.0;
428}
429
430//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
431
433 G4double kineticEnergy)
434{
435 return kineticEnergy;
436}
437
438//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
439
441 const G4Material*, G4double)
442{}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
445
446void
448{
449 if(p != nullptr && pParticleChange != p) { pParticleChange = p; }
450 if(flucModel != f) { flucModel = f; }
451}
452
453//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
454
456{
457 if(p != xSectionTable) {
458 if(xSectionTable != nullptr && localTable) {
460 delete xSectionTable;
461 }
462 xSectionTable = p;
463 }
464 localTable = isLocal;
465}
466
467//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
468
469void G4VEmModel::ModelDescription(std::ostream& outFile) const
470{
471 outFile << "The description for this model has not been written yet.\n";
472}
473
474//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const G4double emax
static const G4double ab
const G4double inveplus
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double TeV
Definition: G4SIunits.hh:204
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
G4int GetN() const
Definition: G4Isotope.hh:93
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
G4double GetPDGCharge() const
void clearAndDestroy()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:365
G4double highLimit
Definition: G4VEmModel.hh:437
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:455
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:415
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:426
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:428
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391
std::vector< G4double > xsec
Definition: G4VEmModel.hh:466
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:209
G4double inveplus
Definition: G4VEmModel.hh:431
G4bool localTable
Definition: G4VEmModel.hh:459
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:123
G4bool isMaster
Definition: G4VEmModel.hh:457
G4double lowLimit
Definition: G4VEmModel.hh:436
const G4Element * fCurrentElement
Definition: G4VEmModel.hh:417
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:472
G4int nSelectors
Definition: G4VEmModel.hh:443
G4bool localElmSelectors
Definition: G4VEmModel.hh:460
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:223
G4int nsec
Definition: G4VEmModel.hh:444
G4int SelectRandomAtomNumber(const G4Material *)
Definition: G4VEmModel.cc:297
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:406
G4VEmAngularDistribution * anglModel
Definition: G4VEmModel.hh:414
void SetCurrentElement(const G4Element *)
Definition: G4VEmModel.hh:497
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:425
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:580
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:427
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:360
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:469
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.cc:319
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:440
G4double pFactor
Definition: G4VEmModel.hh:432
std::vector< G4EmElementSelector * > * elmSelectors
Definition: G4VEmModel.hh:419
G4ElementData * fElementData
Definition: G4VEmModel.hh:424
const G4Isotope * fCurrentIsotope
Definition: G4VEmModel.hh:418
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:66
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:351
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:399
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:429
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:382
G4VEmFluctuationModel * flucModel
Definition: G4VEmModel.hh:413
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:261
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:228
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:138
G4LossTableManager * fEmManager
Definition: G4VEmModel.hh:420
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4VEmModel * fTripletModel
Definition: G4VEmModel.hh:415
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:84
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:270
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:204
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:424
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:432
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:374
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:108
Definition: DoubConv.h:17
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62