Geant4 11.2.2
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
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 fEmManager->DeRegister(this);
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104
106{
107 G4ParticleChangeForLoss* p = nullptr;
108 if (pParticleChange != nullptr) {
109 p = static_cast<G4ParticleChangeForLoss*>(pParticleChange);
110 } else {
111 p = new G4ParticleChangeForLoss();
112 pParticleChange = p;
113 }
114 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
115 return p;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119
121{
122 G4ParticleChangeForGamma* p = nullptr;
123 if (pParticleChange != nullptr) {
124 p = static_cast<G4ParticleChangeForGamma*>(pParticleChange);
125 } else {
126 p = new G4ParticleChangeForGamma();
127 pParticleChange = p;
128 }
129 if(fTripletModel != nullptr) { fTripletModel->SetParticleChange(p); }
130 return p;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
136 const G4DataVector& cuts)
137{
138 if(highLimit <= lowLimit) { return; }
139 G4EmUtility::InitialiseElementSelectors(this,part,cuts,lowLimit,highLimit);
140 localElmSelectors = true;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
151 const G4Material* material)
152{
153 if(material != nullptr) {
154 G4int n = (G4int)material->GetNumberOfElements();
155 for(G4int i=0; i<n; ++i) {
156 G4int Z = material->GetElement(i)->GetZasInt();
157 InitialiseForElement(part, Z);
158 }
159 }
160}
161
162//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
168
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177
179 const G4ParticleDefinition* p,
180 G4double ekin,
181 G4double emin,
182 G4double emax)
183{
184 SetupForMaterial(p, mat, ekin);
185 const G4double* theAtomNumDensityVector = mat->GetVecNbOfAtomsPerVolume();
186 G4int nelm = (G4int)mat->GetNumberOfElements();
187 if(nelm > nsec) {
188 xsec.resize(nelm);
189 nsec = nelm;
190 }
191 G4double cross = 0.0;
192 for (G4int i=0; i<nelm; ++i) {
193 cross += theAtomNumDensityVector[i]*
194 ComputeCrossSectionPerAtom(p,mat->GetElement(i),ekin,emin,emax);
195 xsec[i] = cross;
196 }
197 return cross;
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
201
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215
217 const G4ParticleDefinition* pd,
218 G4double kinEnergy,
219 G4double tcut,
220 G4double tmax)
221{
222 G4int n = (G4int)mat->GetNumberOfElements();
223 fCurrentElement = mat->GetElement(0);
224 if (n > 1) {
225 const G4double x = G4UniformRand()*
226 G4VEmModel::CrossSectionPerVolume(mat,pd,kinEnergy,tcut,tmax);
227 for(G4int i=0; i<n; ++i) {
228 if (x <= xsec[i]) {
229 fCurrentElement = mat->GetElement(i);
230 break;
231 }
232 }
233 }
234 return fCurrentElement;
235}
236
237//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
238
240{
241 const G4Element* elm = fCurrentElement;
242 if(nullptr == elm && nullptr != mat) {
244 }
245 return elm;
246}
247
248//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249
251{
252 const G4Element* elm = GetCurrentElement(mat);
253 return (nullptr == elm) ? 0 : elm->GetZasInt();
254}
255
256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
257
259{
260 const G4Isotope* iso = nullptr;
261 const G4Element* el = elm;
262 if(nullptr == el && nullptr != fCurrentCouple) {
263 el = GetCurrentElement(fCurrentCouple->GetMaterial());
264 }
265 if(nullptr != el) {
267 }
268 return iso;
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
272
274{
275 auto iso = GetCurrentIsotope(elm);
276 return (nullptr != iso) ? iso->GetN() : 0;
277}
278
279//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
280
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289
297
298//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
299
302
303//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
304
306 G4int& numberOfRecoil)
307{
308 numberOfTriplets = 0;
309 numberOfRecoil = 0;
310}
311
312//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
313
319
320//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
321
323 const G4Material*, G4double)
324{
325 const G4double q = p->GetPDGCharge()*inveplus;
326 return q*q;
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
330
336
337//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
338
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
354
357 G4double)
358{
359 return 0.0;
360}
361
362//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
363
369
370//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
371
373 G4double kineticEnergy)
374{
375 return kineticEnergy;
376}
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
379
381 const G4Material* mat, G4double ekin)
382{
383 GetChargeSquareRatio(p, mat, ekin);
384}
385
386//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
387
388void
390{
391 if(p != nullptr && pParticleChange != p) { pParticleChange = p; }
392 if(flucModel != f) { flucModel = f; }
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
396
398{
399 if(p != xSectionTable) {
400 if(xSectionTable != nullptr && localTable) {
402 delete xSectionTable;
403 }
404 xSectionTable = p;
405 }
406 localTable = isLocal;
407}
408
409//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
410
412{
413 if (G4EmParameters::Instance()->Verbose() > 0) {
415 ed << "The obsolete method SetLPMFlag(..) of the model class " << GetName()
416 << " is called. Please, use G4EmParameters::Instance()->SetLPM(..)"
417 << " instead";
418 G4Exception("G4VEmModel::SetLPMFlag", "em0001", JustWarning, ed);
419 }
420}
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)
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
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
void clearAndDestroy()
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
G4double GetKineticEnergy() const
void SetLPMFlag(G4bool)
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
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:85
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