Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4WentzelOKandVIxSection.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // -------------------------------------------------------------------
 28 //
 29 //
 30 // GEANT4 Class header file
 31 //
 32 //
 33 // File name:     G4WentzelOKandVIxSection
 34 //
 35 // Authors:       V.Ivanchenko and O.Kadri 
 36 //
 37 // Creation date: 21.05.2010 
 38 //
 39 // Modifications:
 40 //
 41 //
 42 // Class Description:
 43 //
 44 // Implementation of the computation of total and transport cross sections,
 45 // sample scattering angle for the single scattering case.
 46 // to be used by single and multiple scattering models. References:
 47 // 1) G.Wentzel, Z. Phys. 40 (1927) 590.
 48 // 2) J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
 49 //
 50 // -------------------------------------------------------------------
 51 //
 52 
 53 #ifndef G4WentzelOKandVIxSection_h
 54 #define G4WentzelOKandVIxSection_h 1
 55 
 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 57 
 58 #include "globals.hh"
 59 #include "G4Material.hh"
 60 #include "G4Element.hh"
 61 #include "G4ElementVector.hh"
 62 #include "G4NistManager.hh"
 63 #include "G4NuclearFormfactorType.hh"
 64 #include "G4ThreeVector.hh"
 65 #include "G4Pow.hh"
 66 
 67 class G4ParticleDefinition;
 68 class G4ScreeningMottCrossSection;
 69 
 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 71 
 72 class G4WentzelOKandVIxSection 
 73 {
 74 
 75 public:
 76 
 77   explicit G4WentzelOKandVIxSection(G4bool combined=true);
 78 
 79   virtual ~G4WentzelOKandVIxSection();
 80 
 81   void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
 82 
 83   void SetupParticle(const G4ParticleDefinition*);
 84 
 85   // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering
 86   // cut = DBL_MAX means no scattering off electrons 
 87   virtual G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
 88   G4double SetupTarget(G4int Z, G4double cut);
 89 
 90   G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax);
 91  
 92   G4ThreeVector& SampleSingleScattering(G4double CosThetaMin,
 93           G4double CosThetaMax,
 94           G4double elecRatio);
 95 
 96   G4double ComputeSecondTransportMoment(G4double CosThetaMax);
 97 
 98   inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
 99                G4double CosThetaMax);
100  
101   inline G4double ComputeElectronCrossSection(G4double CosThetaMin,
102                 G4double CosThetaMax);
103    
104   inline void SetTargetMass(G4double value);
105 
106   inline G4double GetMomentumSquare() const;
107 
108   inline G4double GetCosThetaNuc() const;
109 
110   inline G4double GetCosThetaElec() const;
111 
112   //  hide assignment operator
113   G4WentzelOKandVIxSection & operator=
114   (const G4WentzelOKandVIxSection &right) = delete;
115   G4WentzelOKandVIxSection(const  G4WentzelOKandVIxSection&) = delete;
116 
117 protected:
118 
119   void ComputeMaxElectronScattering(G4double cut);
120 
121   void InitialiseA();
122 
123   inline G4double FlatFormfactor(G4double x);
124 
125   const G4ParticleDefinition* theProton;
126   const G4ParticleDefinition* theElectron;
127   const G4ParticleDefinition* thePositron;
128   const G4ParticleDefinition* particle = nullptr;
129   const G4Material* currentMaterial = nullptr;
130 
131   G4NistManager* fNistManager;
132   G4Pow*         fG4pow;
133 
134   G4ScreeningMottCrossSection* fMottXSection = nullptr;
135 
136   G4ThreeVector temp;
137 
138   // single scattering parameters
139   G4double coeff;
140   G4double cosTetMaxElec = 1.0;
141   G4double cosTetMaxNuc = 1.0;
142 
143   // for the combined mode it is cos(thetaMax)
144   // for single scattering it is cos(thetaMin)
145   G4double cosThetaMax = 1.0;
146 
147   G4double chargeSquare = 0.0;
148   G4double charge3 = 0.0;
149   G4double spin = 0.0;
150   G4double mass = 0.0;
151   G4double tkin = 0.0;
152   G4double mom2 = 0.0;
153   G4double momCM2 = 0.0;
154   G4double invbeta2 = 1.0;
155   G4double kinFactor = 1.0;
156   G4double etag = DBL_MAX;
157   G4double ecut = DBL_MAX;
158 
159   // target
160   G4double targetMass;
161   G4double screenZ = 0.0;
162   G4double formfactA = 0.0;
163   G4double factorA2 = 0.0;
164   G4double factB = 0.0;
165   G4double factD = 0.0;
166   G4double fMottFactor = 1.0;
167   G4double gam0pcmp = 1.0;
168   G4double pcmp2 = 1.0;
169 
170   // integer parameters
171   G4int targetZ = 0;
172   G4int nwarnings = 0;
173 
174   G4NuclearFormfactorType fNucFormfactor = fExponentialNF;
175 
176   G4bool isCombined;
177 
178   static G4double ScreenRSquareElec[100];
179   static G4double ScreenRSquare[100];
180   static G4double FormFactor[100];
181 };
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
185 inline void G4WentzelOKandVIxSection::SetTargetMass(G4double value)
186 {
187   targetMass = value;
188   factD = std::sqrt(mom2)/value;
189 }
190 
191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
192 
193 inline G4double G4WentzelOKandVIxSection::GetMomentumSquare() const
194 {
195   return mom2;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199 
200 inline G4double G4WentzelOKandVIxSection::GetCosThetaNuc() const
201 {
202   return cosTetMaxNuc;
203 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206 
207 inline G4double G4WentzelOKandVIxSection::GetCosThetaElec() const
208 {
209   return cosTetMaxElec;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213 
214 inline G4double 
215 G4WentzelOKandVIxSection::ComputeNuclearCrossSection(G4double cosTMin,
216                  G4double cosTMax)
217 {
218   return targetZ*kinFactor*fMottFactor*(cosTMin - cosTMax)/
219     ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
220 }
221 
222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
223 
224 inline G4double 
225 G4WentzelOKandVIxSection::ComputeElectronCrossSection(G4double cosTMin,
226                   G4double cosTMax)
227 {
228   G4double cost1 = std::max(cosTMin,cosTetMaxElec);
229   G4double cost2 = std::max(cosTMax,cosTetMaxElec);
230   return (cost1 <= cost2) ? 0.0 : kinFactor*fMottFactor*(cost1 - cost2)/
231     ((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
232 }
233 
234 inline G4double G4WentzelOKandVIxSection::FlatFormfactor(G4double x)
235 {
236   return 3.0*(std::sin(x) - x*std::cos(x))/(x*x*x);
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
240 
241 #endif
242 
243