Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LCapture.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//
27// $Id$
28//
29//
30// G4 Model: Low-energy Neutron Capture
31// F.W. Jones, TRIUMF, 03-DEC-96
32//
33// This is a prototype of a low-energy neutron capture process.
34// Currently it is based on the GHEISHA routine CAPTUR,
35// and conforms fairly closely to the original Fortran.
36//
37// HPW Capture using models now. the code comes from the
38// original G4LCapture class.
39//
40// 25-JUN-98 FWJ: replaced missing Initialize for ParticleChange.
41//
42
43#include <iostream>
44
45#include "globals.hh"
46#include "Randomize.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4LCapture.hh"
50
53{
54 SetMinEnergy(0.0*GeV);
56 // Description();
57}
58
59
61{
63}
64
65
66void G4LCapture::ModelDescription(std::ostream& outFile) const
67{
68 outFile << "G4LCapture is one of the Low Energy Parameterized\n"
69 << "(LEP) models used to implement neutron capture on nuclei.\n"
70 << "It is a re-engineered version of the GHEISHA code of\n"
71 << "H. Fesefeldt which simply adds the neutron mass and energy\n"
72 << "to the target nucleus, and emits gammas isotropically as\n"
73 << "long as there is sufficient excitation energy in the\n"
74 << "daughter nucleus. The model is applicable to all incident\n"
75 << "neutron energies.\n";
76}
77
78
81{
84
85 G4double N = targetNucleus.GetA_asInt();
86 G4double Z = targetNucleus.GetZ_asInt();
87
88 const G4LorentzVector theMom = aTrack.Get4Momentum();
89 G4double P = theMom.vect().mag()/GeV;
90 G4double Px = theMom.vect().x()/GeV;
91 G4double Py = theMom.vect().y()/GeV;
92 G4double Pz = theMom.vect().z()/GeV;
93 G4double E = theMom.e()/GeV;
94 G4double E0 = aTrack.GetDefinition()->GetPDGMass()/GeV;
95 G4double Q = aTrack.GetDefinition()->GetPDGCharge();
96 if (verboseLevel > 1) {
97 G4cout << "G4LCapture:ApplyYourself: incident particle:" << G4endl;
98 G4cout << "P " << P << " GeV/c" << G4endl;
99 G4cout << "Px " << Px << " GeV/c" << G4endl;
100 G4cout << "Py " << Py << " GeV/c" << G4endl;
101 G4cout << "Pz " << Pz << " GeV/c" << G4endl;
102 G4cout << "E " << E << " GeV" << G4endl;
103 G4cout << "mass " << E0 << " GeV" << G4endl;
104 G4cout << "charge " << Q << G4endl;
105 }
106
107 // GHEISHA ADD operation to get total energy, mass, charge:
108
109 if (verboseLevel > 1) {
110 G4cout << "G4LCapture:ApplyYourself: material:" << G4endl;
111 G4cout << "A " << N << G4endl;
112 G4cout << "Z " << Z << G4endl;
113 G4cout << "atomic mass " <<
114 Atomas(N, Z) << "GeV" << G4endl;
115 }
116 E = E + Atomas(N, Z);
117 G4double E02 = E*E - P*P;
118 E0 = std::sqrt(std::abs(E02));
119 if (E02 < 0) E0 = -E0;
120 Q = Q + Z;
121 if (verboseLevel > 1) {
122 G4cout << "G4LCapture:ApplyYourself: total:" << G4endl;
123 G4cout << "E " << E << " GeV" << G4endl;
124 G4cout << "mass " << E0 << " GeV" << G4endl;
125 G4cout << "charge " << Q << G4endl;
126 }
127 Px = -Px;
128 Py = -Py;
129 Pz = -Pz;
130
131 // Make a gamma...
132
133 G4double p;
134 if (Z == 1 && N == 1) { // special case for hydrogen
135 p = 0.0022;
136 } else {
137 G4double ran = G4RandGauss::shoot();
138 p = 0.0065 + ran*0.0010;
139 }
140
141 G4double ran1 = G4UniformRand();
142 G4double ran2 = G4UniformRand();
143 G4double cost = -1. + 2.*ran1;
144 G4double sint = std::sqrt(std::abs(1. - cost*cost));
145 G4double phi = ran2*twopi;
146
147 G4double px = p*sint*std::sin(phi);
148 G4double py = p*sint*std::cos(phi);
149 G4double pz = p*cost;
150 G4double e = p;
151
152 G4double a = px*Px + py*Py + pz*Pz;
153 a = (a/(E + E0) - e)/E0;
154
155 px = px + a*Px;
156 py = py + a*Py;
157 pz = pz + a*Pz;
158
159 G4DynamicParticle* aGamma;
161 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
163
164 // Make another gamma if there is sufficient energy left over...
165
166 G4double xp = 0.008 - p;
167 if (xp > 0.) {
168 if (Z > 1 || N > 1) {
169 ran1 = G4UniformRand();
170 ran2 = G4UniformRand();
171 cost = -1. + 2.*ran1;
172 sint = std::sqrt(std::abs(1. - cost*cost));
173 phi = ran2*twopi;
174
175 px = xp*sint*std::sin(phi);
176 py = xp*sint*std::cos(phi);
177 pz = xp*cost;
178 e = xp;
179
180 a = px*Px + py*Py + pz*Pz;
181 a = (a/(E + E0) - e)/E0;
182
183 px = px + a*Px;
184 py = py + a*Py;
185 pz = pz + a*Pz;
186
188 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
190 }
191 }
192 return &theParticleChange;
193}
194
195const std::pair<G4double, G4double> G4LCapture::GetFatalEnergyCheckLevels() const
196{
197 // max energy non-conservation is mass of heavy nucleus
198 return std::pair<G4double, G4double>(5*perCent,250*GeV);
199}
@ stopAndKill
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double y() const
double mag() const
Hep3Vector vect() const
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
virtual void ModelDescription(std::ostream &outFile) const
Definition: G4LCapture.cc:66
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
Definition: G4LCapture.cc:195
G4LCapture(const G4String &name="G4LCapture")
Definition: G4LCapture.cc:51
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
Definition: G4LCapture.cc:80
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83