source: trunk/source/processes/electromagnetic/standard/include/G4UrbanMscModel2.hh @ 1006

Last change on this file since 1006 was 1006, checked in by garnier, 15 years ago

fichiers oublies

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