source: trunk/source/processes/electromagnetic/standard/include/G4UrbanMscModel92.hh @ 1315

Last change on this file since 1315 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

File size: 6.5 KB
Line 
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: G4UrbanMscModel92.hh,v 1.1 2009/11/01 13:04:12 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29// -------------------------------------------------------------------
30//
31//
32// GEANT4 Class header file
33//
34//
35// File name:     G4UrbanMscModel92
36//
37// Author:        Laszlo Urban
38//
39// Creation date: 06.03.2008
40//
41// Modifications:
42//
43// 28-10-2009  V.Ivanchenko moved G4UrbanMscModel to G4UrbanMscModel92,
44//             now it is a frozen version of the Urban model corresponding
45//             to g4 9.2
46//
47// Class Description:
48//
49// Implementation of the model of multiple scattering based on
50// H.W.Lewis Phys Rev 78 (1950) 526 and L.Urban model
51
52// -------------------------------------------------------------------
53//
54
55#ifndef G4UrbanMscModel92_h
56#define G4UrbanMscModel92_h 1
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
60#include "G4VMscModel.hh"
61#include "G4PhysicsTable.hh"
62#include "G4MscStepLimitType.hh"
63
64class G4ParticleChangeForMSC;
65class G4SafetyHelper;
66class G4LossTableManager;
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70class G4UrbanMscModel92 : public G4VMscModel
71{
72
73public:
74
75  G4UrbanMscModel92(const G4String& nam = "UrbanMsc92");
76
77  virtual ~G4UrbanMscModel92();
78
79  void Initialise(const G4ParticleDefinition*, const G4DataVector&);
80
81  G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle,
82                                      G4double KineticEnergy,
83                                      G4double AtomicNumber,
84                                      G4double AtomicWeight=0., 
85                                      G4double cut =0.,
86                                      G4double emax=DBL_MAX);
87
88  void SampleScattering(const G4DynamicParticle*,
89                        G4double safety);
90
91  G4double ComputeTruePathLengthLimit(const G4Track& track,
92                                      G4PhysicsTable* theLambdaTable,
93                                      G4double currentMinimalStep);
94
95  G4double ComputeGeomPathLength(G4double truePathLength);
96
97  G4double ComputeTrueStepLength(G4double geomStepLength);
98
99  G4double ComputeTheta0(G4double truePathLength,
100                         G4double KineticEnergy);
101
102private:
103
104  G4double SimpleScattering(G4double xmeanth, G4double x2meanth);
105
106  G4double SampleCosineTheta(G4double trueStepLength, G4double KineticEnergy);
107
108  G4double SampleDisplacement();
109
110  G4double LatCorrelation();
111
112  inline G4double GetLambda(G4double kinEnergy);
113
114  inline void SetParticle(const G4ParticleDefinition*);
115
116  inline void UpdateCache();
117
118  //  hide assignment operator
119  G4UrbanMscModel92 & operator=(const  G4UrbanMscModel92 &right);
120  G4UrbanMscModel92(const  G4UrbanMscModel92&);
121
122  const G4ParticleDefinition* particle;
123  G4ParticleChangeForMSC*     fParticleChange;
124
125  G4PhysicsTable*             theLambdaTable;
126  const G4MaterialCutsCouple* couple;
127  G4LossTableManager*         theManager;
128
129  G4double mass;
130  G4double charge,ChargeSquare;
131  G4double masslimite,lambdalimit,fr;
132
133  G4double taubig;
134  G4double tausmall;
135  G4double taulim;
136  G4double currentTau;
137  G4double tlimit;
138  G4double tlimitmin;
139  G4double tlimitminfix;
140  G4double tgeom;
141
142  G4double geombig;
143  G4double geommin;
144  G4double geomlimit;
145  G4double skindepth;
146  G4double smallstep;
147
148  G4double presafety;
149
150  G4double lambda0;
151  G4double lambdaeff;
152  G4double tPathLength;
153  G4double zPathLength;
154  G4double par1,par2,par3;
155
156  G4double stepmin;
157
158  G4double currentKinEnergy;
159  G4double currentRange; 
160  G4double rangeinit;
161  G4double currentRadLength;
162
163  G4double theta0max,rellossmax;
164  G4double third;
165
166  G4int    currentMaterialIndex;
167
168  G4double y;
169  G4double Zold;
170  G4double Zeff,Z2,Z23,lnZ;
171  G4double coeffth1,coeffth2;
172  G4double coeffc1,coeffc2;
173  G4double scr1ini,scr2ini,scr1,scr2;
174
175  G4bool   isInitialized;
176  G4bool   inside;
177  G4bool   insideskin;
178};
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
182
183inline
184G4double G4UrbanMscModel92::GetLambda(G4double e)
185{
186  G4double x;
187  if(theLambdaTable) {
188    G4bool b;
189    x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
190  } else {
191    x = CrossSection(couple,particle,e);
192  }
193  if(x > DBL_MIN) x = 1./x;
194  else            x = DBL_MAX;
195  return x;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199
200inline
201void G4UrbanMscModel92::SetParticle(const G4ParticleDefinition* p)
202{
203  if (p != particle) {
204    particle = p;
205    mass = p->GetPDGMass();
206    charge = p->GetPDGCharge()/eplus;
207    ChargeSquare = charge*charge;
208  }
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
212
213inline
214void G4UrbanMscModel92::UpdateCache()                                   
215{
216    lnZ = std::log(Zeff);
217    coeffth1 = 0.885+lnZ*(0.104-0.0170*lnZ);
218    coeffth2 = 0.028+lnZ*(0.012-0.00125*lnZ);
219    coeffc1  = 2.134-lnZ*(0.1045-0.00602*lnZ);
220    coeffc2  = 0.001126-lnZ*(0.0001089+0.0000247*lnZ);
221    Z2 = Zeff*Zeff;
222    Z23 = std::exp(2.*lnZ/3.);
223    scr1     = scr1ini*Z23;
224    scr2     = scr2ini*Z2*ChargeSquare;
225  //  lastMaterial = couple->GetMaterial();
226    Zold = Zeff;
227}
228
229#endif
230
Note: See TracBrowser for help on using the repository browser.