Geant4 11.1.1
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 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{
67 dedx = 2.0*CLHEP::MeV*CLHEP::cm2/CLHEP::g;
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
78{
79 // recomputed for each new run
80 if(nullptr == safetyHelper) {
83 safetyHelper->InitialiseHelper();
84 }
85 G4ParticleChangeForMSC* change = nullptr;
86 if (nullptr != pParticleChange) {
87 change = static_cast<G4ParticleChangeForMSC*>(pParticleChange);
88 } else {
89 change = new G4ParticleChangeForMSC();
90 }
91 if(IsMaster() && nullptr != p) {
92
93 // table is always built for low mass particles
94 if(p->GetParticleName() != "GenericIon" &&
95 (p->GetPDGMass() < CLHEP::GeV || ForceBuildTableFlag()) ) {
96
98 G4LossTableBuilder* builder =
102 emin = std::max(emin, param->MinKinEnergy());
103 emax = std::min(emax, param->MaxKinEnergy());
104 if(emin < emax) {
106 emin, emax, useSpline);
107 }
108 }
109 }
110 return change;
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114
116{
117 if(IsLocked()) { return; }
119 if(std::abs(part->GetPDGEncoding()) == 11) {
121 facrange = param->MscRangeFactor();
123 } else {
125 facrange = param->MscMuHadRangeFactor();
127 }
128 skin = param->MscSkin();
129 facgeom = param->MscGeomFactor();
130 facsafety = param->MscSafetyFactor();
131 lambdalimit = param->MscLambdaLimit();
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
136void G4VMscModel::DumpParameters(std::ostream& out) const
137{
138 G4String alg = "UseSafety";
139 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
140 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
141 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
142
143 out << std::setw(18) << "StepLim=" << alg << " Rfact=" << facrange
144 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
145 << " Skin=" << skin << " Llim=" << lambdalimit/CLHEP::mm << " mm" << G4endl;
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
150void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
152 const G4DynamicParticle*,
154{}
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
157
159 const G4MaterialCutsCouple* couple)
160{
161 G4double x;
162 if (nullptr != ionisation) {
163 x = ionisation->GetDEDX(kinEnergy, couple);
164 } else {
165 const G4double q = part->GetPDGCharge()*inveplus;
166 x = dedx*q*q;
167 }
168 return x;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172
174 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
175{
176 G4double x;
177 if (nullptr != ionisation) {
178 x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
179 } else {
180 const G4double q = part->GetPDGCharge()*inveplus;
181 x = dedx*q*q;
182 }
183 return x;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
187
189 const G4MaterialCutsCouple* couple)
190{
191 // << ionisation << " " << part->GetParticleName() << G4endl;
192 localtkin = kinEnergy;
193 if (nullptr != ionisation) {
194 localrange = ionisation->GetRange(kinEnergy, couple);
195 } else {
196 const G4double q = part->GetPDGCharge()*inveplus;
197 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
198 }
199 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
200 return localrange;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204
206 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
207{
208 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
209 // << ionisation << " " << part->GetParticleName() << G4endl;
210 localtkin = kinEnergy;
211 if (nullptr != ionisation) {
212 localrange = ionisation->GetRange(kinEnergy, couple, logKinEnergy);
213 } else {
214 const G4double q = part->GetPDGCharge()*inveplus;
215 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
216 }
217 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
218 return localrange;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
222
224 G4double range, const G4MaterialCutsCouple* couple)
225{
226 G4double e;
227 //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
228 // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin << G4endl;
229 if(nullptr != ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
230 else {
231 e = localtkin;
232 if(localrange > range) {
233 G4double q = part->GetPDGCharge()*inveplus;
234 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
235 }
236 }
237 return e;
238}
239
240//....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 G4Material * GetMaterial() const
G4double GetDensity() const
Definition: G4Material.hh:175
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4PhysicsTable * xSectionTable
Definition: G4VEmModel.hh:422
G4double inveplus
Definition: G4VEmModel.hh:427
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:421
G4double HighEnergyActivationLimit() const
Definition: G4VEmModel.hh:648
G4double LowEnergyActivationLimit() const
Definition: G4VEmModel.hh:655
G4bool IsLocked() const
Definition: G4VEmModel.hh:856
G4bool ForceBuildTableFlag() const
Definition: G4VEmModel.hh:690
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:136
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:158
G4double facrange
Definition: G4VMscModel.hh:199
G4double skin
Definition: G4VMscModel.hh:202
G4VMscModel(const G4String &nam)
Definition: G4VMscModel.cc:59
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
Definition: G4VMscModel.cc:77
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:223
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:188
G4double lambdalimit
Definition: G4VMscModel.hh:204
G4MscStepLimitType steppingAlgorithm
Definition: G4VMscModel.hh:209
~G4VMscModel() override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
Definition: G4VMscModel.cc:150
G4bool latDisplasment
Definition: G4VMscModel.hh:212
G4double facsafety
Definition: G4VMscModel.hh:201
void InitialiseParameters(const G4ParticleDefinition *)
Definition: G4VMscModel.cc:115
G4double facgeom
Definition: G4VMscModel.hh:200
Definition: DoubConv.h:17