G4WentzelVIRelXSection.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 //
00031 // GEANT4 Class header file
00032 //
00033 //
00034 // File name:     G4WentzelVIRelXSection
00035 //
00036 // Authors:       V.Ivanchenko  
00037 //
00038 // Creation date: 08.06.2012 from G4WentzelOKandVIxSection 
00039 //
00040 // Modifications:
00041 //
00042 //
00043 // Class Description:
00044 //
00045 // Implementation of the computation of total and transport cross sections,
00046 // sample scattering angle for the single scattering case.
00047 // to be used by single and multiple scattering models. References:
00048 // 1) G.Wentzel, Z. Phys. 40 (1927) 590.
00049 // 2) J.M. Fernandez-Varea et al., NIM B73 (1993) 447.
00050 //
00051 // -------------------------------------------------------------------
00052 //
00053 
00054 #ifndef G4WentzelVIRelXSection_h
00055 #define G4WentzelVIRelXSection_h 1
00056 
00057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00058 
00059 #include "globals.hh"
00060 #include "G4Material.hh"
00061 #include "G4Element.hh"
00062 #include "G4ElementVector.hh"
00063 #include "G4NistManager.hh"
00064 #include "G4ThreeVector.hh"
00065 #include "G4Pow.hh"
00066 
00067 class G4ParticleDefinition;
00068 
00069 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00070 
00071 class G4WentzelVIRelXSection 
00072 {
00073 
00074 public:
00075 
00076   G4WentzelVIRelXSection();
00077 
00078   virtual ~G4WentzelVIRelXSection();
00079 
00080   void Initialise(const G4ParticleDefinition*, G4double CosThetaLim);
00081 
00082   void SetupParticle(const G4ParticleDefinition*);
00083 
00084   // return cos(ThetaMax) for msc and cos(thetaMin) for single scattering
00085   // cut = DBL_MAX means no scattering off electrons 
00086   G4double SetupTarget(G4int Z, G4double cut = DBL_MAX);
00087 
00088   G4double ComputeTransportCrossSectionPerAtom(G4double CosThetaMax);
00089  
00090   G4ThreeVector SampleSingleScattering(G4double CosThetaMin,
00091                                        G4double CosThetaMax,
00092                                        G4double elecRatio = 0.0);
00093 
00094   inline G4double ComputeNuclearCrossSection(G4double CosThetaMin,
00095                                              G4double CosThetaMax);
00096  
00097   inline G4double ComputeElectronCrossSection(G4double CosThetaMin,
00098                                               G4double CosThetaMax);
00099  
00100   inline G4double SetupKinematic(G4double kinEnergy, const G4Material* mat);
00101   
00102   inline void SetTargetMass(G4double value);
00103 
00104   //obsolete method
00105   inline void SetRelativisticMass(G4double value);
00106 
00107   inline G4double GetMomentumSquare() const;
00108 
00109   inline G4double GetCosThetaNuc() const;
00110 
00111   inline G4double GetCosThetaElec() const;
00112 
00113 private:
00114 
00115   void ComputeMaxElectronScattering(G4double cut);
00116 
00117   //  hide assignment operator
00118   G4WentzelVIRelXSection & operator=(const  G4WentzelVIRelXSection &right);
00119   G4WentzelVIRelXSection(const  G4WentzelVIRelXSection&);
00120 
00121   const G4ParticleDefinition* theProton;
00122   const G4ParticleDefinition* theElectron;
00123   const G4ParticleDefinition* thePositron;
00124   const G4Material* currentMaterial;
00125 
00126   G4NistManager*  fNistManager;
00127   G4Pow*          fG4pow;
00128 
00129   G4double numlimit;
00130 
00131   G4double elecXSRatio;
00132 
00133   // integer parameters
00134   G4int    nwarnings;
00135   G4int    nwarnlimit;
00136 
00137   // single scattering parameters
00138   G4double coeff;
00139   G4double cosTetMaxElec;
00140   G4double cosTetMaxNuc;
00141   G4double cosThetaMax;
00142   G4double alpha2;
00143 
00144   // projectile
00145   const G4ParticleDefinition* particle;
00146 
00147   G4double chargeSquare;
00148   G4double charge3;
00149   G4double spin;
00150   G4double mass;
00151   G4double tkin;
00152   G4double mom2;
00153   G4double momCM2;
00154   G4double invbeta2;
00155   G4double kinFactor;
00156   G4double etag;
00157   G4double ecut;
00158   G4double lowEnergyLimit;
00159 
00160   // target
00161   G4int    targetZ;
00162   G4double targetMass;
00163   G4double screenZ;
00164   G4double formfactA;
00165   G4double factorA2;
00166   G4double factB;
00167   G4double factB1;
00168   G4double factD;
00169   G4double gam0pcmp;
00170   G4double pcmp2;
00171 
00172   static G4double ScreenRSquare[100];
00173   static G4double FormFactor[100];
00174 };
00175 
00176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00177 
00178 inline G4double 
00179 G4WentzelVIRelXSection::SetupKinematic(G4double ekin, const G4Material* mat)
00180 {
00181   if(ekin != tkin || mat != currentMaterial) { 
00182     currentMaterial = mat;
00183     tkin  = ekin;
00184     mom2  = tkin*(tkin + 2.0*mass);
00185     invbeta2 = 1.0 +  mass*mass/mom2;
00186     factB = spin/invbeta2;
00187     cosTetMaxNuc = 
00188         std::max(cosThetaMax,1.-factorA2*mat->GetIonisation()->GetInvA23()/mom2);
00189   } 
00190   return cosTetMaxNuc;
00191 }
00192 
00193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00194 
00195 inline void G4WentzelVIRelXSection::SetTargetMass(G4double value)
00196 {
00197   targetMass = value;
00198   factD = std::sqrt(mom2)/value;
00199 }
00200 
00201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00202 
00203 inline void G4WentzelVIRelXSection::SetRelativisticMass(G4double value)
00204 {
00205   mass = value;
00206 }
00207 
00208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00209 
00210 inline G4double G4WentzelVIRelXSection::GetMomentumSquare() const
00211 {
00212   return mom2;
00213 }
00214 
00215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00216 
00217 inline G4double G4WentzelVIRelXSection::GetCosThetaNuc() const
00218 {
00219   return cosTetMaxNuc;
00220 }
00221 
00222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00223 
00224 inline G4double G4WentzelVIRelXSection::GetCosThetaElec() const
00225 {
00226   return cosTetMaxElec;
00227 }
00228 
00229 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00230 
00231 inline G4double 
00232 G4WentzelVIRelXSection::ComputeNuclearCrossSection(G4double cosTMin,
00233                                                      G4double cosTMax)
00234 {
00235   G4double xsec = 0.0;
00236   if(cosTMax < cosTMin) {
00237     xsec = targetZ*kinFactor*(cosTMin - cosTMax)/
00238       ((1.0 - cosTMin + screenZ)*(1.0 - cosTMax + screenZ));
00239   }
00240   return xsec;
00241 }
00242 
00243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00244 
00245 inline G4double 
00246 G4WentzelVIRelXSection::ComputeElectronCrossSection(G4double cosTMin,
00247                                                       G4double cosTMax)
00248 {
00249   G4double xsec = 0.0;
00250   G4double cost1 = std::max(cosTMin,cosTetMaxElec);
00251   G4double cost2 = std::max(cosTMax,cosTetMaxElec);
00252   if(cost1 > cost2) {
00253     xsec = kinFactor*(cost1 - cost2)/((1.0 - cost1 + screenZ)*(1.0 - cost2 + screenZ));
00254   }
00255   return xsec;
00256 }
00257 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00258 
00259 #endif
00260 

Generated on Mon May 27 17:50:25 2013 for Geant4 by  doxygen 1.4.7