Geant4-11
G4NeutronElectronElXsc.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// 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"
39#include "G4DynamicParticle.hh"
40
41// class G4ParticleDefinition;
43class G4PhysicsTable;
44
46{
47public:
48
51
52 void Initialise();
53
54 virtual
56
57
58 virtual
60 G4int Z, const G4Material*);
61
63 G4int Z, const G4Material*);
64
66
68 G4int Z, const G4Material*);
69
70 G4double CalculateAm( G4double momentum);
71
72 inline G4double GetAm(){return fAm;};
73
76
78
79protected:
84 G4double fMinEnergy, fMaxEnergy, fCutEnergy; // minimal recoil electron energy detected
85 G4double fBiasingFactor; // biasing xsc up
86
88 static const G4double fXscArray[200];
89};
90
91
92
94//
95// return Wentzel atom screening correction for neutron-electron scattering
96
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
109//
110// Slow electron (Tkin << me_c2) in the neutron rest frame
111
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
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4double GetKineticEnergy() const
G4double CalculateAm(G4double momentum)
G4double GetRosenbluthXsc(const G4DynamicParticle *, G4int Z, const G4Material *)
static const G4double fXscArray[200]
G4double GetElementNonRelXsc(const G4DynamicParticle *, G4int Z, const G4Material *)
G4PhysicsLogVector * fEnergyXscVector
G4double XscIntegrand(G4double x)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
void SetBiasingFactor(G4double bf)
static constexpr double Bohr_radius
static constexpr double hbarc