Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// Geant4 Header : G4HadronElastic
28//
29// Author : V.Ivanchenko 29 June 2009 (redesign old elastic model)
30//
31
32#include "G4HadronElastic.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4ParticleTable.hh"
36#include "G4IonTable.hh"
37#include "Randomize.hh"
38#include "G4Proton.hh"
39#include "G4Neutron.hh"
40#include "G4Deuteron.hh"
41#include "G4Alpha.hh"
42#include "G4Pow.hh"
43#include "G4Exp.hh"
44#include "G4Log.hh"
46
47
50{
51 SetMinEnergy( 0.0*GeV );
53 lowestEnergyLimit= 1.e-6*eV;
54 pLocalTmax = 0.0;
55 nwarn = 0;
56
57 theProton = G4Proton::Proton();
58 theNeutron = G4Neutron::Neutron();
59 theDeuteron = G4Deuteron::Deuteron();
60 theAlpha = G4Alpha::Alpha();
61}
62
64{}
65
66
67void G4HadronElastic::ModelDescription(std::ostream& outFile) const
68{
69 outFile << "G4HadronElastic is the base class for all hadron-nucleus\n"
70 << "elastic scattering models except HP.\n"
71 << "By default it 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 but fit the data only for relativistic scattering.\n";
76}
77
79 const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
80{
82
83 const G4HadProjectile* aParticle = &aTrack;
84 G4double ekin = aParticle->GetKineticEnergy();
85
86 // no scattering below the limit
87 if(ekin <= lowestEnergyLimit) {
90 return &theParticleChange;
91 }
92
93 G4int A = targetNucleus.GetA_asInt();
94 G4int Z = targetNucleus.GetZ_asInt();
95
96 // Scattered particle referred to axis of incident particle
97 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
98 G4double m1 = theParticle->GetPDGMass();
99 G4double plab = std::sqrt(ekin*(ekin + 2.0*m1));
100
101 if (verboseLevel>1) {
102 G4cout << "G4HadronElastic: "
103 << aParticle->GetDefinition()->GetParticleName()
104 << " Plab(GeV/c)= " << plab/GeV
105 << " Ekin(MeV) = " << ekin/MeV
106 << " scattered off Z= " << Z
107 << " A= " << A
108 << G4endl;
109 }
110
112 G4double e1 = m1 + ekin;
113 G4LorentzVector lv(0.0,0.0,plab,e1+mass2);
114 G4ThreeVector bst = lv.boostVector();
115 G4double momentumCMS = plab*mass2/std::sqrt(m1*m1 + mass2*mass2 + 2.*mass2*e1);
116
117 pLocalTmax = 4.0*momentumCMS*momentumCMS;
118
119 // Sampling in CM system
120 G4double t = SampleInvariantT(theParticle, plab, Z, A);
121
122 if(t < 0.0 || t > pLocalTmax) {
123 // For the very rare cases where cos(theta) is greater than 1 or smaller than -1,
124 // print some debugging information via a "JustWarning" exception, and resample
125 // using the default algorithm
126#ifdef G4VERBOSE
127 if(nwarn < 2) {
129 ed << GetModelName() << " wrong sampling t= " << t << " tmax= " << pLocalTmax
130 << " for " << aParticle->GetDefinition()->GetParticleName()
131 << " ekin=" << ekin << " MeV"
132 << " off (Z,A)=(" << Z << "," << A << ") - will be resampled" << G4endl;
133 G4Exception( "G4HadronElastic::ApplyYourself", "hadEla001", JustWarning, ed);
134 ++nwarn;
135 }
136#endif
137 t = G4HadronElastic::SampleInvariantT(theParticle, plab, Z, A);
138 }
139
140 G4double phi = G4UniformRand()*CLHEP::twopi;
141 G4double cost = 1. - 2.0*t/pLocalTmax;
142
143 if (cost > 1.0) { cost = 1.0; }
144 else if(cost < -1.0) { cost = -1.0; }
145
146 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
147
148 if (verboseLevel>1) {
149 G4cout << " t= " << t << " tmax(GeV^2)= " << pLocalTmax/(GeV*GeV)
150 << " Pcms(GeV)= " << momentumCMS/GeV << " cos(t)=" << cost
151 << " sin(t)=" << sint << G4endl;
152 }
153 G4LorentzVector nlv1(momentumCMS*sint*std::cos(phi),
154 momentumCMS*sint*std::sin(phi),
155 momentumCMS*cost,
156 std::sqrt(momentumCMS*momentumCMS + m1*m1));
157
158 nlv1.boost(bst);
159
160 G4double eFinal = nlv1.e() - m1;
161 if (verboseLevel > 1) {
162 G4cout <<"G4HadronElastic: m= " << m1 << " Efin(MeV)= " << eFinal
163 << " 4-M Final: " << nlv1
164 << G4endl;
165 }
166
167 if(eFinal <= 0.0) {
170 } else {
173 }
174 lv -= nlv1;
175 G4double erec = std::max(lv.e() - mass2, 0.0);
176 if (verboseLevel > 1) {
177 G4cout << "Recoil: " <<" m= " << mass2 << " Erec(MeV)= " << erec
178 << " 4-mom: " << lv
179 << G4endl;
180 }
181
182 // the recoil is created if kinetic energy above the threshold
183 if(erec > GetRecoilEnergyThreshold()) {
184 G4ParticleDefinition * theDef = nullptr;
185 if(Z == 1 && A == 1) { theDef = theProton; }
186 else if (Z == 1 && A == 2) { theDef = theDeuteron; }
187 else if (Z == 1 && A == 3) { theDef = G4Triton::Triton(); }
188 else if (Z == 2 && A == 3) { theDef = G4He3::He3(); }
189 else if (Z == 2 && A == 4) { theDef = theAlpha; }
190 else {
191 theDef =
193 }
194 G4DynamicParticle * aSec = new G4DynamicParticle(theDef, lv.vect().unit(), erec);
196 } else {
198 }
199
200 return &theParticleChange;
201}
202
203// sample momentum transfer in the CMS system
206 G4double mom, G4int, G4int A)
207{
208 const G4double plabLowLimit = 400.0*CLHEP::MeV;
209 const G4double GeV2 = GeV*GeV;
210 const G4double z07in13 = std::pow(0.7, 0.3333333333);
211 const G4double numLimit = 18.;
212
213 G4int pdg = std::abs(part->GetPDGEncoding());
214 G4double tmax = pLocalTmax/GeV2;
215
216 G4double aa, bb, cc, dd;
217 G4Pow* g4pow = G4Pow::GetInstance();
218 if (A <= 62) {
219 if (pdg == 211){ //Pions
220 if(mom >= plabLowLimit){ //High energy
221 bb = 14.5*g4pow->Z23(A);/*14.5*/
222 dd = 10.;
223 cc = 0.075*g4pow->Z13(A)/dd;//1.4
224 //aa = g4pow->powZ(A, 1.93)/bb;//1.63
225 aa = (A*A)/bb;//1.63
226 } else { //Low energy
227 bb = 29.*z07in13*z07in13*g4pow->Z23(A);
228 dd = 15.;
229 cc = 0.04*g4pow->Z13(A)/dd;//1.4
230 aa = g4pow->powZ(A, 1.63)/bb;//1.63
231 }
232 } else { //Other particles
233 bb = 14.5*g4pow->Z23(A);
234 dd = 20.;
235 aa = (A*A)/bb;//1.63
236 cc = 1.4*g4pow->Z13(A)/dd;
237 }
238 //===========================
239 } else { //(A>62)
240 if (pdg == 211) {
241 if(mom >= plabLowLimit){ //high
242 bb = 60.*z07in13*g4pow->Z13(A);//60
243 dd = 30.;
244 aa = 0.5*(A*A)/bb;//1.33
245 cc = 4.*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
246 } else { //low
247 bb = 120.*z07in13*g4pow->Z13(A);//60
248 dd = 30.;
249 aa = 2.*g4pow->powZ(A,1.33)/bb;
250 cc = 4.*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
251 }
252 } else {
253 bb = 60.*g4pow->Z13(A);
254 dd = 25.;
255 aa = g4pow->powZ(A,1.33)/bb;//1.33
256 cc = 0.2*g4pow->powZ(A,0.4)/dd;//1:0.4 --- 2: 0.4
257 }
258 }
259 G4double q1 = 1.0 - G4Exp(-std::min(bb*tmax, numLimit));
260 G4double q2 = 1.0 - G4Exp(-std::min(dd*tmax, numLimit));
261 G4double s1 = q1*aa;
262 G4double s2 = q2*cc;
263 if((s1 + s2)*G4UniformRand() < s2) {
264 q1 = q2;
265 bb = dd;
266 }
267 return -GeV2*G4Log(1.0 - G4UniformRand()*q1)/bb;
268}
269
270//////////////////////////////////////////////
271//
272// Cofs for s-,c-,b-particles ds/dt slopes
273
275{
276 // The input parameter "pdg" should be the absolute value of the PDG code
277 // (i.e. the same value for a particle and its antiparticle).
278
279 G4double coeff = 1.0;
280
281 // heavy barions
282
283 static const G4double lBarCof1S = 0.88;
284 static const G4double lBarCof2S = 0.76;
285 static const G4double lBarCof3S = 0.64;
286 static const G4double lBarCof1C = 0.784378;
287 static const G4double lBarCofSC = 0.664378;
288 static const G4double lBarCof2SC = 0.544378;
289 static const G4double lBarCof1B = 0.740659;
290 static const G4double lBarCofSB = 0.620659;
291 static const G4double lBarCof2SB = 0.500659;
292
293 if( pdg == 3122 || pdg == 3222 || pdg == 3112 || pdg == 3212 )
294 {
295 coeff = lBarCof1S; // Lambda, Sigma+, Sigma-, Sigma0
296
297 } else if( pdg == 3322 || pdg == 3312 )
298 {
299 coeff = lBarCof2S; // Xi-, Xi0
300 }
301 else if( pdg == 3324)
302 {
303 coeff = lBarCof3S; // Omega
304 }
305 else if( pdg == 4122 || pdg == 4212 || pdg == 4222 || pdg == 4112 )
306 {
307 coeff = lBarCof1C; // LambdaC+, SigmaC+, SigmaC++, SigmaC0
308 }
309 else if( pdg == 4332 )
310 {
311 coeff = lBarCof2SC; // OmegaC
312 }
313 else if( pdg == 4232 || pdg == 4132 )
314 {
315 coeff = lBarCofSC; // XiC+, XiC0
316 }
317 else if( pdg == 5122 || pdg == 5222 || pdg == 5112 || pdg == 5212 )
318 {
319 coeff = lBarCof1B; // LambdaB, SigmaB+, SigmaB-, SigmaB0
320 }
321 else if( pdg == 5332 )
322 {
323 coeff = lBarCof2SB; // OmegaB-
324 }
325 else if( pdg == 5132 || pdg == 5232 ) // XiB-, XiB0
326 {
327 coeff = lBarCofSB;
328 }
329 // heavy mesons Kaons?
330 static const G4double lMesCof1S = 0.82; // Kp/piP kaons?
331 static const G4double llMesCof1C = 0.676568;
332 static const G4double llMesCof1B = 0.610989;
333 static const G4double llMesCof2C = 0.353135;
334 static const G4double llMesCof2B = 0.221978;
335 static const G4double llMesCofSC = 0.496568;
336 static const G4double llMesCofSB = 0.430989;
337 static const G4double llMesCofCB = 0.287557;
338 static const G4double llMesCofEtaP = 0.88;
339 static const G4double llMesCofEta = 0.76;
340
341 if( pdg == 321 || pdg == 311 || pdg == 310 )
342 {
343 coeff = lMesCof1S; //K+-0
344 }
345 else if( pdg == 511 || pdg == 521 )
346 {
347 coeff = llMesCof1B; // BMeson0, BMeson+
348 }
349 else if(pdg == 421 || pdg == 411 )
350 {
351 coeff = llMesCof1C; // DMeson+, DMeson0
352 }
353 else if( pdg == 531 )
354 {
355 coeff = llMesCofSB; // BSMeson0
356 }
357 else if( pdg == 541 )
358 {
359 coeff = llMesCofCB; // BCMeson+-
360 }
361 else if(pdg == 431 )
362 {
363 coeff = llMesCofSC; // DSMeson+-
364 }
365 else if(pdg == 441 || pdg == 443 )
366 {
367 coeff = llMesCof2C; // Etac, JPsi
368 }
369 else if(pdg == 553 )
370 {
371 coeff = llMesCof2B; // Upsilon
372 }
373 else if(pdg == 221 )
374 {
375 coeff = llMesCofEta; // Eta
376 }
377 else if(pdg == 331 )
378 {
379 coeff = llMesCofEtaP; // Eta'
380 }
381 return coeff;
382}
383
384
double A(double temperature)
const G4double GeV2
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
void SetLocalEnergyDeposit(G4double aE)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
~G4HadronElastic() override
void ModelDescription(std::ostream &) const override
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A) override
G4double GetSlopeCof(const G4int pdg)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
G4HadronElastic(const G4String &name="hElasticLHEP")
void SetMinEnergy(G4double anEnergy)
G4double GetRecoilEnergyThreshold() const
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:94