Geant4-11
source
processes
hadronic
cross_sections
include
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
"
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
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
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
G4DynamicParticle.hh
G4double
double G4double
Definition:
G4Types.hh:83
G4bool
bool G4bool
Definition:
G4Types.hh:86
G4int
int G4int
Definition:
G4Types.hh:85
G4VCrossSectionDataSet.hh
Z
const G4int Z[17]
Definition:
G4WaterStopping.cc:51
G4DynamicParticle
Definition:
G4DynamicParticle.hh:65
G4DynamicParticle::GetKineticEnergy
G4double GetKineticEnergy() const
G4Material
Definition:
G4Material.hh:118
G4NeutronElectronElXsc
Definition:
G4NeutronElectronElXsc.hh:46
G4NeutronElectronElXsc::fCutEnergy
G4double fCutEnergy
Definition:
G4NeutronElectronElXsc.hh:84
G4NeutronElectronElXsc::fee2
G4double fee2
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::~G4NeutronElectronElXsc
~G4NeutronElectronElXsc()
Definition:
G4NeutronElectronElXsc.cc:87
G4NeutronElectronElXsc::GetCutEnergy
G4double GetCutEnergy()
Definition:
G4NeutronElectronElXsc.hh:75
G4NeutronElectronElXsc::fAm
G4double fAm
Definition:
G4NeutronElectronElXsc.hh:82
G4NeutronElectronElXsc::fCofXsc
G4double fCofXsc
Definition:
G4NeutronElectronElXsc.hh:81
G4NeutronElectronElXsc::CalculateAm
G4double CalculateAm(G4double momentum)
Definition:
G4NeutronElectronElXsc.hh:97
G4NeutronElectronElXsc::fme
G4double fme
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::GetRosenbluthXsc
G4double GetRosenbluthXsc(const G4DynamicParticle *, G4int Z, const G4Material *)
Definition:
G4NeutronElectronElXsc.cc:163
G4NeutronElectronElXsc::GetAm
G4double GetAm()
Definition:
G4NeutronElectronElXsc.hh:72
G4NeutronElectronElXsc::fXscArray
static const G4double fXscArray[200]
Definition:
G4NeutronElectronElXsc.hh:88
G4NeutronElectronElXsc::fM
G4double fM
Definition:
G4NeutronElectronElXsc.hh:77
G4NeutronElectronElXsc::fMv2
G4double fMv2
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::fme2
G4double fme2
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::GetElementNonRelXsc
G4double GetElementNonRelXsc(const G4DynamicParticle *, G4int Z, const G4Material *)
Definition:
G4NeutronElectronElXsc.hh:113
G4NeutronElectronElXsc::fBiasingFactor
G4double fBiasingFactor
Definition:
G4NeutronElectronElXsc.hh:85
G4NeutronElectronElXsc::fEnergyXscVector
G4PhysicsLogVector * fEnergyXscVector
Definition:
G4NeutronElectronElXsc.hh:87
G4NeutronElectronElXsc::XscIntegrand
G4double XscIntegrand(G4double x)
Definition:
G4NeutronElectronElXsc.cc:190
G4NeutronElectronElXsc::fEnergyBin
G4int fEnergyBin
Definition:
G4NeutronElectronElXsc.hh:83
G4NeutronElectronElXsc::G4NeutronElectronElXsc
G4NeutronElectronElXsc()
Definition:
G4NeutronElectronElXsc.cc:50
G4NeutronElectronElXsc::fee
G4double fee
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::SetCutEnergy
void SetCutEnergy(G4double ec)
Definition:
G4NeutronElectronElXsc.hh:74
G4NeutronElectronElXsc::IsElementApplicable
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
Definition:
G4NeutronElectronElXsc.cc:101
G4NeutronElectronElXsc::fMaxEnergy
G4double fMaxEnergy
Definition:
G4NeutronElectronElXsc.hh:84
G4NeutronElectronElXsc::fMinEnergy
G4double fMinEnergy
Definition:
G4NeutronElectronElXsc.hh:84
G4NeutronElectronElXsc::fM2
G4double fM2
Definition:
G4NeutronElectronElXsc.hh:80
G4NeutronElectronElXsc::GetElementCrossSection
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
Definition:
G4NeutronElectronElXsc.cc:142
G4NeutronElectronElXsc::Initialise
void Initialise()
Definition:
G4NeutronElectronElXsc.cc:116
G4NeutronElectronElXsc::SetBiasingFactor
void SetBiasingFactor(G4double bf)
Definition:
G4NeutronElectronElXsc.hh:77
G4PhysicsLogVector
Definition:
G4PhysicsLogVector.hh:48
G4PhysicsTable
Definition:
G4PhysicsTable.hh:56
G4VCrossSectionDataSet
Definition:
G4VCrossSectionDataSet.hh:70
globals.hh
CLHEP::Bohr_radius
static constexpr double Bohr_radius
Definition:
PhysicalConstants.h:100
CLHEP::hbarc
static constexpr double hbarc
Definition:
PhysicalConstants.h:66
Generated by
1.9.3