Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4RPGProtonInelastic.cc
이 파일의 문서화 페이지로 가기
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: G4RPGProtonInelastic.cc 79697 2014-03-12 13:10:09Z gcosmo $
27 //
28 
29 #include "G4RPGProtonInelastic.hh"
30 #include "G4SystemOfUnits.hh"
31 #include "Randomize.hh"
32 
35  G4Nucleus& targetNucleus )
36 {
38  const G4HadProjectile *originalIncident = &aTrack;
39  if (originalIncident->GetKineticEnergy()<= 0.1)
40  {
44  return &theParticleChange;
45  }
46 
47  //
48  // create the target particle
49  //
50  G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
51 
52  if (originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9. )
53  {
54  SlowProton( originalIncident, targetNucleus );
55  delete originalTarget;
56  return &theParticleChange;
57  }
58 
59  // Fermi motion and evaporation
60  // As of Geant3, the Fermi energy calculation had not been Done
61 
62  G4double ek = originalIncident->GetKineticEnergy();
63  G4double amas = originalIncident->GetDefinition()->GetPDGMass();
64  G4ReactionProduct modifiedOriginal;
65  modifiedOriginal = *originalIncident;
66 
67  G4double tkin = targetNucleus.Cinema( ek );
68  ek += tkin;
69  modifiedOriginal.SetKineticEnergy(ek);
70  G4double et = ek + amas;
71  G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
72  G4double pp = modifiedOriginal.GetMomentum().mag();
73  if (pp > 0.0) {
74  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
75  modifiedOriginal.SetMomentum( momentum * (p/pp) );
76  }
77  //
78  // calculate black track energies
79  //
80  tkin = targetNucleus.EvaporationEffects(ek);
81  ek -= tkin;
82  modifiedOriginal.SetKineticEnergy(ek);
83  et = ek + amas;
84  p = std::sqrt( std::abs((et-amas)*(et+amas)) );
85  pp = modifiedOriginal.GetMomentum().mag();
86  if (pp > 0.0) {
87  G4ThreeVector momentum = modifiedOriginal.GetMomentum();
88  modifiedOriginal.SetMomentum( momentum * (p/pp) );
89  }
90  const G4double cutOff = 0.1;
91  if (modifiedOriginal.GetKineticEnergy() < cutOff) {
92  SlowProton( originalIncident, targetNucleus );
93  delete originalTarget;
94  return &theParticleChange;
95  }
96 
97  G4ReactionProduct currentParticle = modifiedOriginal;
98  G4ReactionProduct targetParticle;
99  targetParticle = *originalTarget;
100  currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
101  targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
102  G4bool incidentHasChanged = false;
103  G4bool targetHasChanged = false;
104  G4bool quasiElastic = false;
105  G4FastVector<G4ReactionProduct,256> vec; // vec will contain the sec. particles
106  G4int vecLen = 0;
107  vec.Initialize( 0 );
108 
109  InitialCollision(vec, vecLen, currentParticle, targetParticle,
110  incidentHasChanged, targetHasChanged);
111 
112  CalculateMomenta(vec, vecLen,
113  originalIncident, originalTarget, modifiedOriginal,
114  targetNucleus, currentParticle, targetParticle,
115  incidentHasChanged, targetHasChanged, quasiElastic);
116 
117  SetUpChange( vec, vecLen,
118  currentParticle, targetParticle,
119  incidentHasChanged );
120 
121  delete originalTarget;
122  return &theParticleChange;
123 }
124 
125 
126 void
128  G4Nucleus &targetNucleus )
129 {
130  const G4double A = targetNucleus.GetA_asInt(); // atomic weight
131  const G4double Z = targetNucleus.GetZ_asInt(); // atomic number
132  //
133  // calculate Q-value of reactions
134  //
135  G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
136  G4double massVec[9];
137  massVec[0] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
138  massVec[1] = 0.;
139  if (A > Z+1.0)
140  massVec[1] = targetNucleus.AtomicMass( A , Z+1.0 );
141  massVec[2] = theAtomicMass;
142  massVec[3] = 0.;
143  if (A > 1.0 && A-1.0 > Z)
144  massVec[3] = targetNucleus.AtomicMass( A-1.0, Z );
145  massVec[4] = 0.;
146  if (A > 2.0 && A-2.0 > Z)
147  massVec[4] = targetNucleus.AtomicMass( A-2.0, Z );
148  massVec[5] = 0.;
149  if (A > 3.0 && Z > 1.0 && A-3.0 > Z-1.0)
150  massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-1.0 );
151  massVec[6] = 0.;
152  if (A > 1.0 && A-1.0 > Z+1.0)
153  massVec[6] = targetNucleus.AtomicMass( A-1.0, Z+1.0 );
154  massVec[7] = massVec[3];
155  massVec[8] = 0.;
156  if (A > 1.0 && Z > 1.0)
157  massVec[8] = targetNucleus.AtomicMass( A-1.0, Z-1.0 );
158 
159  G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
160  G4int vecLen = 0;
161  vec.Initialize( 0 );
162 
163  twoBody.NuclearReaction( vec, vecLen, originalIncident,
164  targetNucleus, theAtomicMass, massVec );
165 
168 
169  G4DynamicParticle *pd;
170  for( G4int i=0; i<vecLen; ++i )
171  {
172  pd = new G4DynamicParticle();
173  pd->SetDefinition( vec[i]->GetDefinition() );
174  pd->SetMomentum( vec[i]->GetMomentum() );
176  delete vec[i];
177  }
178 }
179 
180 
181 // Initial Collision
182 // selects the particle types arising from the initial collision of
183 // the proton and target nucleon. Secondaries are assigned to forward
184 // and backward reaction hemispheres, but final state energies and
185 // momenta are not calculated here.
186 
187 void
189  G4int& vecLen,
190  G4ReactionProduct& currentParticle,
191  G4ReactionProduct& targetParticle,
192  G4bool& incidentHasChanged,
193  G4bool& targetHasChanged)
194 {
195  G4double KE = currentParticle.GetKineticEnergy()/GeV;
196 
197  G4int mult;
198  G4int partType;
199  std::vector<G4int> fsTypes;
200  G4int part1;
201  G4int part2;
202 
203  G4double testCharge;
204  G4double testBaryon;
205  G4double testStrange;
206 
207  // Get particle types according to incident and target types
208 
209  if (targetParticle.GetDefinition() == particleDef[pro]) {
210  mult = GetMultiplicityT1(KE);
211  fsTypes = GetFSPartTypesForPP(mult, KE);
212 
213  part1 = fsTypes[0];
214  part2 = fsTypes[1];
215  currentParticle.SetDefinition(particleDef[part1]);
216  targetParticle.SetDefinition(particleDef[part2]);
217  if (part1 == pro) {
218  if (part2 == neu) {
219  if (G4UniformRand() > 0.5) {
220  incidentHasChanged = true;
221  targetParticle.SetDefinition(particleDef[part1]);
222  currentParticle.SetDefinition(particleDef[part2]);
223  } else {
224  targetHasChanged = true;
225  }
226  } else if (part2 > neu && part2 < xi0) {
227  targetHasChanged = true;
228  }
229 
230  } else { // neutron
231  targetHasChanged = true;
232  incidentHasChanged = true;
233  }
234 
235  testCharge = 2.0;
236  testBaryon = 2.0;
237  testStrange = 0.0;
238 
239  } else { // target was a neutron
240  mult = GetMultiplicityT0(KE);
241  fsTypes = GetFSPartTypesForPN(mult, KE);
242 
243  part1 = fsTypes[0];
244  part2 = fsTypes[1];
245  currentParticle.SetDefinition(particleDef[part1]);
246  targetParticle.SetDefinition(particleDef[part2]);
247  if (part1 == pro) {
248  if (part2 == pro) {
249  targetHasChanged = true;
250  } else if (part2 == neu) {
251  if (G4UniformRand() > 0.5) {
252  incidentHasChanged = true;
253  targetHasChanged = true;
254  targetParticle.SetDefinition(particleDef[part1]);
255  currentParticle.SetDefinition(particleDef[part2]);
256  }
257  } else { // hyperon
258  targetHasChanged = true;
259  }
260 
261  } else { // neutron
262  incidentHasChanged = true;
263  if (part2 > neu && part2 < xi0) targetHasChanged = true;
264  }
265 
266  testCharge = 1.0;
267  testBaryon = 2.0;
268  testStrange = 0.0;
269  }
270 
271  // Remove incident and target from fsTypes
272 
273  fsTypes.erase(fsTypes.begin());
274  fsTypes.erase(fsTypes.begin());
275 
276  // Remaining particles are secondaries. Put them into vec.
277 
278  G4ReactionProduct* rp(0);
279  for(G4int i=0; i < mult-2; ++i ) {
280  partType = fsTypes[i];
281  rp = new G4ReactionProduct();
282  rp->SetDefinition(particleDef[partType]);
283  (G4UniformRand() < 0.5) ? rp->SetSide(-1) : rp->SetSide(1);
284  vec.SetElement(vecLen++, rp);
285  }
286 
287  // Check conservation of charge, strangeness, baryon number
288 
289  CheckQnums(vec, vecLen, currentParticle, targetParticle,
290  testCharge, testBaryon, testStrange);
291 
292  return;
293 }
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
G4double ek
G4RPGTwoBody twoBody
std::vector< G4int > GetFSPartTypesForPP(G4int mult, G4double KE) const
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:382
void SetMomentumChange(const G4ThreeVector &aV)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4double x, const G4double y, const G4double z)
const char * p
Definition: xmltok.h:285
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * GetDefinition() const
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)
G4double GetPDGMass() const
void SlowProton(const G4HadProjectile *originalIncident, G4Nucleus &targetNucleus)
void CheckQnums(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4double Q, G4double B, G4double S)
G4DynamicParticle * ReturnTargetParticle() const
Definition: G4Nucleus.cc:241
Float_t Z
G4int GetMultiplicityT0(G4double KE) const
void SetEnergyChange(G4double anEnergy)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetKineticEnergy(const G4double en)
G4double GetKineticEnergy() const
G4int GetMultiplicityT1(G4double KE) const
void NuclearReaction(G4FastVector< G4ReactionProduct, 4 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4Nucleus &aNucleus, const G4double theAtomicMass, const G4double *massVec)
#define G4UniformRand()
Definition: Randomize.hh:53
double A(double temperature)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
const G4LorentzVector & Get4Momentum() const
Hep3Vector unit() const
G4double GetKineticEnergy() const
G4ParticleDefinition * particleDef[18]
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
const G4ParticleDefinition * GetDefinition() const
double mag() const
int G4int
Definition: G4Types.hh:78
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
void SetMomentum(const G4ThreeVector &momentum)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
Hep3Vector vect() const
void InitialCollision(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged)
G4double EvaporationEffects(G4double kineticEnergy)
Definition: G4Nucleus.cc:278
static constexpr double GeV
Definition: G4SIunits.hh:217
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
std::vector< G4int > GetFSPartTypesForPN(G4int mult, G4double KE) const
void SetStatusChange(G4HadFinalStateStatus aS)
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:254