Geant4 10.7.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4VMscModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 07.03.2008
37//
38// Modifications:
39//
40//
41// Class Description:
42//
43// General interface to msc models
44
45// -------------------------------------------------------------------
46//
47
48#include "G4VMscModel.hh"
49#include "G4SystemOfUnits.hh"
52#include "G4LossTableManager.hh"
53#include "G4LossTableBuilder.hh"
54#include "G4EmParameters.hh"
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
60 G4VEmModel(nam),
61 safetyHelper(nullptr),
62 ionisation(nullptr),
63 facrange(0.04),
64 facgeom(2.5),
65 facsafety(0.6),
66 skin(1.0),
67 dtrl(0.05),
68 lambdalimit(1.*CLHEP::mm),
69 geomMin(1.e-6*CLHEP::mm),
70 geomMax(1.e50*CLHEP::mm),
71 fDisplacement(0.,0.,0.),
72 steppingAlgorithm(fUseSafety),
73 samplez(false),
74 latDisplasment(true)
75{
76 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
77 localrange = DBL_MAX;
78 localtkin = 0.0;
79 currentPart = nullptr;
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86{}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
92{
93 // recomputed for each new run
94 if(!safetyHelper) {
97 safetyHelper->InitialiseHelper();
98 }
99 G4ParticleChangeForMSC* change = nullptr;
100 if (pParticleChange) {
101 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
102 } else {
103 change = new G4ParticleChangeForMSC();
104 }
105 if(p) {
106
107 // table is never built for GenericIon
108 if(p->GetParticleName() == "GenericIon") {
109 if(xSectionTable) {
111 delete xSectionTable;
112 xSectionTable = nullptr;
113 }
114
115 // table is always built for low mass particles
116 } else if(p->GetPDGMass() < 4.5*GeV || ForceBuildTableFlag()) {
117
119 idxTable = 0;
120 G4LossTableBuilder* builder =
122 if(IsMaster()) {
125 emin = std::max(emin, param->MinKinEnergy());
126 emax = std::min(emax, param->MaxKinEnergy());
127 if(emin < emax) {
129 emin, emax, true);
130 }
131 }
132 }
133 }
134 return change;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
140{
141 if(IsLocked()) { return; }
143 if(std::abs(part->GetPDGEncoding()) == 11) {
145 facrange = param->MscRangeFactor();
147 } else {
149 facrange = param->MscMuHadRangeFactor();
151 }
152 skin = param->MscSkin();
153 facgeom = param->MscGeomFactor();
154 facsafety = param->MscSafetyFactor();
155 lambdalimit = param->MscLambdaLimit();
156}
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159
160void G4VMscModel::DumpParameters(std::ostream& out) const
161{
162 G4String alg = "UseSafety";
163 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
164 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
165 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
166
167 out << std::setw(22) << "StepLim=" << alg << " Rfact=" << facrange
168 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
169 << " Skin=" << skin << " Llimit=" << lambdalimit << G4endl;
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173
174void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
176 const G4DynamicParticle*,
178{}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fMinimal
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
static G4EmParameters * Instance()
G4double MscMuHadRangeFactor() const
G4double MinKinEnergy() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscSafetyFactor() const
G4MscStepLimitType MscStepLimitType() const
G4double MscGeomFactor() const
G4double MaxKinEnergy() const
G4bool LateralDisplacement() const
G4bool MuHadLateralDisplacement() const
G4double MscLambdaLimit() const
G4double MscRangeFactor() const
G4double MscSkin() const
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
const G4String & GetParticleName() const
void clearAndDestroy()
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:440
size_t idxTable
Definition: G4VEmModel.hh:444
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:439
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:659
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:743
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:666
G4bool IsLocked() const
Definition: G4VEmModel.hh:867
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:701
void DumpParameters(std::ostream &out) const
Definition: G4VMscModel.cc:160
G4double facrange
Definition: G4VMscModel.hh:193
G4double skin
Definition: G4VMscModel.hh:196
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:59
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:91
G4double lambdalimit
Definition: G4VMscModel.hh:198
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:203
~G4VMscModel() override
Definition: G4VMscModel.cc:85
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
Definition: G4VMscModel.cc:174
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4double facsafety
Definition: G4VMscModel.hh:195
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:139
G4double facgeom
Definition: G4VMscModel.hh:194
Definition: DoubConv.h:17
#define DBL_MAX
Definition: templates.hh:62