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

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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