source: trunk/source/processes/hadronic/models/rpg/src/G4RPGNeutronInelastic.cc@ 1103

Last change on this file since 1103 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 12.5 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// $Id: G4RPGNeutronInelastic.cc,v 1.4 2008/05/05 21:21:55 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29
30#include "G4RPGNeutronInelastic.hh"
31#include "Randomize.hh"
32
33
34G4HadFinalState*
35G4RPGNeutronInelastic::ApplyYourself(const G4HadProjectile& aTrack,
36 G4Nucleus& targetNucleus)
37{
38 theParticleChange.Clear();
39 const G4HadProjectile* originalIncident = &aTrack;
40
41 //
42 // create the target particle
43 //
44 G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
45
46 G4ReactionProduct modifiedOriginal;
47 modifiedOriginal = *originalIncident;
48 G4ReactionProduct targetParticle;
49 targetParticle = *originalTarget;
50 if( originalIncident->GetKineticEnergy()/GeV < 0.01 + 2.*G4UniformRand()/9. )
51 {
52 SlowNeutron(originalIncident,modifiedOriginal,targetParticle,targetNucleus );
53 delete originalTarget;
54 return &theParticleChange;
55 }
56
57 //
58 // Fermi motion and evaporation
59 // As of Geant3, the Fermi energy calculation had not been Done
60 //
61 G4double ek = originalIncident->GetKineticEnergy()/MeV;
62 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
63
64 G4double tkin = targetNucleus.Cinema( ek );
65 ek += tkin;
66 modifiedOriginal.SetKineticEnergy( ek*MeV );
67 G4double et = ek + amas;
68 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
69 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
70 if( pp > 0.0 )
71 {
72 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
73 modifiedOriginal.SetMomentum( momentum * (p/pp) );
74 }
75 //
76 // calculate black track energies
77 //
78 tkin = targetNucleus.EvaporationEffects( ek );
79 ek -= tkin;
80 modifiedOriginal.SetKineticEnergy(ek);
81 et = ek + amas;
82 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
83 pp = modifiedOriginal.GetMomentum().mag();
84 if( pp > 0.0 )
85 {
86 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
87 modifiedOriginal.SetMomentum( momentum * (p/pp) );
88 }
89 const G4double cutOff = 0.1;
90 if( modifiedOriginal.GetKineticEnergy()/MeV <= cutOff )
91 {
92 SlowNeutron( originalIncident, modifiedOriginal, targetParticle, targetNucleus );
93 delete originalTarget;
94 return &theParticleChange;
95 }
96
97 G4ReactionProduct currentParticle = modifiedOriginal;
98 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
99 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
100 G4bool incidentHasChanged = false;
101 G4bool targetHasChanged = false;
102 G4bool quasiElastic = false;
103 G4FastVector<G4ReactionProduct,256> vec; // vec will contain sec. particles
104 G4int vecLen = 0;
105 vec.Initialize( 0 );
106
107 InitialCollision(vec, vecLen, currentParticle, targetParticle,
108 incidentHasChanged, targetHasChanged);
109
110 CalculateMomenta(vec, vecLen,
111 originalIncident, originalTarget, modifiedOriginal,
112 targetNucleus, currentParticle, targetParticle,
113 incidentHasChanged, targetHasChanged, quasiElastic);
114
115 SetUpChange(vec, vecLen,
116 currentParticle, targetParticle,
117 incidentHasChanged);
118
119 delete originalTarget;
120 return &theParticleChange;
121}
122
123void
124G4RPGNeutronInelastic::SlowNeutron(const G4HadProjectile* originalIncident,
125 G4ReactionProduct& modifiedOriginal,
126 G4ReactionProduct& targetParticle,
127 G4Nucleus& targetNucleus)
128{
129 const G4double A = targetNucleus.GetN(); // atomic weight
130 const G4double Z = targetNucleus.GetZ(); // atomic number
131
132 G4double currentKinetic = modifiedOriginal.GetKineticEnergy()/MeV;
133 G4double currentMass = modifiedOriginal.GetMass()/MeV;
134 if( A < 1.5 ) // Hydrogen
135 {
136 //
137 // very simple simulation of scattering angle and energy
138 // nonrelativistic approximation with isotropic angular
139 // distribution in the cms system
140 //
141 G4double cost1, eka = 0.0;
142 while (eka <= 0.0)
143 {
144 cost1 = -1.0 + 2.0*G4UniformRand();
145 eka = 1.0 + 2.0*cost1*A + A*A;
146 }
147 G4double cost = std::min( 1.0, std::max( -1.0, (A*cost1+1.0)/std::sqrt(eka) ) );
148 eka /= (1.0+A)*(1.0+A);
149 G4double ek = currentKinetic*MeV/GeV;
150 G4double amas = currentMass*MeV/GeV;
151 ek *= eka;
152 G4double en = ek + amas;
153 G4double p = std::sqrt(std::abs(en*en-amas*amas));
154 G4double sint = std::sqrt(std::abs(1.0-cost*cost));
155 G4double phi = G4UniformRand()*twopi;
156 G4double px = sint*std::sin(phi);
157 G4double py = sint*std::cos(phi);
158 G4double pz = cost;
159 targetParticle.SetMomentum( px*GeV, py*GeV, pz*GeV );
160 G4double pxO = originalIncident->Get4Momentum().x()/GeV;
161 G4double pyO = originalIncident->Get4Momentum().y()/GeV;
162 G4double pzO = originalIncident->Get4Momentum().z()/GeV;
163 G4double ptO = pxO*pxO + pyO+pyO;
164 if( ptO > 0.0 )
165 {
166 G4double pO = std::sqrt(pxO*pxO+pyO*pyO+pzO*pzO);
167 cost = pzO/pO;
168 sint = 0.5*(std::sqrt(std::abs((1.0-cost)*(1.0+cost)))+std::sqrt(ptO)/pO);
169 G4double ph = pi/2.0;
170 if( pyO < 0.0 )ph = ph*1.5;
171 if( std::abs(pxO) > 0.000001 )ph = std::atan2(pyO,pxO);
172 G4double cosp = std::cos(ph);
173 G4double sinp = std::sin(ph);
174 px = cost*cosp*px - sinp*py+sint*cosp*pz;
175 py = cost*sinp*px + cosp*py+sint*sinp*pz;
176 pz = -sint*px + cost*pz;
177 }
178 else
179 {
180 if( pz < 0.0 )pz *= -1.0;
181 }
182 G4double pu = std::sqrt(px*px+py*py+pz*pz);
183 modifiedOriginal.SetMomentum( targetParticle.GetMomentum() * (p/pu) );
184 modifiedOriginal.SetKineticEnergy( ek*GeV );
185
186 targetParticle.SetMomentum(
187 originalIncident->Get4Momentum().vect() - modifiedOriginal.GetMomentum() );
188 G4double pp = targetParticle.GetMomentum().mag();
189 G4double tarmas = targetParticle.GetMass();
190 targetParticle.SetTotalEnergy( std::sqrt( pp*pp + tarmas*tarmas ) );
191
192 theParticleChange.SetEnergyChange( modifiedOriginal.GetKineticEnergy() );
193 G4DynamicParticle *pd = new G4DynamicParticle;
194 pd->SetDefinition( targetParticle.GetDefinition() );
195 pd->SetMomentum( targetParticle.GetMomentum() );
196 theParticleChange.AddSecondary( pd );
197 return;
198 }
199
200 G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
201 G4int vecLen = 0;
202 vec.Initialize( 0 );
203
204 G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
205 G4double massVec[9];
206 massVec[0] = targetNucleus.AtomicMass( A+1.0, Z );
207 massVec[1] = theAtomicMass;
208 massVec[2] = 0.;
209 if (Z > 1.0) massVec[2] = targetNucleus.AtomicMass(A, Z-1.0);
210 massVec[3] = 0.;
211 if (Z > 1.0 && A > 1.0) massVec[3] = targetNucleus.AtomicMass(A-1.0, Z-1.0 );
212 massVec[4] = 0.;
213 if (Z > 1.0 && A > 2.0 && A-2.0 > Z-1.0)
214 massVec[4] = targetNucleus.AtomicMass( A-2.0, Z-1.0 );
215 massVec[5] = 0.;
216 if (Z > 2.0 && A > 3.0 && A-3.0 > Z-2.0)
217 massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-2.0 );
218 massVec[6] = 0.;
219 if (A > 1.0 && A-1.0 > Z) massVec[6] = targetNucleus.AtomicMass(A-1.0, Z);
220 massVec[7] = massVec[3];
221 massVec[8] = 0.;
222 if (Z > 2.0 && A > 1.0) massVec[8] = targetNucleus.AtomicMass( A-1.0,Z-2.0 );
223
224 twoBody.NuclearReaction(vec, vecLen, originalIncident,
225 targetNucleus, theAtomicMass, massVec );
226
227 theParticleChange.SetStatusChange( stopAndKill );
228 theParticleChange.SetEnergyChange( 0.0 );
229
230 G4DynamicParticle* pd;
231 for( G4int i=0; i<vecLen; ++i ) {
232 pd = new G4DynamicParticle();
233 pd->SetDefinition( vec[i]->GetDefinition() );
234 pd->SetMomentum( vec[i]->GetMomentum() );
235 theParticleChange.AddSecondary( pd );
236 delete vec[i];
237 }
238}
239
240
241// Initial Collision
242// selects the particle types arising from the initial collision of
243// the neutron and target nucleon. Secondaries are assigned to
244// forward and backward reaction hemispheres, but final state energies
245// and momenta are not calculated here.
246
247void
248G4RPGNeutronInelastic::InitialCollision(G4FastVector<G4ReactionProduct,256>& vec,
249 G4int& vecLen,
250 G4ReactionProduct& currentParticle,
251 G4ReactionProduct& targetParticle,
252 G4bool& incidentHasChanged,
253 G4bool& targetHasChanged)
254{
255 G4double KE = currentParticle.GetKineticEnergy()/GeV;
256
257 G4int mult;
258 G4int partType;
259 std::vector<G4int> fsTypes;
260 G4int part1;
261 G4int part2;
262
263 G4double testCharge;
264 G4double testBaryon;
265 G4double testStrange;
266
267 // Get particle types according to incident and target types
268
269 if (targetParticle.GetDefinition() == particleDef[neu]) {
270 mult = GetMultiplicityT1(KE);
271 fsTypes = GetFSPartTypesForNN(mult, KE);
272
273 part1 = fsTypes[0];
274 part2 = fsTypes[1];
275 currentParticle.SetDefinition(particleDef[part1]);
276 targetParticle.SetDefinition(particleDef[part2]);
277 if (part1 == pro) {
278 if (part2 == neu) {
279 if (G4UniformRand() > 0.5) {
280 incidentHasChanged = true;
281 } else {
282 targetHasChanged = true;
283 currentParticle.SetDefinition(particleDef[part2]);
284 targetParticle.SetDefinition(particleDef[part1]);
285 }
286 } else {
287 targetHasChanged = true;
288 incidentHasChanged = true;
289 }
290
291 } else { // neutron
292 if (part2 > neu && part2 < xi0) targetHasChanged = true;
293 }
294
295 testCharge = 0.0;
296 testBaryon = 2.0;
297 testStrange = 0.0;
298
299 } else { // target was a proton
300 mult = GetMultiplicityT0(KE);
301 fsTypes = GetFSPartTypesForNP(mult, KE);
302
303 part1 = fsTypes[0];
304 part2 = fsTypes[1];
305 currentParticle.SetDefinition(particleDef[part1]);
306 targetParticle.SetDefinition(particleDef[part2]);
307 if (part1 == pro) {
308 if (part2 == pro) {
309 incidentHasChanged = true;
310 } else if (part2 == neu) {
311 if (G4UniformRand() > 0.5) {
312 incidentHasChanged = true;
313 targetHasChanged = true;
314 } else {
315 currentParticle.SetDefinition(particleDef[part2]);
316 targetParticle.SetDefinition(particleDef[part1]);
317 }
318
319 } else if (part2 > neu && part2 < xi0) {
320 incidentHasChanged = true;
321 targetHasChanged = true;
322 }
323
324 } else { // neutron
325 targetHasChanged = true;
326 }
327
328 testCharge = 1.0;
329 testBaryon = 2.0;
330 testStrange = 0.0;
331 }
332
333 // if (mult == 2 && !incidentHasChanged && !targetHasChanged)
334 // quasiElastic = true;
335
336 // Remove incident and target from fsTypes
337
338 fsTypes.erase(fsTypes.begin());
339 fsTypes.erase(fsTypes.begin());
340
341 // Remaining particles are secondaries. Put them into vec.
342
343 G4ReactionProduct* rp(0);
344 for(G4int i=0; i < mult-2; ++i ) {
345 partType = fsTypes[i];
346 rp = new G4ReactionProduct();
347 rp->SetDefinition(particleDef[partType]);
348 (G4UniformRand() < 0.5) ? rp->SetSide(-1) : rp->SetSide(1);
349 vec.SetElement(vecLen++, rp);
350 }
351
352 // Check conservation of charge, strangeness, baryon number
353
354 CheckQnums(vec, vecLen, currentParticle, targetParticle,
355 testCharge, testBaryon, testStrange);
356
357 return;
358}
Note: See TracBrowser for help on using the repository browser.