Geant4  10.02.p02
G4HadronElastic.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: G4HadronElastic.cc 91806 2015-08-06 12:20:45Z gcosmo $
27 //
28 // Geant4 Header : G4HadronElastic
29 //
30 // Author : V.Ivanchenko 29 June 2009 (redesign old elastic model)
31 //
32 
33 #include "G4HadronElastic.hh"
34 #include "G4SystemOfUnits.hh"
35 #include "G4ParticleTable.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4IonTable.hh"
38 #include "Randomize.hh"
39 #include "G4Proton.hh"
40 #include "G4Neutron.hh"
41 #include "G4Deuteron.hh"
42 #include "G4Alpha.hh"
43 #include "G4Pow.hh"
44 #include "G4Exp.hh"
45 #include "G4Log.hh"
46 
47 
49  : G4HadronicInteraction(name)
50 {
51  SetMinEnergy( 0.0*GeV );
52  SetMaxEnergy( 100.*TeV );
53  lowestEnergyLimit= 1.e-6*eV;
54 
59  //Description();
60 }
61 
62 
64 {}
65 
66 
67 void G4HadronElastic::ModelDescription(std::ostream& outFile) const
68 {
69 
70  outFile << "G4HadronElastic is a hadron-nucleus elastic scattering\n"
71  << "model which uses the Gheisha two-exponential momentum\n"
72  << "transfer parameterization. The model is fully relativistic\n"
73  << "as opposed to the original Gheisha model which was not.\n"
74  << "This model may be used for all long-lived hadrons at all\n"
75  << "incident energies.\n";
76 
77 }
78 
79 
81  const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
82 {
84 
85  const G4HadProjectile* aParticle = &aTrack;
86  G4double ekin = aParticle->GetKineticEnergy();
87  if(ekin <= lowestEnergyLimit) {
89  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
90  return &theParticleChange;
91  }
92 
93  G4int A = targetNucleus.GetA_asInt();
94  G4int Z = targetNucleus.GetZ_asInt();
95 
96  G4double plab = aParticle->GetTotalMomentum();
97 
98  // Scattered particle referred to axis of incident particle
99  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
100  G4double m1 = theParticle->GetPDGMass();
101 
102  if (verboseLevel>1) {
103  G4cout << "G4HadronElastic: "
104  << aParticle->GetDefinition()->GetParticleName()
105  << " Plab(GeV/c)= " << plab/GeV
106  << " Ekin(MeV) = " << ekin/MeV
107  << " scattered off Z= " << Z
108  << " A= " << A
109  << G4endl;
110  }
111 
113  G4LorentzVector lv1 = aParticle->Get4Momentum();
114  G4LorentzVector lv(0.0,0.0,0.0,mass2);
115  lv += lv1;
116 
117  G4ThreeVector bst = lv.boostVector();
118  lv1.boost(-bst);
119 
120  G4ThreeVector p1 = lv1.vect();
121  G4double momentumCMS = p1.mag();
122  G4double tmax = 4.0*momentumCMS*momentumCMS;
123 
124  // Sampling in CM system
125  G4double t = SampleInvariantT(theParticle, plab, Z, A);
127  G4double cost = 1. - 2.0*t/tmax;
128  G4double sint;
129 
130  // problem in sampling
131  if(cost > 1.0 || cost < -1.0) {
132  //if(verboseLevel > 0) {
133  G4cout << "G4HadronElastic WARNING (1 - cost)= " << 1 - cost
134  << " after scattering of "
135  << aParticle->GetDefinition()->GetParticleName()
136  << " p(GeV/c)= " << plab/GeV
137  << " on an ion Z= " << Z << " A= " << A
138  << G4endl;
139  //}
140  cost = 1.0;
141  sint = 0.0;
142 
143  // normal situation
144  } else {
145  sint = std::sqrt((1.0-cost)*(1.0+cost));
146  }
147  if (verboseLevel>1) {
148  G4cout << " t= " << t << " tmax(GeV^2)= " << tmax/(GeV*GeV)
149  << " Pcms(GeV)= " << momentumCMS/GeV << " cos(t)=" << cost
150  << " sin(t)=" << sint << G4endl;
151  }
152  G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost);
153  v1 *= momentumCMS;
154  G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(),
155  std::sqrt(momentumCMS*momentumCMS + m1*m1));
156 
157  nlv1.boost(bst);
158 
159  G4double eFinal = nlv1.e() - m1;
160  if (verboseLevel > 1) {
161  G4cout <<" m= " << m1 << " Efin(MeV)= " << eFinal
162  << " Proj: 4-mom " << lv1 << " Final: " << nlv1
163  << G4endl;
164  }
165  if(eFinal <= lowestEnergyLimit) {
166  if(eFinal < 0.0 && verboseLevel > 0) {
167  G4cout << "G4HadronElastic WARNING Efinal= " << eFinal
168  << " after scattering of "
169  << aParticle->GetDefinition()->GetParticleName()
170  << " p(GeV/c)= " << plab/GeV
171  << " on an ion Z= " << Z << " A= " << A
172  << G4endl;
173  }
175  nlv1 = G4LorentzVector(0.0,0.0,0.0,m1);
176 
177  } else {
178  theParticleChange.SetMomentumChange(nlv1.vect().unit());
180  }
181 
182  lv -= nlv1;
183  G4double erec = lv.e() - mass2;
184  if (verboseLevel > 1) {
185  G4cout << "Recoil: " <<" m= " << mass2 << " Erec(MeV)= " << erec
186  << " 4-mom: " << lv
187  << G4endl;
188  }
189 
190  if(erec > GetRecoilEnergyThreshold()) {
191  G4ParticleDefinition * theDef = 0;
192  if(Z == 1 && A == 1) { theDef = theProton; }
193  else if (Z == 1 && A == 2) { theDef = theDeuteron; }
194  else if (Z == 1 && A == 3) { theDef = G4Triton::Triton(); }
195  else if (Z == 2 && A == 3) { theDef = G4He3::He3(); }
196  else if (Z == 2 && A == 4) { theDef = theAlpha; }
197  else {
198  theDef =
200  }
201  G4DynamicParticle * aSec = new G4DynamicParticle(theDef, lv);
203  } else if(erec > 0.0) {
205  }
206 
207  return &theParticleChange;
208 }
209 
210 // sample momentum transfer in the CMS system
211 G4double
213  G4double plab,
214  G4int Z, G4int A)
215 {
216  static const G4double GeV2 = GeV*GeV;
217  G4double momentumCMS = ComputeMomentumCMS(p,plab,Z,A);
218  G4double tmax = 4.0*momentumCMS*momentumCMS/GeV2;
219  G4double aa, bb, cc;
220  G4double dd = 10.;
221  G4Pow* g4pow = G4Pow::GetInstance();
222  if (A <= 62) {
223  bb = 14.5*g4pow->Z23(A);
224  aa = g4pow->powZ(A, 1.63)/bb;
225  cc = 1.4*g4pow->Z13(A)/dd;
226  } else {
227  bb = 60.*g4pow->Z13(A);
228  aa = g4pow->powZ(A, 1.33)/bb;
229  cc = 0.4*g4pow->powZ(A, 0.4)/dd;
230  }
231  G4double q1 = 1.0 - G4Exp(-bb*tmax);
232  G4double q2 = 1.0 - G4Exp(-dd*tmax);
233  G4double s1 = q1*aa;
234  G4double s2 = q2*cc;
235  if((s1 + s2)*G4UniformRand() < s2) {
236  q1 = q2;
237  bb = dd;
238  }
239  return -GeV2*G4Log(1.0 - G4UniformRand()*q1)/bb;
240 }
241 
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual void ModelDescription(std::ostream &) const
virtual G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Pow.hh:56
G4HadronElastic(const G4String &name="hElasticLHEP")
G4String name
Definition: TRTMaterials.hh:40
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4ParticleDefinition * theNeutron
G4ParticleDefinition * theProton
int G4int
Definition: G4Types.hh:78
G4double GetRecoilEnergyThreshold() const
const G4String & GetParticleName() const
static ulg bb
Definition: csz_inflate.cc:344
void SetMinEnergy(G4double anEnergy)
G4ParticleDefinition * theAlpha
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
double A(double temperature)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static const double twopi
Definition: G4SIunits.hh:75
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:214
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4LorentzVector & Get4Momentum() const
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const double eV
Definition: G4SIunits.hh:212
void SetEnergyChange(G4double anEnergy)
G4double ComputeMomentumCMS(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
G4ParticleDefinition * theDeuteron
G4double lowestEnergyLimit
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double Z23(G4int Z) const
Definition: G4Pow.hh:154
virtual ~G4HadronElastic()
void SetMaxEnergy(const G4double anEnergy)
void SetLocalEnergyDeposit(G4double aE)
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:215
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:254
double G4double
Definition: G4Types.hh:76
static G4He3 * He3()
Definition: G4He3.cc:94
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetTotalMomentum() const
CLHEP::HepLorentzVector G4LorentzVector