Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UrbanMscModel90.hh
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//
31// GEANT4 Class header file
32//
33//
34// File name: G4UrbanMscModel90
35//
36// Author: V.Ivanchenko clone Laszlo Urban model
37//
38// Creation date: 07.12.2007
39//
40// Modifications:
41//
42//
43// Class Description:
44//
45// Implementation of the model of multiple scattering based on
46// H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
47
48// -------------------------------------------------------------------
49//
50
51#ifndef G4UrbanMscModel90_h
52#define G4UrbanMscModel90_h 1
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
57
58#include "G4VMscModel.hh"
59#include "G4MscStepLimitType.hh"
60
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
67{
68
69public:
70
71 G4UrbanMscModel90(const G4String& nam = "UrbanMsc90");
72
73 virtual ~G4UrbanMscModel90();
74
75 void Initialise(const G4ParticleDefinition*, const G4DataVector&);
76
78
80 G4double KineticEnergy,
81 G4double AtomicNumber,
82 G4double AtomicWeight=0.,
83 G4double cut =0.,
84 G4double emax=DBL_MAX);
85
86 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
88 const G4DynamicParticle*,
90 G4double);
91
93
95 G4double& currentMinimalStep);
96
98
100
101 G4double ComputeTheta0(G4double truePathLength,
102 G4double KineticEnergy);
103
104private:
105
106 G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
107
108 G4double SampleDisplacement();
109
110 G4double LatCorrelation();
111
112 inline void SetParticle(const G4ParticleDefinition*);
113
114 // hide assignment operator
115 G4UrbanMscModel90 & operator=(const G4UrbanMscModel90 &right);
117
118 const G4ParticleDefinition* particle;
119 G4ParticleChangeForMSC* fParticleChange;
120
121 const G4MaterialCutsCouple* couple;
122 G4LossTableManager* theManager;
123
124 G4double mass;
125 G4double charge;
126
127 G4double masslimite,masslimitmu;
128
129 G4double taubig;
130 G4double tausmall;
131 G4double taulim;
132 G4double currentTau;
133 G4double frscaling1,frscaling2;
134 G4double tlimit;
135 G4double tlimitmin;
136 G4double tlimitminfix;
137
138 G4double nstepmax;
139 G4double geombig;
140 G4double geommin;
141 G4double geomlimit;
142 G4double skindepth;
143 G4double smallstep;
144
145 G4double presafety;
146
147 G4double lambda0;
148 G4double lambdaeff;
149 G4double tPathLength;
150 G4double zPathLength;
151 G4double par1,par2,par3;
152
153 G4double stepmin;
154
155 G4double currentKinEnergy;
156 G4double currentRange;
157 G4double currentRadLength;
158
159 G4double Zeff;
160
161 G4int currentMaterialIndex;
162
163 G4bool firstStep;
164 G4bool inside;
165 G4bool insideskin;
166
167};
168
169//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
170
171inline
172void G4UrbanMscModel90::SetParticle(const G4ParticleDefinition* p)
173{
174 if (p != particle) {
175 particle = p;
176 mass = p->GetPDGMass();
177 charge = p->GetPDGCharge()/CLHEP::eplus;
178 }
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183#endif
184
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetPDGCharge() const
G4double ComputeTheta0(G4double truePathLength, G4double KineticEnergy)
G4double ComputeGeomPathLength(G4double truePathLength)
void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
G4double ComputeTrueStepLength(G4double geomStepLength)
G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double AtomicWeight=0., G4double cut=0., G4double emax=DBL_MAX)
void StartTracking(G4Track *)
#define DBL_MAX
Definition: templates.hh:83