source: trunk/source/processes/electromagnetic/standard/include/G4GoudsmitSaundersonMscModel.hh@ 1190

Last change on this file since 1190 was 1058, checked in by garnier, 17 years ago

file release beta

  • Property svn:executable set to *
File size: 6.5 KB
RevLine 
[1058]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: G4GoudsmitSaundersonMscModel.hh,v 1.1 2009/03/05 18:48:30 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31//
32// GEANT4 Class header file
33//
34//
35// File name: G4GoudsmitSaundersonMscModel
36//
37// Author: Omrane Kadri
38//
39// Creation date: 20.02.2009
40//
41// Modifications:
42// 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
43//
44// Class description:
45//
46// Multiple scattering model using classical Goudsmit-Saunderson model
47//
48//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
49//REFERENCES:
50//Ref.1:E. Benedito et al.,"Mixed simulation ... cross-sections", NIMB 174 (2001) pp 91-110;
51//Ref.2:I. Kawrakow et al.,"On the condensed ... transport",NIMB 142 (1998) pp 253-280;
52//Ref.3:I. Kawrakow et al.,"On the representation ... calculations",NIMB 134 (1998) pp 325-336;
53//Ref.4:I. Kawrakow et al.,"The EGSnrc code ... Transport",NRCC Report PIRS-701, Sept. 21, 2006;
54//Ref.5:F. Salvat et al.,"ELSEPA--Dirac partial ...molecules", Comp. Phys. Comm. 165 (2005) pp 157-190;
55//Ref.6:G4UrbanMscModel G4_v9.1Ref09;
56//Ref.7:G4eCoulombScatteringModel G4_v9.1Ref09.
57//@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
58
59#ifndef G4GoudsmitSaundersonMscModel_h
60#define G4GoudsmitSaundersonMscModel_h 1
61
62#include "G4VMscModel.hh"
63#include "G4PhysicsTable.hh"
64#include "globals.hh"
65#include "G4DataInterpolation.hh"
66
67class G4DataVector;
68class G4ParticleChangeForMSC;
69class G4LossTableManager;
70class G4GoudsmitSaundersonTable;
71class G4PhysicsTable;
72
73class G4GoudsmitSaundersonMscModel : public G4VMscModel
74{
75public:
76
77 G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson");
78
79 virtual ~G4GoudsmitSaundersonMscModel();
80
81 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82
83 virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition* particle,
84 G4double KineticEnergy,
85 G4double AtomicNumber, G4double,
86 G4double, G4double);
87
88 virtual void SampleScattering(const G4DynamicParticle* ,
89 G4double );
90
91 virtual G4double ComputeTruePathLengthLimit(
92 const G4Track& track,
93 G4PhysicsTable* theLambdaTable,
94 G4double currentMinimalStep);
95
96 virtual G4double ComputeGeomPathLength(G4double truePathLength);
97
98 virtual G4double ComputeTrueStepLength(G4double geomStepLength);
99
100private:
101 void SampleCosineTheta(G4double,G4double,G4double &,G4double &);
102
103 void CalculateIntegrals(const G4ParticleDefinition* ,G4double,G4double,
104 G4double&,G4double&);
105
106 void LoadELSEPAXSections();
107
108 inline void SetParticle(const G4ParticleDefinition* p);
109
110 inline G4double GetLambda(G4double);
111
112 // hide assignment operator
113 G4GoudsmitSaundersonMscModel & operator=(const G4GoudsmitSaundersonMscModel &right);
114 G4GoudsmitSaundersonMscModel(const G4GoudsmitSaundersonMscModel&);
115
116 G4double lowKEnergy;
117 G4double highKEnergy;
118 G4double currentKinEnergy;
119 G4double currentRange;
120
121 G4double smallstep,tlimitminfix,skindepth;
122 G4double fr,rangeinit,masslimite,tgeom;
123 G4double par1,par2,par3,zPathLength,truePathLength;
124 G4double tausmall,taulim,tlimit,tlimitmin,geommin,geombig;
125 G4double charge,lambdalimit;
126 G4double tPathLength,stepmin ;
127 G4double lambda1,lambda11;
128 G4double mass;
129 G4double lambda0;
130 G4int currentMaterialIndex;
131
132 G4bool inside;
133 G4bool insideskin;
134 G4bool isInitialized;
135
136 G4PhysicsTable* theLambdaTable;
137 G4GoudsmitSaundersonTable* GSTable;
138 G4LossTableManager* theManager;
139 const G4ParticleDefinition* particle;
140 G4ParticleChangeForMSC* fParticleChange;
141 G4DataInterpolation* MyValue;
142 const G4MaterialCutsCouple* currentCouple;
143
144 static G4double ener[106];
145 static G4double TCSE[103][106];
146 static G4double FTCSE[103][106];
147 static G4double TCSP[103][106];
148 static G4double FTCSP[103][106];
149
150};
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
153inline
154void G4GoudsmitSaundersonMscModel::SetParticle(const G4ParticleDefinition* p)
155{
156 if (p != particle) {
157 particle = p;
158 mass = p->GetPDGMass();
159 charge = p->GetPDGCharge()/eplus;
160 }
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164
165inline G4double G4GoudsmitSaundersonMscModel::GetLambda(G4double e)
166{
167 G4double x;
168 if(theLambdaTable) {
169 G4bool b;
170 x = ((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b);
171 } else {
172 x = CrossSection(currentCouple,particle,e);
173 }
174
175 if(x > DBL_MIN) x = 1./x;
176 else x = DBL_MAX;
177 return x;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
182#endif
183
Note: See TracBrowser for help on using the repository browser.