Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RPGPiPlusInelastic.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 "G4RPGPiPlusInelastic.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "Randomize.hh"
32 
35  G4Nucleus& targetNucleus)
36 {
37  const G4HadProjectile *originalIncident = &aTrack;
38  if (originalIncident->GetKineticEnergy()<= 0.1) {
42  return &theParticleChange;
43  }
44 
45  // create the target particle
46 
47  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
48  G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
49 
50  G4ReactionProduct currentParticle(
51  const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
52  currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
53  currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
54 
55  // Fermi motion and evaporation
56  // As of Geant3, the Fermi energy calculation had not been Done
57 
58  G4double ek = originalIncident->GetKineticEnergy();
59  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
60 
61  G4double tkin = targetNucleus.Cinema( ek );
62  ek += tkin;
63  currentParticle.SetKineticEnergy( ek );
64  G4double et = ek + amas;
65  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
66  G4double pp = currentParticle.GetMomentum().mag();
67  if( pp > 0.0 )
68  {
69  G4ThreeVector momentum = currentParticle.GetMomentum();
70  currentParticle.SetMomentum( momentum * (p/pp) );
71  }
72 
73  // calculate black track energies
74 
75  tkin = targetNucleus.EvaporationEffects( ek );
76  ek -= tkin;
77  currentParticle.SetKineticEnergy( ek );
78  et = ek + amas;
79  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
80  pp = currentParticle.GetMomentum().mag();
81  if( pp > 0.0 )
82  {
83  G4ThreeVector momentum = currentParticle.GetMomentum();
84  currentParticle.SetMomentum( momentum * (p/pp) );
85  }
86 
87  G4ReactionProduct modifiedOriginal = currentParticle;
88 
89  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
90  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
91  G4bool incidentHasChanged = false;
92  G4bool targetHasChanged = false;
93  G4bool quasiElastic = false;
94  G4FastVector<G4ReactionProduct,256> vec; // vec will contain the secondary particles
95  G4int vecLen = 0;
96  vec.Initialize( 0 );
97 
98  const G4double cutOff = 0.1;
99  if( currentParticle.GetKineticEnergy() > cutOff )
100  InitialCollision(vec, vecLen, currentParticle, targetParticle,
101  incidentHasChanged, targetHasChanged);
102 
103  CalculateMomenta( vec, vecLen,
104  originalIncident, originalTarget, modifiedOriginal,
105  targetNucleus, currentParticle, targetParticle,
106  incidentHasChanged, targetHasChanged, quasiElastic );
107 
108  SetUpChange( vec, vecLen,
109  currentParticle, targetParticle,
110  incidentHasChanged );
111 
112  delete originalTarget;
113  return &theParticleChange;
114 }
115 
116 
117 // Initial Collision
118 // selects the particle types arising from the initial collision of
119 // the projectile and target nucleon. Secondaries are assigned to
120 // forward and backward reaction hemispheres, but final state energies
121 // and momenta are not calculated here.
122 
123 void
124 G4RPGPiPlusInelastic::InitialCollision(G4FastVector<G4ReactionProduct,256>& vec,
125  G4int& vecLen,
126  G4ReactionProduct& currentParticle,
127  G4ReactionProduct& targetParticle,
128  G4bool& incidentHasChanged,
129  G4bool& targetHasChanged)
130 {
131  G4double KE = currentParticle.GetKineticEnergy()/GeV;
132 
133  G4int mult;
134  G4int partType;
135  std::vector<G4int> fsTypes;
136 
137  G4double testCharge;
138  G4double testBaryon;
139  G4double testStrange;
140 
141  // Get particle types according to incident and target types
142 
143  if (targetParticle.GetDefinition() == particleDef[pro]) {
144  mult = GetMultiplicityT32(KE);
145  fsTypes = GetFSPartTypesForPipP(mult, KE);
146  partType = fsTypes[0];
147  if (partType != pro) {
148  targetHasChanged = true;
149  targetParticle.SetDefinition(particleDef[partType]);
150  }
151 
152  testCharge = 2.0;
153  testBaryon = 1.0;
154  testStrange = 0.0;
155 
156  } else { // target was a neutron
157  mult = GetMultiplicityT12(KE);
158  fsTypes = GetFSPartTypesForPipN(mult, KE);
159  partType = fsTypes[0];
160  if (partType != neu) {
161  targetHasChanged = true;
162  targetParticle.SetDefinition(particleDef[partType]);
163  }
164 
165  testCharge = 1.0;
166  testBaryon = 1.0;
167  testStrange = 0.0;
168  }
169 
170  // Remove target particle from list
171 
172  fsTypes.erase(fsTypes.begin());
173 
174  // See if the incident particle changed type
175 
176  G4int choose = -1;
177  for(G4int i=0; i < mult-1; ++i ) {
178  partType = fsTypes[i];
179  if (partType == pip) {
180  choose = i;
181  break;
182  }
183  }
184  if (choose == -1) {
185  incidentHasChanged = true;
186  choose = G4int(G4UniformRand()*(mult-1) );
187  partType = fsTypes[choose];
188  currentParticle.SetDefinition(particleDef[partType]);
189  }
190  fsTypes.erase(fsTypes.begin()+choose);
191 
192  // Remaining particles are secondaries. Put them into vec.
193  // Improve this by randomizing secondary order, then alternate
194  // which secondary is put into forward or backward hemisphere
195 
196  G4ReactionProduct* rp(0);
197  for(G4int i=0; i < mult-2; ++i ) {
198  partType = fsTypes[i];
199  rp = new G4ReactionProduct();
200  rp->SetDefinition(particleDef[partType]);
201  (G4UniformRand() < 0.5) ? rp->SetSide(-1) : rp->SetSide(1);
202  if (partType > pim && partType < pro) rp->SetMayBeKilled(false); // kaons
203  vec.SetElement(vecLen++, rp);
204  }
205 
206  // if (mult == 2 && !incidentHasChanged && !targetHasChanged)
207  // quasiElastic = true;
208 
209  // Check conservation of charge, strangeness, baryon number
210 
211  CheckQnums(vec, vecLen, currentParticle, targetParticle,
212  testCharge, testBaryon, testStrange);
213 
214  return;
215 }
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:264
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
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)
G4ParticleDefinition * GetDefinition() const
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)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ParticleDefinition * GetDefinition() const
G4ParticleDefinition * particleDef[18]
std::vector< G4int > GetFSPartTypesForPipP(G4int mult, G4double KE) const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetKineticEnergy() const
void SetEnergyChange(G4double anEnergy)
std::vector< G4int > GetFSPartTypesForPipN(G4int mult, G4double KE) const
G4double GetPDGMass() const
Hep3Vector unit() const
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4int GetMultiplicityT12(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 SetMomentumChange(const G4ThreeVector &aV)
G4int GetMultiplicityT32(G4double KE) const