Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/include/G4NeutronElectronElXsc.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 // Neutron-electron elastic cross section base on the integration of
 27 // the Rosenbluth differential xsc
 28 //
 29 // 16.05.17 V. Grichine
 30 //
 31 //
 32 
 33 #ifndef G4NeutronElectronElXsc_h
 34 #define G4NeutronElectronElXsc_h
 35 
 36 
 37 #include "globals.hh"
 38 #include "G4VCrossSectionDataSet.hh"
 39 #include "G4DynamicParticle.hh"
 40 
 41 // class G4ParticleDefinition;
 42 class G4PhysicsLogVector;
 43 class G4PhysicsTable;
 44 
 45 class G4NeutronElectronElXsc : public G4VCrossSectionDataSet
 46 {
 47 public:
 48    
 49   G4NeutronElectronElXsc();
 50   ~G4NeutronElectronElXsc();
 51 
 52   void Initialise();
 53 
 54   virtual
 55   G4bool IsElementApplicable(const G4DynamicParticle*, G4int Z, const G4Material*);
 56 
 57 
 58   virtual
 59   G4double GetElementCrossSection(const G4DynamicParticle*, 
 60           G4int Z, const G4Material*);
 61 
 62   G4double GetRosenbluthXsc(const G4DynamicParticle*, 
 63           G4int Z, const G4Material*);
 64 
 65   G4double XscIntegrand(G4double x);
 66 
 67    G4double GetElementNonRelXsc(const G4DynamicParticle*, 
 68           G4int Z, const G4Material*);
 69 
 70   G4double CalculateAm( G4double momentum);
 71 
 72   inline G4double GetAm(){return fAm;};
 73 
 74   void SetCutEnergy(G4double ec){fCutEnergy=ec;};
 75   G4double GetCutEnergy(){return fCutEnergy;};
 76 
 77   void SetBiasingFactor(G4double bf){fBiasingFactor=bf;};
 78 
 79 protected:
 80   G4double fM, fM2, fMv2, fme, fme2, fee, fee2;
 81   G4double fCofXsc;    // 
 82   G4double fAm;    //
 83   G4int fEnergyBin; 
 84   G4double fMinEnergy, fMaxEnergy, fCutEnergy; // minimal recoil electron energy detected
 85   G4double fBiasingFactor; // biasing xsc up
 86 
 87   G4PhysicsLogVector* fEnergyXscVector;
 88   static const G4double fXscArray[200];
 89 };
 90 
 91 
 92 
 93 ////////////////////////////////////////////////////////////////////
 94 //
 95 // return Wentzel atom screening correction for neutron-electron scattering
 96 
 97 inline  G4double G4NeutronElectronElXsc::CalculateAm( G4double momentum)
 98 {
 99   G4double k   = momentum/CLHEP::hbarc;
100   G4double ch  = 1.13;
101   G4double zn  = 1.77*k*CLHEP::Bohr_radius;
102   G4double zn2 = zn*zn;
103   fAm          = ch/zn2;
104 
105   return fAm;
106 }
107 
108 ////////////////////////////////////////////////////
109 //
110 // Slow electron (Tkin << me_c2) in the neutron rest frame
111 
112 inline G4double G4NeutronElectronElXsc::
113 GetElementNonRelXsc(const G4DynamicParticle* aPart, G4int ZZ,  
114            const G4Material*) 
115 {
116   G4double result(0.), te(0.), momentum(0.);
117 
118   te = aPart->GetKineticEnergy()*fme/fM;
119   momentum = std::sqrt( te*(te + 2.*fme) );
120   fAm = CalculateAm(momentum);
121 
122   result = 1. + std::log(1. +1./fAm);
123   result *= fCofXsc; //*energy;
124   result *= ZZ;  // incoherent sum over  all element electrons
125 
126   return result;
127 }
128 
129 #endif
130