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