Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "G4EmUtility.hh"
60#include "G4Log.hh"
61#include "Randomize.hh"
62#include <iostream>
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
68 inveplus(1.0/CLHEP::eplus),
69 lowLimit(0.1*CLHEP::keV),
70 highLimit(100.0*CLHEP::TeV),
71 polarAngleLimit(CLHEP::pi),
72 name(nam)
73{
74 xsec.resize(nsec);
75 fEmManager = G4LossTableManager::Instance();
76 fEmManager->Register(this);
77
78 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86{
87 if(localElmSelectors) {
88 for(G4int i=0; i<nSelectors; ++i) {
89 delete (*elmSelectors)[i];
90 }
91 delete elmSelectors;
92 }
93 delete anglModel;
94
95 if(localTable && xSectionTable != nullptr) {
97 delete xSectionTable;
98 xSectionTable = nullptr;
99 }
100 if(isMaster && fElementData != nullptr) {
101 delete fElementData;
102 fElementData = nullptr;
103 }
104 fEmManager->DeRegister(this);
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
110{
111 G4ParticleChangeForLoss* p = nullptr;
112 if (pParticleChange != nullptr) {
113 p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
114 } else {
115 p = new G4ParticleChangeForLoss();
116 pParticleChange = p;
117 }
118 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
119 return p;
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123
125{
126 G4ParticleChangeForGamma* p = nullptr;
127 if (pParticleChange != nullptr) {
128 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
129 } else {
130 p = new G4ParticleChangeForGamma();
131 pParticleChange = p;
132 }
133 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
134 return p;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
140 const G4DataVector& cuts)
141{
142 if(highLimit <= lowLimit) { return; }
143 G4EmUtility::InitialiseElementSelectors(this,part,cuts,lowLimit,highLimit);
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
149{}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
154 const G4Material* material)
155{
156 if(material != nullptr) {
157 G4int n = (G4int)material->GetNumberOfElements();
158 for(G4int i=0; i<n; ++i) {
159 G4int Z = material->GetElement(i)->GetZasInt();
160 InitialiseForElement(part, Z);
161 }
162 }
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166
168{}
169
170//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171
175{
176 return 0.0;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
180
182 const G4ParticleDefinition* p,
183 G4double ekin,
184 G4double emin,
185 G4double emax)
186{
187 SetupForMaterial(p, mat, ekin);
188 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
189 G4int nelm = (G4int)mat->GetNumberOfElements();
190 if(nelm > nsec) {
191 xsec.resize(nelm);
192 nsec = nelm;
193 }
194 G4double cross = 0.0;
195 for (G4int i=0; i<nelm; ++i) {
196 cross += theAtomNumDensityVector[i]*
197 ComputeCrossSectionPerAtom(p,mat->GetElement(i),ekin,emin,emax);
198 xsec[i] = cross;
199 }
200 return cross;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204
207 G4double)
208{
209 return 0.0;
210}
211
212//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213
215{}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
218
220 const G4ParticleDefinition* pd,
221 G4double kinEnergy,
222 G4double tcut,
223 G4double tmax)
224{
225 G4int n = (G4int)mat->GetNumberOfElements();
226 fCurrentElement = mat->GetElement(0);
227 if (n > 1) {
228 const G4double x = G4UniformRand()*
229 G4VEmModel::CrossSectionPerVolume(mat,pd,kinEnergy,tcut,tmax);
230 for(G4int i=0; i<n; ++i) {
231 if (x <= xsec[i]) {
232 fCurrentElement = mat->GetElement(i);
233 break;
234 }
235 }
236 }
237 return fCurrentElement;
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241
243{
244 const G4Element* elm = fCurrentElement;
245 if(nullptr == elm && nullptr != mat) {
247 }
248 return elm;
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252
254{
255 const G4Element* elm = GetCurrentElement(mat);
256 return (nullptr == elm) ? 0 : elm->GetZasInt();
257}
258
259//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
260
262{
263 const G4Isotope* iso = nullptr;
264 const G4Element* el = elm;
265 if(nullptr == el && nullptr != fCurrentCouple) {
266 el = GetCurrentElement(fCurrentCouple->GetMaterial());
267 }
268 if(nullptr != el) {
270 }
271 return iso;
272}
273
274//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
275
277{
278 auto iso = GetCurrentIsotope(elm);
279 return (nullptr != iso) ? iso->GetN() : 0;
280}
281
282//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
283
287{
288 return 0.0;
289}
290
291//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
292
295 G4int, G4int,
297{
298 return 0.0;
299}
300
301//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
302
304{}
305
306//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
307
309 G4int& numberOfRecoil)
310{
311 numberOfTriplets = 0;
312 numberOfRecoil = 0;
313}
314
315//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
316
318{
320 track.GetMaterial(), track.GetKineticEnergy());
321}
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
324
326 const G4Material*, G4double)
327{
328 const G4double q = p->GetPDGCharge()*inveplus;
329 return q*q;
330}
331
332//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
333
335 const G4Material*, G4double)
336{
337 return p->GetPDGCharge();
338}
339
340//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
341
343 const G4DynamicParticle*,
344 const G4double&,G4double&)
345{}
346
347//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
348
350 const G4ParticleDefinition* p, G4double e)
351{
352 SetCurrentCouple(couple);
354}
355
356//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
357
360 G4double)
361{
362 return 0.0;
363}
364
365//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
366
369{
370 return 0.0;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
374
376 G4double kineticEnergy)
377{
378 return kineticEnergy;
379}
380
381//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
382
384 const G4Material*, G4double)
385{}
386
387//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
388
389void
391{
392 if(p != nullptr && pParticleChange != p) { pParticleChange = p; }
393 if(flucModel != f) { flucModel = f; }
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397
399{
400 if(p != xSectionTable) {
401 if(xSectionTable != nullptr && localTable) {
403 delete xSectionTable;
404 }
405 xSectionTable = p;
406 }
407 localTable = isLocal;
408}
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
411
412void G4VEmModel::ModelDescription(std::ostream& outFile) const
413{
414 outFile << "The description for this model has not been written yet.\n";
415}
416
417//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
const G4double inveplus
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
G4int GetZasInt() const
Definition: G4Element.hh:132
static const G4Element * SampleRandomElement(const G4Material *)
Definition: G4EmUtility.cc:66
static const G4Isotope * SampleRandomIsotope(const G4Element *)
Definition: G4EmUtility.cc:84
static void InitialiseElementSelectors(G4VEmModel *, const G4ParticleDefinition *, const G4DataVector &cuts, const G4double emin, const G4double emax)
Definition: G4EmUtility.cc:332
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)
const G4Material * GetMaterial() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
G4double GetPDGCharge() const
void clearAndDestroy()
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:308
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:398
G4int SelectIsotopeNumber(const G4Element *) const
Definition: G4VEmModel.cc:276
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:358
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:422
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:284
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:424
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:334
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
Definition: G4VEmModel.cc:153
G4double inveplus
Definition: G4VEmModel.hh:427
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4int SelectRandomAtomNumber(const G4Material *) const
Definition: G4VEmModel.cc:253
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:181
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:390
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
Definition: G4VEmModel.cc:167
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:421
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:561
const G4Material * pBaseMaterial
Definition: G4VEmModel.hh:423
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
Definition: G4VEmModel.cc:261
virtual void DefineForRegion(const G4Region *)
Definition: G4VEmModel.cc:303
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4VEmModel.cc:412
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:383
G4double pFactor
Definition: G4VEmModel.hh:428
G4ElementData * fElementData
Definition: G4VEmModel.hh:420
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:67
virtual G4double ComputeCrossSectionPerShell(const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:294
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:342
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:425
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:325
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
Definition: G4VEmModel.cc:242
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:205
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
Definition: G4VEmModel.cc:172
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
virtual ~G4VEmModel()
Definition: G4VEmModel.cc:85
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:214
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:148
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
Definition: G4VEmModel.cc:367
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:375
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:317
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:62