Geant4-11
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 lambdalimit(1.*CLHEP::mm),
62 geomMin(1.e-6*CLHEP::mm),
63 geomMax(1.e50*CLHEP::mm),
64 fDisplacement(0.,0.,0.),
65 steppingAlgorithm(fUseSafety)
66{
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
73{}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
76
79{
80 // recomputed for each new run
81 if(nullptr == safetyHelper) {
85 }
86 G4ParticleChangeForMSC* change = nullptr;
87 if (nullptr != pParticleChange) {
88 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
89 } else {
90 change = new G4ParticleChangeForMSC();
91 }
92 if(IsMaster() && nullptr != p) {
93
94 // table is always built for low mass particles
95 if(p->GetParticleName() != "GenericIon" &&
97
99 G4LossTableBuilder* builder =
103 emin = std::max(emin, param->MinKinEnergy());
104 emax = std::min(emax, param->MaxKinEnergy());
105 if(emin < emax) {
107 emin, emax, true);
108 }
109 }
110 }
111 return change;
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
117{
118 if(IsLocked()) { return; }
120 if(std::abs(part->GetPDGEncoding()) == 11) {
122 facrange = param->MscRangeFactor();
124 } else {
126 facrange = param->MscMuHadRangeFactor();
128 }
129 skin = param->MscSkin();
130 facgeom = param->MscGeomFactor();
131 facsafety = param->MscSafetyFactor();
132 lambdalimit = param->MscLambdaLimit();
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136
137void G4VMscModel::DumpParameters(std::ostream& out) const
138{
139 G4String alg = "UseSafety";
140 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
141 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
142 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
143
144 out << std::setw(18) << "StepLim=" << alg << " Rfact=" << facrange
145 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
146 << " Skin=" << skin << " Llim=" << lambdalimit/CLHEP::mm << " mm" << G4endl;
147}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
151void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
153 const G4DynamicParticle*,
155{}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
160 const G4MaterialCutsCouple* couple)
161{
162 G4double x;
163 if (nullptr != ionisation) {
164 x = ionisation->GetDEDX(kinEnergy, couple);
165 } else {
166 const G4double q = part->GetPDGCharge()*inveplus;
167 x = dedx*q*q;
168 }
169 return x;
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173
175 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
176{
177 G4double x;
178 if (nullptr != ionisation) {
179 x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
180 } else {
181 const G4double q = part->GetPDGCharge()*inveplus;
182 x = dedx*q*q;
183 }
184 return x;
185}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
188
190 const G4MaterialCutsCouple* couple)
191{
192 // << ionisation << " " << part->GetParticleName() << G4endl;
193 localtkin = kinEnergy;
194 if (nullptr != ionisation) {
195 localrange = ionisation->GetRange(kinEnergy, couple);
196 } else {
197 const G4double q = part->GetPDGCharge()*inveplus;
198 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
199 }
200 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
201 return localrange;
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205
207 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
208{
209 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
210 // << ionisation << " " << part->GetParticleName() << G4endl;
211 localtkin = kinEnergy;
212 if (nullptr != ionisation) {
213 localrange = ionisation->GetRange(kinEnergy, couple, logKinEnergy);
214 } else {
215 const G4double q = part->GetPDGCharge()*inveplus;
216 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
217 }
218 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
219 return localrange;
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223
225 G4double range, const G4MaterialCutsCouple* couple)
226{
227 G4double e;
228 //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
229 // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin << G4endl;
230 if(nullptr != ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
231 else {
232 e = localtkin;
233 if(localrange > range) {
234 G4double q = part->GetPDGCharge()*inveplus;
235 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
236 }
237 }
238 return e;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static const G4double emax
@ fMinimal
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
static constexpr double mm
Definition: G4SIunits.hh:95
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 G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:176
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:426
G4double inveplus
Definition: G4VEmModel.hh:431
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:662
G4bool IsMaster() const
Definition: G4VEmModel.hh:746
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:425
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:669
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:676
G4bool IsLocked() const
Definition: G4VEmModel.hh:877
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:711
G4double GetKineticEnergy(G4double range, const G4MaterialCutsCouple *)
G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple *)
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
void DumpParameters(std::ostream &out) const
Definition: G4VMscModel.cc:137
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:159
G4double facrange
Definition: G4VMscModel.hh:196
G4double skin
Definition: G4VMscModel.hh:199
G4double localtkin
Definition: G4VMscModel.hh:191
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:59
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:78
G4VEnergyLossProcess * ionisation
Definition: G4VMscModel.hh:187
G4double dedx
Definition: G4VMscModel.hh:190
G4double localrange
Definition: G4VMscModel.hh:192
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:224
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:189
G4double lambdalimit
Definition: G4VMscModel.hh:201
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:209
G4SafetyHelper * safetyHelper
Definition: G4VMscModel.hh:186
~G4VMscModel() override
Definition: G4VMscModel.cc:72
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
Definition: G4VMscModel.cc:151
G4bool latDisplasment
Definition: G4VMscModel.hh:206
G4double facsafety
Definition: G4VMscModel.hh:198
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:116
G4double facgeom
Definition: G4VMscModel.hh:197
Definition: DoubConv.h:17
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double GeV
static constexpr double MeV
static constexpr double g
static constexpr double cm2
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments