source: trunk/source/processes/hadronic/models/low_energy/src/G4LCapture.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 5.8 KB
Line 
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: G4LCapture.cc,v 1.14 2007/02/24 05:17:29 dennis Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31// G4 Model: Low-energy Neutron Capture
32// F.W. Jones, TRIUMF, 03-DEC-96
33//
34// This is a prototype of a low-energy neutron capture process.
35// Currently it is based on the GHEISHA routine CAPTUR,
36// and conforms fairly closely to the original Fortran.
37//
38// HPW Capture using models now. the code comes from the
39// original G4LCapture class.
40//
41// 25-JUN-98 FWJ: replaced missing Initialize for ParticleChange.
42//
43
44#include "globals.hh"
45#include "G4LCapture.hh"
46#include "Randomize.hh"
47
48G4LCapture::G4LCapture() :
49 G4HadronicInteraction("G4LCapture")
50{
51 SetMinEnergy( 0.0*GeV );
52 SetMaxEnergy( DBL_MAX );
53}
54
55G4LCapture::~G4LCapture()
56{
57 theParticleChange.Clear();
58}
59
60G4HadFinalState*
61G4LCapture::ApplyYourself(const G4HadProjectile & aTrack, G4Nucleus& targetNucleus)
62{
63
64 theParticleChange.Clear();
65 theParticleChange.SetStatusChange(stopAndKill);
66
67 G4double N = targetNucleus.GetN();
68 G4double Z = targetNucleus.GetZ();
69
70 const G4LorentzVector theMom = aTrack.Get4Momentum();
71 G4double P = theMom.vect().mag()/GeV;
72 G4double Px = theMom.vect().x()/GeV;
73 G4double Py = theMom.vect().y()/GeV;
74 G4double Pz = theMom.vect().z()/GeV;
75 G4double E = theMom.e()/GeV;
76 G4double E0 = aTrack.GetDefinition()->GetPDGMass()/GeV;
77 G4double Q = aTrack.GetDefinition()->GetPDGCharge();
78 if (verboseLevel > 1) {
79 G4cout << "G4LCapture:ApplyYourself: incident particle:" << G4endl;
80 G4cout << "P " << P << " GeV/c" << G4endl;
81 G4cout << "Px " << Px << " GeV/c" << G4endl;
82 G4cout << "Py " << Py << " GeV/c" << G4endl;
83 G4cout << "Pz " << Pz << " GeV/c" << G4endl;
84 G4cout << "E " << E << " GeV" << G4endl;
85 G4cout << "mass " << E0 << " GeV" << G4endl;
86 G4cout << "charge " << Q << G4endl;
87 }
88
89 // GHEISHA ADD operation to get total energy, mass, charge:
90
91 if (verboseLevel > 1) {
92 G4cout << "G4LCapture:ApplyYourself: material:" << G4endl;
93 G4cout << "A " << N << G4endl;
94 G4cout << "Z " << Z << G4endl;
95 G4cout << "atomic mass " <<
96 Atomas(N, Z) << "GeV" << G4endl;
97 }
98 E = E + Atomas(N, Z);
99 G4double E02 = E*E - P*P;
100 E0 = std::sqrt(std::abs(E02));
101 if (E02 < 0) E0 = -E0;
102 Q = Q + Z;
103 if (verboseLevel > 1) {
104 G4cout << "G4LCapture:ApplyYourself: total:" << G4endl;
105 G4cout << "E " << E << " GeV" << G4endl;
106 G4cout << "mass " << E0 << " GeV" << G4endl;
107 G4cout << "charge " << Q << G4endl;
108 }
109 Px = -Px;
110 Py = -Py;
111 Pz = -Pz;
112
113 // Make a gamma...
114
115 G4double p;
116 if (Z == 1 && N == 1) { // special case for hydrogen
117 p = 0.0022;
118 } else {
119 G4double ran = G4RandGauss::shoot();
120 p = 0.0065 + ran*0.0010;
121 }
122
123 G4double ran1 = G4UniformRand();
124 G4double ran2 = G4UniformRand();
125 G4double cost = -1. + 2.*ran1;
126 G4double sint = std::sqrt(std::abs(1. - cost*cost));
127 G4double phi = ran2*twopi;
128
129 G4double px = p*sint*std::sin(phi);
130 G4double py = p*sint*std::cos(phi);
131 G4double pz = p*cost;
132 G4double e = p;
133 G4double e0 = 0.;
134
135 G4double a = px*Px + py*Py + pz*Pz;
136 a = (a/(E + E0) - e)/E0;
137
138 px = px + a*Px;
139 py = py + a*Py;
140 pz = pz + a*Pz;
141
142 G4DynamicParticle* aGamma;
143 aGamma = new G4DynamicParticle(G4Gamma::GammaDefinition(),
144 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
145 theParticleChange.AddSecondary(aGamma);
146
147 // Make another gamma if there is sufficient energy left over...
148
149 G4double xp = 0.008 - p;
150 if (xp > 0.) {
151 if (Z > 1 || N > 1) {
152 ran1 = G4UniformRand();
153 ran2 = G4UniformRand();
154 cost = -1. + 2.*ran1;
155 sint = std::sqrt(std::abs(1. - cost*cost));
156 phi = ran2*twopi;
157
158 px = xp*sint*std::sin(phi);
159 py = xp*sint*std::cos(phi);
160 pz = xp*cost;
161 e = xp;
162 e0 = 0.;
163
164 a = px*Px + py*Py + pz*Pz;
165 a = (a/(E + E0) - e)/E0;
166
167 px = px + a*Px;
168 py = py + a*Py;
169 pz = pz + a*Pz;
170
171 aGamma = new G4DynamicParticle(G4Gamma::GammaDefinition(),
172 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
173 theParticleChange.AddSecondary(aGamma);
174 }
175 }
176 return &theParticleChange;
177}
Note: See TracBrowser for help on using the repository browser.