Geant4_10
G4GoudsmitSaundersonMscModel.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: G4GoudsmitSaundersonMscModel.hh 67990 2013-03-13 10:56:28Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 //
31 // GEANT4 Class header file
32 //
33 //
34 // File name: G4GoudsmitSaundersonMscModel
35 //
36 // Author: Omrane Kadri
37 //
38 // Creation date: 20.02.2009
39 //
40 // Modifications:
41 // 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
42 // 12.05.2010 O.Kadri: adding Qn1 and Qn12 as private doubles
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:G4WentzelVIModel G4_v9.3.
57 //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
58 
59 #ifndef G4GoudsmitSaundersonMscModel_h
60 #define G4GoudsmitSaundersonMscModel_h 1
61 
63 
64 #include "G4VMscModel.hh"
65 #include "G4PhysicsTable.hh"
66 #include "globals.hh"
67 
68 class G4DataVector;
70 class G4LossTableManager;
72 
74 {
75 public:
76 
77  G4GoudsmitSaundersonMscModel(const G4String& nam = "GoudsmitSaunderson");
78 
80 
81  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
82 
83  void StartTracking(G4Track*);
84 
86  G4double KineticEnergy,
87  G4double AtomicNumber, G4double,
89 
91  G4double safety);
92 
93  virtual G4double ComputeTruePathLengthLimit(const G4Track& track,
94  G4double& currentMinimalStep);
95 
96  virtual G4double ComputeGeomPathLength(G4double truePathLength);
97 
98  virtual G4double ComputeTrueStepLength(G4double geomStepLength);
99 
100 private:
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
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 lambda0,lambda1,lambda11;
128  G4double mass;
129  G4int currentMaterialIndex;
130 
131  G4bool firstStep;
132  G4bool inside;
133  G4bool insideskin;
134 
135  G4GoudsmitSaundersonTable* GSTable;
136  G4LossTableManager* theManager;
137  const G4ParticleDefinition* particle;
138  G4ParticleChangeForMSC* fParticleChange;
139  const G4MaterialCutsCouple* currentCouple;
140 
141  static G4double ener[106];
142  static G4double TCSE[103][106];
143  static G4double FTCSE[103][106];
144  static G4double TCSP[103][106];
145  static G4double FTCSP[103][106];
146 
147 };
148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149 
150 inline
151 void G4GoudsmitSaundersonMscModel::SetParticle(const G4ParticleDefinition* p)
152 {
153  if (p != particle) {
154  particle = p;
155  mass = p->GetPDGMass();
156  charge = p->GetPDGCharge()/CLHEP::eplus;
157  }
158 }
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161 
162 #endif
163 
G4GoudsmitSaundersonMscModel(const G4String &nam="GoudsmitSaunderson")
virtual G4double ComputeTrueStepLength(G4double geomStepLength)
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &currentMinimalStep)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
bool G4bool
Definition: G4Types.hh:79
G4double GetPDGMass() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *particle, G4double KineticEnergy, G4double AtomicNumber, G4double, G4double, G4double)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
virtual G4double ComputeGeomPathLength(G4double truePathLength)