Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 isMaster = fEmManager->IsMaster();
78
79 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87{
88 if(localElmSelectors) {
89 for(G4int i=0; i<nSelectors; ++i) {
90 delete (*elmSelectors)[i];
91 }
92 delete elmSelectors;
93 }
94 delete anglModel;
95
96 if(localTable && xSectionTable != nullptr) {
97 xSectionTable->clearAndDestroy();
98 delete xSectionTable;
99 xSectionTable = nullptr;
100 }
101 fEmManager->DeRegister(this);
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105
107{
108 G4ParticleChangeForLoss* p = nullptr;
109 if (pParticleChange != nullptr) {
110 p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
111 } else {
112 p = new G4ParticleChangeForLoss();
113 pParticleChange = p;
114 }
115 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
116 return p;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
122{
123 G4ParticleChangeForGamma* p = nullptr;
124 if (pParticleChange != nullptr) {
125 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
126 } else {
127 p = new G4ParticleChangeForGamma();
128 pParticleChange = p;
129 }
130 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
131 return p;
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
137 const G4DataVector& cuts)
138{
139 if(highLimit <= lowLimit) { return; }
140 G4EmUtility::InitialiseElementSelectors(this,part,cuts,lowLimit,highLimit);
141 localElmSelectors = true;
142}
143
144//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
145
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
152 const G4Material* material)
153{
154 if(material != nullptr) {
155 G4int n = (G4int)material->GetNumberOfElements();
156 for(G4int i=0; i<n; ++i) {
157 G4int Z = material->GetElement(i)->GetZasInt();
158 InitialiseForElement(part, Z);
159 }
160 }
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
169
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
178
180 const G4ParticleDefinition* p,
181 G4double ekin,
182 G4double emin,
183 G4double emax)
184{
185 SetupForMaterial(p, mat, ekin);
186 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
187 G4int nelm = (G4int)mat->GetNumberOfElements();
188 if(nelm > nsec) {
189 xsec.resize(nelm);
190 nsec = nelm;
191 }
192 G4double cross = 0.0;
193 for (G4int i=0; i<nelm; ++i) {
194 cross += theAtomNumDensityVector[i]*
195 ComputeCrossSectionPerAtom(p,mat->GetElement(i),ekin,emin,emax);
196 xsec[i] = cross;
197 }
198 return cross;
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
211
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216
218 const G4ParticleDefinition* pd,
219 G4double kinEnergy,
220 G4double tcut,
221 G4double tmax)
222{
223 G4int n = (G4int)mat->GetNumberOfElements();
224 fCurrentElement = mat->GetElement(0);
225 if (n > 1) {
226 const G4double x = G4UniformRand()*
227 G4VEmModel::CrossSectionPerVolume(mat,pd,kinEnergy,tcut,tmax);
228 for(G4int i=0; i<n; ++i) {
229 if (x <= xsec[i]) {
230 fCurrentElement = mat->GetElement(i);
231 break;
232 }
233 }
234 }
235 return fCurrentElement;
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
241{
242 const G4Element* elm = fCurrentElement;
243 if(nullptr == elm && nullptr != mat) {
245 }
246 return elm;
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
250
252{
253 const G4Element* elm = GetCurrentElement(mat);
254 return (nullptr == elm) ? 0 : elm->GetZasInt();
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
258
260{
261 const G4Isotope* iso = nullptr;
262 const G4Element* el = elm;
263 if(nullptr == el && nullptr != fCurrentCouple) {
264 el = GetCurrentElement(fCurrentCouple->GetMaterial());
265 }
266 if(nullptr != el) {
268 }
269 return iso;
270}
271
272//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
273
275{
276 auto iso = GetCurrentIsotope(elm);
277 return (nullptr != iso) ? iso->GetN() : 0;
278}
279
280//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
281
288
289//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
290
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
300
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305
307 G4int& numberOfRecoil)
308{
309 numberOfTriplets = 0;
310 numberOfRecoil = 0;
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
314
320
321//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
322
324 const G4Material*, G4double)
325{
326 const G4double q = p->GetPDGCharge()*inveplus;
327 return q*q;
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
331
337
338//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
339
344
345//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
346
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
355
358 G4double)
359{
360 return 0.0;
361}
362
363//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
364
370
371//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
372
374 G4double kineticEnergy)
375{
376 return kineticEnergy;
377}
378
379//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
380
382 const G4Material* mat, G4double ekin)
383{
384 GetChargeSquareRatio(p, mat, ekin);
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 xSectionTable = p;
401 localTable = isLocal;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
405
407{
408 if (G4EmParameters::Instance()->Verbose() > 0) {
410 ed << "The obsolete method SetLPMFlag(..) of the model class " << GetName()
411 << " is called. Please, use G4EmParameters::Instance()->SetLPM(..)"
412 << " instead";
413 G4Exception("G4VEmModel::SetLPMFlag", "em0001", JustWarning, ed);
414 }
415}
416
417//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
418
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
423
424void G4VEmModel::ModelDescription(std::ostream& outFile) const
425{
426 outFile << "The description for this model has not been written yet.\n";
427}
428
429//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
G4int GetZasInt() const
Definition G4Element.hh:120
static G4EmParameters * Instance()
static const G4Element * SampleRandomElement(const G4Material *)
static const G4Isotope * SampleRandomIsotope(const G4Element *)
static void InitialiseElementSelectors(G4VEmModel *, const G4ParticleDefinition *, const G4DataVector &cuts, const G4double emin, const G4double emax)
static const std::vector< G4double > * GetDensityFactors()
static const std::vector< G4int > * GetCoupleIndexes()
static G4LossTableManager * Instance()
const G4Element * GetElement(G4int iel) const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
void SetLPMFlag(G4bool)
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
void SetMasterThread(G4bool)
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
G4int SelectIsotopeNumber(const G4Element *) const
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
G4PhysicsTable * xSectionTable
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const std::vector< G4double > * theDensityFactor
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
virtual void InitialiseForMaterial(const G4ParticleDefinition *, const G4Material *)
G4double inveplus
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4int SelectRandomAtomNumber(const G4Material *) const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
G4VParticleChange * pParticleChange
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const G4Material * pBaseMaterial
const G4Isotope * GetCurrentIsotope(const G4Element *elm=nullptr) const
virtual void DefineForRegion(const G4Region *)
virtual void ModelDescription(std::ostream &outFile) const
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double pFactor
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)
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
const std::vector< G4int > * theDensityIdx
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4String & GetName() const
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
virtual G4double GetPartialCrossSection(const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
virtual ~G4VEmModel()
Definition G4VEmModel.cc:86
virtual void StartTracking(G4Track *)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double ChargeSquareRatio(const G4Track &)
G4ParticleChangeForLoss * GetParticleChangeForLoss()
#define DBL_MAX
Definition templates.hh:62