Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VMscModel.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4VMscModel
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 07.03.2008
38//
39// Modifications:
40//
41//
42// Class Description:
43//
44// General interface to msc models
45
46// -------------------------------------------------------------------
47//
48
49#include "G4VMscModel.hh"
50#include "G4SystemOfUnits.hh"
53#include "G4LossTableBuilder.hh"
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
59 G4VEmModel(nam),
60 safetyHelper(0),
61 ionisation(0),
62 facrange(0.04),
63 facgeom(2.5),
64 facsafety(0.3),
65 skin(1.0),
66 dtrl(0.05),
67 lambdalimit(mm),
68 geomMin(1.e-6*CLHEP::mm),
69 geomMax(1.e50*CLHEP::mm),
70 steppingAlgorithm(fUseSafety),
71 samplez(false),
72 latDisplasment(true)
73{
74 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
75 localrange = DBL_MAX;
76 localtkin = 0.0;
78 currentPart = 0;
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
84{}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
90{
91 if(!safetyHelper) {
94 safetyHelper->InitialiseHelper();
95 }
96 G4ParticleChangeForMSC* change = 0;
97 if (pParticleChange) {
98 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
99 } else {
100 change = new G4ParticleChangeForMSC();
101 }
102 if(p) {
103
104 // table is never built for GenericIon
105 if(p->GetParticleName() == "GenericIon") {
106 if(xSectionTable) {
108 delete xSectionTable;
109 xSectionTable = 0;
110 }
111
112 // table is always built for low mass particles
113 } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) {
116 emin = std::max(emin, man->MinKinEnergy());
117 emax = std::min(emax, man->MaxKinEnergy());
118 G4LossTableBuilder* builder = man->GetTableBuilder();
120 emin, emax, true);
122 theDensityIdx = builder->GetCoupleIndexes();
123 }
124 }
125 return change;
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
129
132{
133 return fDisplacement;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
139{
140 return DBL_MAX;
141}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
146{
147 return truePathLength;
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
153{
154 return geomPathLength;
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
159void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
161 const G4DynamicParticle*,
163{}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fUseSafety
double G4double
Definition: G4Types.hh:64
const std::vector< G4double > * GetDensityFactors()
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
const std::vector< G4int > * GetCoupleIndexes()
G4double MinKinEnergy() const
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
G4double MaxKinEnergy() const
const G4String & GetParticleName() const
void clearAndDestroy()
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:352
const std::vector< G4double > * theDensityFactor
Definition: G4VEmModel.hh:353
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:351
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:536
const std::vector< G4int > * theDensityIdx
Definition: G4VEmModel.hh:354
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:543
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:578
virtual G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
Definition: G4VMscModel.cc:131
virtual G4double ComputeGeomPathLength(G4double truePathLength)
Definition: G4VMscModel.cc:145
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:58
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=0)
Definition: G4VMscModel.cc:89
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax)
Definition: G4VMscModel.cc:159
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:152
virtual ~G4VMscModel()
Definition: G4VMscModel.cc:83
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:138
G4ThreeVector fDisplacement
Definition: G4VMscModel.hh:185
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:83