Geant4 11.3.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 lambdalimit(1.*CLHEP::mm),
62 geomMin(1.e-6*CLHEP::mm),
63 geomMax(1.e50*CLHEP::mm),
64 fDisplacement(0.,0.,0.),
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 return change;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95
97{
98 if(IsLocked()) { return; }
100 if(std::abs(part->GetPDGEncoding()) == 11) {
102 facrange = param->MscRangeFactor();
104 } else {
106 facrange = param->MscMuHadRangeFactor();
108 }
109 skin = param->MscSkin();
110 facgeom = param->MscGeomFactor();
111 facsafety = param->MscSafetyFactor();
112 lambdalimit = param->MscLambdaLimit();
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116
117void G4VMscModel::DumpParameters(std::ostream& out) const
118{
119 G4String alg = "UseSafety";
120 if (steppingAlgorithm == fUseDistanceToBoundary) alg = "DistanceToBoundary";
121 else if (steppingAlgorithm == fMinimal) alg = "Minimal";
122 else if (steppingAlgorithm == fUseSafetyPlus) alg = "SafetyPlus";
123
124 out << std::setw(18) << "StepLim=" << alg << " Rfact=" << facrange
125 << " Gfact=" << facgeom << " Sfact=" << facsafety << " DispFlag:" << latDisplasment
126 << " Skin=" << skin << " Llim=" << lambdalimit/CLHEP::mm << " mm" << G4endl;
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130
131void G4VMscModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
133 const G4DynamicParticle*,
135{}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
138
140 const G4MaterialCutsCouple* couple)
141{
142 G4double x;
143 if (nullptr != ionisation) {
144 x = ionisation->GetDEDX(kinEnergy, couple);
145 } else {
146 const G4double q = part->GetPDGCharge()*inveplus;
147 x = dedx*q*q;
148 }
149 return x;
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153
155 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
156{
157 G4double x;
158 if (nullptr != ionisation) {
159 x = ionisation->GetDEDX(kinEnergy, couple, logKinEnergy);
160 } else {
161 const G4double q = part->GetPDGCharge()*inveplus;
162 x = dedx*q*q;
163 }
164 return x;
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
168
170 const G4MaterialCutsCouple* couple)
171{
172 // << ionisation << " " << part->GetParticleName() << G4endl;
173 localtkin = kinEnergy;
174 if (nullptr != ionisation) {
175 localrange = ionisation->GetRange(kinEnergy, couple);
176 } else {
177 const G4double q = part->GetPDGCharge()*inveplus;
178 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
179 }
180 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
181 return localrange;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185
187 const G4MaterialCutsCouple* couple, G4double logKinEnergy)
188{
189 //G4cout << "G4VMscModel::GetRange E(MeV)= " << kinEnergy << " "
190 // << ionisation << " " << part->GetParticleName() << G4endl;
191 localtkin = kinEnergy;
192 if (nullptr != ionisation) {
193 localrange = ionisation->GetRange(kinEnergy, couple, logKinEnergy);
194 } else {
195 const G4double q = part->GetPDGCharge()*inveplus;
196 localrange = kinEnergy/(dedx*q*q*couple->GetMaterial()->GetDensity());
197 }
198 //G4cout << "R(mm)= " << localrange << " " << ionisation << G4endl;
199 return localrange;
200}
201
202//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
203
205 G4double range, const G4MaterialCutsCouple* couple)
206{
207 G4double e;
208 //G4cout << "G4VMscModel::GetEnergy R(mm)= " << range << " " << ionisation
209 // << " Rlocal(mm)= " << localrange << " Elocal(MeV)= " << localtkin << G4endl;
210 if(nullptr != ionisation) { e = ionisation->GetKineticEnergy(range, couple); }
211 else {
212 e = localtkin;
213 if(localrange > range) {
214 G4double q = part->GetPDGCharge()*inveplus;
215 e -= (localrange - range)*dedx*q*q*couple->GetMaterial()->GetDensity();
216 }
217 }
218 return e;
219}
220
221//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
static G4EmParameters * Instance()
G4double MscMuHadRangeFactor() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4double MscSafetyFactor() const
G4MscStepLimitType MscStepLimitType() const
G4double MscGeomFactor() const
G4bool LateralDisplacement() const
G4bool MuHadLateralDisplacement() const
G4double MscLambdaLimit() const
G4double MscRangeFactor() const
G4double MscSkin() const
const G4Material * GetMaterial() const
G4double GetDensity() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
G4double inveplus
G4VParticleChange * pParticleChange
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
G4bool IsLocked() const
void DumpParameters(std::ostream &out) const
G4double GetDEDX(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double facrange
G4double skin
G4double geomMin
G4VMscModel(const G4String &nam)
G4ParticleChangeForMSC * GetParticleChangeForMSC(const G4ParticleDefinition *p=nullptr)
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
G4double geomMax
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4double lambdalimit
G4MscStepLimitType steppingAlgorithm
~G4VMscModel() override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
G4bool latDisplasment
G4double facsafety
G4ThreeVector fDisplacement
void InitialiseParameters(const G4ParticleDefinition *)
G4double facgeom