Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
hadronic
models
neutron_hp
src
G4NeutronHPIsotropic.cc
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_hp -- source file
27
// J.P. Wellisch, Nov-1996
28
// A prototype of the low energy neutron transport model.
29
//
30
#include "
G4NeutronHPIsotropic.hh
"
31
#include "
G4PhysicalConstants.hh
"
32
#include "
G4SystemOfUnits.hh
"
33
#include "
Randomize.hh
"
34
#include "
G4Gamma.hh
"
35
#include "
G4Electron.hh
"
36
#include "
G4Positron.hh
"
37
#include "
G4Neutron.hh
"
38
#include "
G4Proton.hh
"
39
#include "
G4Deuteron.hh
"
40
#include "
G4Triton.hh
"
41
#include "
G4He3.hh
"
42
#include "
G4Alpha.hh
"
43
#include "
G4ParticleTable.hh
"
44
45
46
void
G4NeutronHPIsotropic::Init
(std::ifstream & )
47
{
48
}
49
50
G4ReactionProduct
*
G4NeutronHPIsotropic::Sample
(
G4double
anEnergy,
G4double
massCode,
G4double
)
51
{
52
G4ReactionProduct
* result =
new
G4ReactionProduct
;
53
G4int
Z
=
static_cast<
G4int
>
(massCode/1000);
54
G4int
A =
static_cast<
G4int
>
(massCode-1000*
Z
);
55
56
if
(massCode==0)
57
{
58
result->
SetDefinition
(
G4Gamma::Gamma
());
59
}
60
else
if
(A==0)
61
{
62
result->
SetDefinition
(
G4Electron::Electron
());
63
if
(Z==1) result->
SetDefinition
(
G4Positron::Positron
());
64
}
65
else
if
(A==1)
66
{
67
result->
SetDefinition
(
G4Neutron::Neutron
());
68
if
(Z==1) result->
SetDefinition
(
G4Proton::Proton
());
69
}
70
else
if
(A==2)
71
{
72
result->
SetDefinition
(
G4Deuteron::Deuteron
());
73
}
74
else
if
(A==3)
75
{
76
result->
SetDefinition
(
G4Triton::Triton
());
77
if
(Z==2) result->
SetDefinition
(
G4He3::He3
());
78
}
79
else
if
(A==4)
80
{
81
result->
SetDefinition
(
G4Alpha::Alpha
());
82
//110607 TK modified following parts for migration to G4NDL3.15 (ENDF VII.r0)
83
//if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "Unknown ion case 1");
84
if
(Z!=2)
85
{
86
result->
SetDefinition
(
G4ParticleTable::GetParticleTable
()->GetIon ( Z , A , 0.0 ) );
87
}
88
}
89
else
90
{
91
//110607 TK modified following parts for migration to G4NDL3.15 (ENDF VII.r0)
92
result->
SetDefinition
(
G4ParticleTable::GetParticleTable
()->GetIon ( Z , A , 0.0 ) );
93
//throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPIsotropic: Unknown ion case 2");
94
}
95
96
G4double
cosTh =
G4UniformRand
();
97
G4double
phi =
twopi
*
G4UniformRand
();
98
G4double
theta = std::acos(cosTh);
99
G4double
sinth = std::sin(theta);
100
101
// we need the the Q value of the reaction
102
result->
SetKineticEnergy
(std::max(0.001*
MeV
, anEnergy+
GetQValue
()));
103
G4double
mtot = result->
GetTotalMomentum
();
104
G4ThreeVector
tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) );
105
result->
SetMomentum
(tempVector);
106
107
return
result;
108
}
Generated on Sat May 25 2013 14:34:04 for Geant4 by
1.8.4