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: G4VHadronElastic.cc,v 1.4 2009/09/23 14:37:44 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03-ref-09 $ |
---|
28 | // |
---|
29 | // Geant4 Header : G4VHadronElastic |
---|
30 | // |
---|
31 | // Author : V.Ivanchenko 29 June 2009 (redesign old elastic model) |
---|
32 | // |
---|
33 | // Modified: |
---|
34 | // |
---|
35 | // |
---|
36 | |
---|
37 | #include "G4VHadronElastic.hh" |
---|
38 | #include "G4ParticleTable.hh" |
---|
39 | #include "G4ParticleDefinition.hh" |
---|
40 | #include "G4IonTable.hh" |
---|
41 | #include "Randomize.hh" |
---|
42 | #include "G4Proton.hh" |
---|
43 | #include "G4Neutron.hh" |
---|
44 | #include "G4Deuteron.hh" |
---|
45 | #include "G4Alpha.hh" |
---|
46 | #include "G4Pow.hh" |
---|
47 | |
---|
48 | G4VHadronElastic::G4VHadronElastic(const G4String& name) |
---|
49 | : G4HadronicInteraction(name) |
---|
50 | { |
---|
51 | SetMinEnergy( 0.0*GeV ); |
---|
52 | SetMaxEnergy( 100.*TeV ); |
---|
53 | lowestEnergyLimit= 1.e-6*eV; |
---|
54 | |
---|
55 | theProton = G4Proton::Proton(); |
---|
56 | theNeutron = G4Neutron::Neutron(); |
---|
57 | theDeuteron = G4Deuteron::Deuteron(); |
---|
58 | theAlpha = G4Alpha::Alpha(); |
---|
59 | } |
---|
60 | |
---|
61 | G4VHadronElastic::~G4VHadronElastic() |
---|
62 | {} |
---|
63 | |
---|
64 | G4HadFinalState* G4VHadronElastic::ApplyYourself( |
---|
65 | const G4HadProjectile& aTrack, G4Nucleus& targetNucleus) |
---|
66 | { |
---|
67 | theParticleChange.Clear(); |
---|
68 | |
---|
69 | const G4HadProjectile* aParticle = &aTrack; |
---|
70 | G4double ekin = aParticle->GetKineticEnergy(); |
---|
71 | if(ekin <= lowestEnergyLimit) { |
---|
72 | theParticleChange.SetEnergyChange(ekin); |
---|
73 | theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); |
---|
74 | return &theParticleChange; |
---|
75 | } |
---|
76 | |
---|
77 | G4double aTarget = targetNucleus.GetN(); |
---|
78 | G4double zTarget = targetNucleus.GetZ(); |
---|
79 | |
---|
80 | G4double plab = aParticle->GetTotalMomentum(); |
---|
81 | |
---|
82 | // Scattered particle referred to axis of incident particle |
---|
83 | const G4ParticleDefinition* theParticle = aParticle->GetDefinition(); |
---|
84 | G4double m1 = theParticle->GetPDGMass(); |
---|
85 | |
---|
86 | G4int Z = static_cast<G4int>(zTarget+0.5); |
---|
87 | G4int A = static_cast<G4int>(aTarget+0.5); |
---|
88 | if (verboseLevel>1) { |
---|
89 | G4cout << "G4VHadronElastic: " |
---|
90 | << aParticle->GetDefinition()->GetParticleName() |
---|
91 | << " Plab(GeV/c)= " << plab/GeV |
---|
92 | << " Ekin(MeV) = " << ekin/MeV |
---|
93 | << " scattered off Z= " << Z |
---|
94 | << " A= " << A |
---|
95 | << G4endl; |
---|
96 | } |
---|
97 | G4ParticleDefinition * theDef = 0; |
---|
98 | |
---|
99 | if(Z == 1 && A == 1) theDef = theProton; |
---|
100 | else if (Z == 1 && A == 2) theDef = theDeuteron; |
---|
101 | else if (Z == 1 && A == 3) theDef = G4Triton::Triton(); |
---|
102 | else if (Z == 2 && A == 3) theDef = G4He3::He3(); |
---|
103 | else if (Z == 2 && A == 4) theDef = theAlpha; |
---|
104 | else theDef = G4ParticleTable::GetParticleTable()->FindIon(Z,A,0,Z); |
---|
105 | |
---|
106 | G4double m2 = theDef->GetPDGMass(); |
---|
107 | G4LorentzVector lv1 = aParticle->Get4Momentum(); |
---|
108 | G4LorentzVector lv(0.0,0.0,0.0,m2); |
---|
109 | lv += lv1; |
---|
110 | |
---|
111 | G4ThreeVector bst = lv.boostVector(); |
---|
112 | lv1.boost(-bst); |
---|
113 | |
---|
114 | G4ThreeVector p1 = lv1.vect(); |
---|
115 | momentumCMS = p1.mag(); |
---|
116 | G4double tmax = 4.0*momentumCMS*momentumCMS; |
---|
117 | |
---|
118 | // Sampling in CM system |
---|
119 | G4double t = SampleInvariantT(theParticle, plab, Z, A); |
---|
120 | G4double phi = G4UniformRand()*CLHEP::twopi; |
---|
121 | G4double cost = 1. - 2.0*t/tmax; |
---|
122 | G4double sint; |
---|
123 | |
---|
124 | // problem in sampling |
---|
125 | if(cost > 1.0 || cost < -1.0) { |
---|
126 | if(verboseLevel > 0) { |
---|
127 | G4cout << "G4VHadronElastic WARNING cost= " << cost |
---|
128 | << " after scattering of " |
---|
129 | << aParticle->GetDefinition()->GetParticleName() |
---|
130 | << " p(GeV/c)= " << plab |
---|
131 | << " on " << theDef->GetParticleName() |
---|
132 | << G4endl; |
---|
133 | } |
---|
134 | cost = 1.0; |
---|
135 | sint = 0.0; |
---|
136 | |
---|
137 | // normal situation |
---|
138 | } else { |
---|
139 | sint = std::sqrt((1.0-cost)*(1.0+cost)); |
---|
140 | } |
---|
141 | if (verboseLevel>1) { |
---|
142 | G4cout << " t= " << t << " tmax= " << tmax |
---|
143 | << " Pcms= " << momentumCMS << "cos(t)=" << cost |
---|
144 | << " std::sin(t)=" << sint << G4endl; |
---|
145 | } |
---|
146 | G4ThreeVector v1(sint*std::cos(phi),sint*std::sin(phi),cost); |
---|
147 | v1 *= momentumCMS; |
---|
148 | G4LorentzVector nlv1(v1.x(),v1.y(),v1.z(), |
---|
149 | std::sqrt(momentumCMS*momentumCMS + m1*m1)); |
---|
150 | |
---|
151 | nlv1.boost(bst); |
---|
152 | |
---|
153 | G4double eFinal = nlv1.e() - m1; |
---|
154 | if (verboseLevel > 1) { |
---|
155 | G4cout <<" m= " << m1 << " Efin(MeV)= " << eFinal |
---|
156 | << " Proj: 4-mom " << lv1 << " Final: " << nlv1 |
---|
157 | << G4endl; |
---|
158 | } |
---|
159 | if(eFinal <= lowestEnergyLimit) { |
---|
160 | if(eFinal < 0.0 && verboseLevel > 0) { |
---|
161 | G4cout << "G4VHadronElastic WARNING Efinal= " << eFinal |
---|
162 | << " after scattering of " |
---|
163 | << aParticle->GetDefinition()->GetParticleName() |
---|
164 | << " p(GeV/c)= " << plab |
---|
165 | << " on " << theDef->GetParticleName() |
---|
166 | << G4endl; |
---|
167 | } |
---|
168 | theParticleChange.SetEnergyChange(0.0); |
---|
169 | nlv1 = G4LorentzVector(0.0,0.0,0.0,m1); |
---|
170 | |
---|
171 | } else { |
---|
172 | theParticleChange.SetMomentumChange(nlv1.vect().unit()); |
---|
173 | theParticleChange.SetEnergyChange(eFinal); |
---|
174 | } |
---|
175 | |
---|
176 | G4LorentzVector nlv0 = lv - nlv1; |
---|
177 | G4double erec = nlv0.e() - m2; |
---|
178 | if (verboseLevel > 1) { |
---|
179 | G4cout << "Recoil: " <<" m= " << m2 << " Erec(MeV)= " << erec |
---|
180 | << " 4-mom: " << nlv0 |
---|
181 | << G4endl; |
---|
182 | } |
---|
183 | |
---|
184 | if(erec > GetRecoilEnergyThreshold()) { |
---|
185 | G4DynamicParticle * aSec = new G4DynamicParticle(theDef, nlv0); |
---|
186 | theParticleChange.AddSecondary(aSec); |
---|
187 | } else if(erec > 0.0) { |
---|
188 | theParticleChange.SetLocalEnergyDeposit(erec); |
---|
189 | } |
---|
190 | |
---|
191 | return &theParticleChange; |
---|
192 | } |
---|
193 | |
---|
194 | // sample momentum transfer in the CMS system |
---|
195 | G4double |
---|
196 | G4VHadronElastic::SampleInvariantT(const G4ParticleDefinition* /*p*/, |
---|
197 | G4double /*ptot*/, |
---|
198 | G4int /*Z*/, G4int A) |
---|
199 | { |
---|
200 | static const G4double GeV2 = GeV*GeV; |
---|
201 | G4double tmax = 4.0*momentumCMS*momentumCMS/GeV2; |
---|
202 | G4double aa, bb, cc; |
---|
203 | G4double dd = 10.; |
---|
204 | G4Pow* p = G4Pow::GetInstance(); |
---|
205 | if (A <= 62) { |
---|
206 | bb = 14.5*p->Z23(A); |
---|
207 | aa = p->powZ(A, 1.63)/bb; |
---|
208 | cc = 1.4*p->Z13(A)/dd; |
---|
209 | } else { |
---|
210 | bb = 60.*p->Z13(A); |
---|
211 | aa = p->powZ(A, 1.33)/bb; |
---|
212 | cc = 0.4*p->powZ(A, 0.4)/dd; |
---|
213 | } |
---|
214 | G4double q1 = 1.0 - std::exp(-bb*tmax); |
---|
215 | G4double q2 = 1.0 - std::exp(-dd*tmax); |
---|
216 | G4double s1 = q1*aa; |
---|
217 | G4double s2 = q2*cc; |
---|
218 | if((s1 + s2)*G4UniformRand() < s2) { |
---|
219 | q1 = q2; |
---|
220 | bb = dd; |
---|
221 | } |
---|
222 | return -GeV2*std::log(1.0 - G4UniformRand()*q1)/bb; |
---|
223 | } |
---|
224 | |
---|