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