Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4RPGNeutronInelastic.cc 92494 2015-09-02 07:19:42Z gcosmo $
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) { /* Loop checking, 01.09.2015, D.Wright */
139  cost1 = -1.0 + 2.0*G4UniformRand();
140  eka = 1.0 + 2.0*cost1*A + A*A;
141  }
142  G4double cost = std::min( 1.0, std::max( -1.0, (A*cost1+1.0)/std::sqrt(eka) ) );
143  eka /= (1.0+A)*(1.0+A);
144  G4double ek = currentKinetic*MeV/GeV;
145  G4double amas = currentMass*MeV/GeV;
146  ek *= eka;
147  G4double en = ek + amas;
148  G4double p = std::sqrt(std::abs(en*en-amas*amas));
149  G4double sint = std::sqrt(std::abs(1.0-cost*cost));
150  G4double phi = G4UniformRand()*twopi;
151  G4double px = sint*std::sin(phi);
152  G4double py = sint*std::cos(phi);
153  G4double pz = cost;
154  targetParticle.SetMomentum( px*GeV, py*GeV, pz*GeV );
155  G4double pxO = originalIncident->Get4Momentum().x()/GeV;
156  G4double pyO = originalIncident->Get4Momentum().y()/GeV;
157  G4double pzO = originalIncident->Get4Momentum().z()/GeV;
158  G4double ptO = pxO*pxO + pyO+pyO;
159  if( ptO > 0.0 )
160  {
161  G4double pO = std::sqrt(pxO*pxO+pyO*pyO+pzO*pzO);
162  cost = pzO/pO;
163  sint = 0.5*(std::sqrt(std::abs((1.0-cost)*(1.0+cost)))+std::sqrt(ptO)/pO);
164  G4double ph = pi/2.0;
165  if( pyO < 0.0 )ph = ph*1.5;
166  if( std::abs(pxO) > 0.000001 )ph = std::atan2(pyO,pxO);
167  G4double cosp = std::cos(ph);
168  G4double sinp = std::sin(ph);
169  px = cost*cosp*px - sinp*py+sint*cosp*pz;
170  py = cost*sinp*px + cosp*py+sint*sinp*pz;
171  pz = -sint*px + cost*pz;
172  }
173  else
174  {
175  if( pz < 0.0 )pz *= -1.0;
176  }
177  G4double pu = std::sqrt(px*px+py*py+pz*pz);
178  modifiedOriginal.SetMomentum( targetParticle.GetMomentum() * (p/pu) );
179  modifiedOriginal.SetKineticEnergy( ek*GeV );
180 
181  targetParticle.SetMomentum(
182  originalIncident->Get4Momentum().vect() - modifiedOriginal.GetMomentum() );
183  G4double pp = targetParticle.GetMomentum().mag();
184  G4double tarmas = targetParticle.GetMass();
185  targetParticle.SetTotalEnergy( std::sqrt( pp*pp + tarmas*tarmas ) );
186 
187  theParticleChange.SetEnergyChange( modifiedOriginal.GetKineticEnergy() );
189  pd->SetDefinition( targetParticle.GetDefinition() );
190  pd->SetMomentum( targetParticle.GetMomentum() );
192  return;
193  }
194 
195  G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
196  G4int vecLen = 0;
197  vec.Initialize( 0 );
198 
199  G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
200  G4double massVec[9];
201  massVec[0] = targetNucleus.AtomicMass( A+1.0, Z );
202  massVec[1] = theAtomicMass;
203  massVec[2] = 0.;
204  if (Z > 1.0) massVec[2] = targetNucleus.AtomicMass(A, Z-1.0);
205  massVec[3] = 0.;
206  if (Z > 1.0 && A > 1.0) massVec[3] = targetNucleus.AtomicMass(A-1.0, Z-1.0 );
207  massVec[4] = 0.;
208  if (Z > 1.0 && A > 2.0 && A-2.0 > Z-1.0)
209  massVec[4] = targetNucleus.AtomicMass( A-2.0, Z-1.0 );
210  massVec[5] = 0.;
211  if (Z > 2.0 && A > 3.0 && A-3.0 > Z-2.0)
212  massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-2.0 );
213  massVec[6] = 0.;
214  if (A > 1.0 && A-1.0 > Z) massVec[6] = targetNucleus.AtomicMass(A-1.0, Z);
215  massVec[7] = massVec[3];
216  massVec[8] = 0.;
217  if (Z > 2.0 && A > 1.0) massVec[8] = targetNucleus.AtomicMass( A-1.0,Z-2.0 );
218 
219  twoBody.NuclearReaction(vec, vecLen, originalIncident,
220  targetNucleus, theAtomicMass, massVec );
221 
224 
225  G4DynamicParticle* pd;
226  for( G4int i=0; i<vecLen; ++i ) {
227  pd = new G4DynamicParticle();
228  pd->SetDefinition( vec[i]->GetDefinition() );
229  pd->SetMomentum( vec[i]->GetMomentum() );
231  delete vec[i];
232  }
233 }
234 
235 
236 // Initial Collision
237 // selects the particle types arising from the initial collision of
238 // the neutron and target nucleon. Secondaries are assigned to
239 // forward and backward reaction hemispheres, but final state energies
240 // and momenta are not calculated here.
241 
242 void
243 G4RPGNeutronInelastic::InitialCollision(G4FastVector<G4ReactionProduct,256>& vec,
244  G4int& vecLen,
245  G4ReactionProduct& currentParticle,
246  G4ReactionProduct& targetParticle,
247  G4bool& incidentHasChanged,
248  G4bool& targetHasChanged)
249 {
250  G4double KE = currentParticle.GetKineticEnergy()/GeV;
251 
252  G4int mult;
253  G4int partType;
254  std::vector<G4int> fsTypes;
255  G4int part1;
256  G4int part2;
257 
258  G4double testCharge;
259  G4double testBaryon;
260  G4double testStrange;
261 
262  // Get particle types according to incident and target types
263 
264  if (targetParticle.GetDefinition() == particleDef[neu]) {
265  mult = GetMultiplicityT1(KE);
266  fsTypes = GetFSPartTypesForNN(mult, KE);
267 
268  part1 = fsTypes[0];
269  part2 = fsTypes[1];
270  currentParticle.SetDefinition(particleDef[part1]);
271  targetParticle.SetDefinition(particleDef[part2]);
272  if (part1 == pro) {
273  if (part2 == neu) {
274  if (G4UniformRand() > 0.5) {
275  incidentHasChanged = true;
276  } else {
277  targetHasChanged = true;
278  currentParticle.SetDefinition(particleDef[part2]);
279  targetParticle.SetDefinition(particleDef[part1]);
280  }
281  } else {
282  targetHasChanged = true;
283  incidentHasChanged = true;
284  }
285 
286  } else { // neutron
287  if (part2 > neu && part2 < xi0) targetHasChanged = true;
288  }
289 
290  testCharge = 0.0;
291  testBaryon = 2.0;
292  testStrange = 0.0;
293 
294  } else { // target was a proton
295  mult = GetMultiplicityT0(KE);
296  fsTypes = GetFSPartTypesForNP(mult, KE);
297 
298  part1 = fsTypes[0];
299  part2 = fsTypes[1];
300  currentParticle.SetDefinition(particleDef[part1]);
301  targetParticle.SetDefinition(particleDef[part2]);
302  if (part1 == pro) {
303  if (part2 == pro) {
304  incidentHasChanged = true;
305  } else if (part2 == neu) {
306  if (G4UniformRand() > 0.5) {
307  incidentHasChanged = true;
308  targetHasChanged = true;
309  } else {
310  currentParticle.SetDefinition(particleDef[part2]);
311  targetParticle.SetDefinition(particleDef[part1]);
312  }
313 
314  } else if (part2 > neu && part2 < xi0) {
315  incidentHasChanged = true;
316  targetHasChanged = true;
317  }
318 
319  } else { // neutron
320  targetHasChanged = true;
321  }
322 
323  testCharge = 1.0;
324  testBaryon = 2.0;
325  testStrange = 0.0;
326  }
327 
328  // if (mult == 2 && !incidentHasChanged && !targetHasChanged)
329  // quasiElastic = true;
330 
331  // Remove incident and target from fsTypes
332 
333  fsTypes.erase(fsTypes.begin());
334  fsTypes.erase(fsTypes.begin());
335 
336  // Remaining particles are secondaries. Put them into vec.
337 
338  G4ReactionProduct* rp(0);
339  for(G4int i=0; i < mult-2; ++i ) {
340  partType = fsTypes[i];
341  rp = new G4ReactionProduct();
342  rp->SetDefinition(particleDef[partType]);
343  (G4UniformRand() < 0.5) ? rp->SetSide(-1) : rp->SetSide(1);
344  vec.SetElement(vecLen++, rp);
345  }
346 
347  // Check conservation of charge, strangeness, baryon number
348 
349  CheckQnums(vec, vecLen, currentParticle, targetParticle,
350  testCharge, testBaryon, testStrange);
351 
352  return;
353 }
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:254
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
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:241
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static constexpr double twopi
Definition: G4SIunits.hh:76
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetMultiplicityT1(G4double KE) const
const G4ParticleDefinition * GetDefinition() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
const G4ParticleDefinition * GetDefinition() const
G4ParticleDefinition * particleDef[18]
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
void SetTotalEnergy(const G4double en)
G4double ek
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:382
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
static constexpr double GeV
Definition: G4SIunits.hh:217
G4ThreeVector GetMomentum() const
static constexpr double MeV
Definition: G4SIunits.hh:214
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)
static constexpr double pi
Definition: G4SIunits.hh:75
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4int GetMultiplicityT0(G4double KE) const
double mag() const
G4double GetMass() const
G4RPGTwoBody twoBody