Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LElastic.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 // Physics model class G4LElastic
27 // G4 Model: Low-energy Elastic scattering
28 // F.W. Jones, TRIUMF, 04-JUN-96
29 //
30 // use -scheme for elastic scattering: HPW, 20th June 1997
31 // most of the code comes from the old Low-energy Elastic class
32 //
33 // 25-JUN-98 FWJ: replaced missing Initialize for ParticleChange.
34 // 14-DEC-05 V.Ivanchenko: restore 1.19 version (7.0)
35 // 23-JAN-07 V.Ivanchenko: add protection inside sqrt
36 
37 #include <iostream>
38 
39 #include "G4LElastic.hh"
40 #include "globals.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "Randomize.hh"
44 #include "G4ParticleTable.hh"
45 #include "G4IonTable.hh"
46 
49 {
50  SetMinEnergy(0.0);
52 }
53 
54 
55 void G4LElastic::ModelDescription(std::ostream& outFile) const
56 {
57  outFile << "G4LElastic is one of the Low Energy Parameterized (LEP)\n"
58  << "models used to implement elastic hadron scattering from nuclei.\n"
59  << "It is a re-engineered version of the GHEISHA code of\n"
60  << "H. Fesefeldt. It performs simplified two-body elastic\n"
61  << "scattering for all long-lived hadronic projectiles by using\n"
62  << "a two-exponential parameterization in momentum transfer.\n"
63  << "It is valid for incident hadrons of all energies.\n";
64 }
65 
66 
68 G4LElastic::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
69 {
70  if(getenv("debug_LElastic")) verboseLevel = 5;
72  const G4HadProjectile* aParticle = &aTrack;
73  G4double atno2 = targetNucleus.GetA_asInt();
74  G4double zTarget = targetNucleus.GetZ_asInt();
77 
78  // Elastic scattering off Hydrogen
79 
80  G4DynamicParticle* aSecondary = 0;
81  if (atno2 < 1.5) {
82  const G4ParticleDefinition* aParticleType = aParticle->GetDefinition();
83  if (aParticleType == G4PionPlus::PionPlus())
84  aSecondary = LightMedia.PionPlusExchange(aParticle, targetNucleus);
85  else if (aParticleType == G4PionMinus::PionMinus())
86  aSecondary = LightMedia.PionMinusExchange(aParticle, targetNucleus);
87  else if (aParticleType == G4KaonPlus::KaonPlus())
88  aSecondary = LightMedia.KaonPlusExchange(aParticle, targetNucleus);
89  else if (aParticleType == G4KaonZeroShort::KaonZeroShort())
90  aSecondary = LightMedia.KaonZeroShortExchange(aParticle,targetNucleus);
91  else if (aParticleType == G4KaonZeroLong::KaonZeroLong())
92  aSecondary = LightMedia.KaonZeroLongExchange(aParticle, targetNucleus);
93  else if (aParticleType == G4KaonMinus::KaonMinus())
94  aSecondary = LightMedia.KaonMinusExchange(aParticle, targetNucleus);
95  else if (aParticleType == G4Proton::Proton())
96  aSecondary = LightMedia.ProtonExchange(aParticle, targetNucleus);
97  else if (aParticleType == G4AntiProton::AntiProton())
98  aSecondary = LightMedia.AntiProtonExchange(aParticle, targetNucleus);
99  else if (aParticleType == G4Neutron::Neutron())
100  aSecondary = LightMedia.NeutronExchange(aParticle, targetNucleus);
101  else if (aParticleType == G4AntiNeutron::AntiNeutron())
102  aSecondary = LightMedia.AntiNeutronExchange(aParticle, targetNucleus);
103  else if (aParticleType == G4Lambda::Lambda())
104  aSecondary = LightMedia.LambdaExchange(aParticle, targetNucleus);
105  else if (aParticleType == G4AntiLambda::AntiLambda())
106  aSecondary = LightMedia.AntiLambdaExchange(aParticle, targetNucleus);
107  else if (aParticleType == G4SigmaPlus::SigmaPlus())
108  aSecondary = LightMedia.SigmaPlusExchange(aParticle, targetNucleus);
109  else if (aParticleType == G4SigmaMinus::SigmaMinus())
110  aSecondary = LightMedia.SigmaMinusExchange(aParticle, targetNucleus);
111  else if (aParticleType == G4AntiSigmaPlus::AntiSigmaPlus())
112  aSecondary = LightMedia.AntiSigmaPlusExchange(aParticle,targetNucleus);
113  else if (aParticleType == G4AntiSigmaMinus::AntiSigmaMinus())
114  aSecondary= LightMedia.AntiSigmaMinusExchange(aParticle,targetNucleus);
115  else if (aParticleType == G4XiZero::XiZero())
116  aSecondary = LightMedia.XiZeroExchange(aParticle, targetNucleus);
117  else if (aParticleType == G4XiMinus::XiMinus())
118  aSecondary = LightMedia.XiMinusExchange(aParticle, targetNucleus);
119  else if (aParticleType == G4AntiXiZero::AntiXiZero())
120  aSecondary = LightMedia.AntiXiZeroExchange(aParticle, targetNucleus);
121  else if (aParticleType == G4AntiXiMinus::AntiXiMinus())
122  aSecondary = LightMedia.AntiXiMinusExchange(aParticle, targetNucleus);
123  else if (aParticleType == G4OmegaMinus::OmegaMinus())
124  aSecondary = LightMedia.OmegaMinusExchange(aParticle, targetNucleus);
125  else if (aParticleType == G4AntiOmegaMinus::AntiOmegaMinus())
126  aSecondary= LightMedia.AntiOmegaMinusExchange(aParticle,targetNucleus);
127  else if (aParticleType == G4KaonPlus::KaonPlus())
128  aSecondary = LightMedia.KaonPlusExchange(aParticle, targetNucleus);
129  }
130 
131  // Has a charge or strangeness exchange occurred?
132  if (aSecondary) {
133  aSecondary->SetMomentum(aParticle->Get4Momentum().vect());
135  theParticleChange.AddSecondary(aSecondary);
136  }
137  // G4cout << "Entering elastic scattering 1"<<G4endl;
138 
139  G4double p = aParticle->GetTotalMomentum()/GeV;
140  if (verboseLevel > 1)
141  G4cout << "G4LElastic::DoIt: Incident particle p=" << p << " GeV" << G4endl;
142 
143  if (p < 0.01) return &theParticleChange;
144 
145  // G4cout << "Entering elastic scattering 2"<<G4endl;
146 // Compute the direction of elastic scattering.
147 // It is planned to replace this code with a method based on
148 // parameterized functions and a Monte Carlo method to invert the CDF.
149 
150  G4double ran = G4UniformRand();
151  G4double aa, bb, cc, dd, rr;
152  if (atno2 <= 62.) {
153  aa = std::pow(atno2, 1.63);
154  bb = 14.5*std::pow(atno2, 0.66);
155  cc = 1.4*std::pow(atno2, 0.33);
156  dd = 10.;
157  }
158  else {
159  aa = std::pow(atno2, 1.33);
160  bb = 60.*std::pow(atno2, 0.33);
161  cc = 0.4*std::pow(atno2, 0.40);
162  dd = 10.;
163  }
164  aa = aa/bb;
165  cc = cc/dd;
166  rr = (aa + cc)*ran;
167  if (verboseLevel > 1) {
168  G4cout << "DoIt: aa,bb,cc,dd,rr" << G4endl;
169  G4cout << aa << " " << bb << " " << cc << " " << dd << " " << rr << G4endl;
170  }
171  G4double t1 = -std::log(ran)/bb;
172  G4double t2 = -std::log(ran)/dd;
173  if (verboseLevel > 1) {
174  G4cout << "t1,Fctcos " << t1 << " " << Fctcos(t1, aa, bb, cc, dd, rr) <<
175  G4endl;
176  G4cout << "t2,Fctcos " << t2 << " " << Fctcos(t2, aa, bb, cc, dd, rr) <<
177  G4endl;
178  }
179  G4double eps = 0.001;
180  G4int ind1 = 10;
181  G4double t;
182  G4int ier1;
183  ier1 = Rtmi(&t, t1, t2, eps, ind1,
184  aa, bb, cc, dd, rr);
185  if (verboseLevel > 1) {
186  G4cout << "From Rtmi, ier1=" << ier1 << G4endl;
187  G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) <<
188  G4endl;
189  }
190  if (ier1 != 0) t = 0.25*(3.*t1 + t2);
191  if (verboseLevel > 1) {
192  G4cout << "t, Fctcos " << t << " " << Fctcos(t, aa, bb, cc, dd, rr) <<
193  G4endl;
194  }
195  G4double phi = G4UniformRand()*twopi;
196  rr = 0.5*t/(p*p);
197  if (rr > 1.) rr = 0.;
198  if (verboseLevel > 1)
199  G4cout << "rr=" << rr << G4endl;
200  G4double cost = 1. - rr;
201  G4double sint = std::sqrt(std::max(rr*(2. - rr), 0.));
202  if (sint == 0.) return &theParticleChange;
203  // G4cout << "Entering elastic scattering 3"<<G4endl;
204  if (verboseLevel > 1)
205  G4cout << "cos(t)=" << cost << " std::sin(t)=" << sint << G4endl;
206 
207  // Scattered particle referred to axis of incident particle
208  G4double mass1 = aParticle->GetDefinition()->GetPDGMass();
209  G4int Z=static_cast<G4int>(zTarget+.5);
210  G4int A=static_cast<G4int>(atno2);
211  if(G4UniformRand()<atno2-A) A++;
212  //G4cout << " ion info "<<atno2 << " "<<A<<" "<<Z<<" "<<zTarget<<G4endl;
213  G4double mass2 =
215 
216  // non relativistic approximation
217  G4double a = 1 + mass2/mass1;
218  G4double b = -2.*p*cost;
219  G4double c = p*p*(1 - mass2/mass1);
220  G4double p1 = (-b+std::sqrt(std::max(0.0,b*b-4.*a*c)))/(2.*a);
221  G4double px = p1*sint*std::sin(phi);
222  G4double py = p1*sint*std::cos(phi);
223  G4double pz = p1*cost;
224 
225 // relativistic calculation
226  G4double etot = std::sqrt(mass1*mass1+p*p)+mass2;
227  a = etot*etot-p*p*cost*cost;
228  b = 2*p*p*(mass1*cost*cost-etot);
229  c = p*p*p*p*sint*sint;
230 
231  G4double de = (-b-std::sqrt(std::max(0.0,b*b-4.*a*c)))/(2.*a);
232  G4double e1 = std::sqrt(p*p+mass1*mass1)-de;
233  G4double p12=e1*e1-mass1*mass1;
234  p1 = std::sqrt(std::max(1.*eV*eV,p12));
235  px = p1*sint*std::sin(phi);
236  py = p1*sint*std::cos(phi);
237  pz = p1*cost;
238 
239  if (verboseLevel > 1)
240  {
241  G4cout << "Relevant test "<<p<<" "<<p1<<" "<<cost<<" "<<de<<G4endl;
242  G4cout << "p1/p = "<<p1/p<<" "<<mass1<<" "<<mass2<<" "<<a<<" "<<b<<" "<<c<<G4endl;
243  G4cout << "rest = "<< b*b<<" "<<4.*a*c<<" "<<G4endl;
244  G4cout << "make p1 = "<< p12<<" "<<e1*e1<<" "<<mass1*mass1<<" "<<G4endl;
245  }
246 // Incident particle
247  G4double pxinc = p*aParticle->Get4Momentum().vect().unit().x();
248  G4double pyinc = p*aParticle->Get4Momentum().vect().unit().y();
249  G4double pzinc = p*aParticle->Get4Momentum().vect().unit().z();
250  if (verboseLevel > 1) {
251  G4cout << "NOM SCAT " << px << " " << py << " " << pz << G4endl;
252  G4cout << "INCIDENT " << pxinc << " " << pyinc << " " << pzinc << G4endl;
253  }
254 
255 // Transform scattered particle to reflect direction of incident particle
256  G4double pxnew, pynew, pznew;
257  Defs1(p, px, py, pz, pxinc, pyinc, pzinc, &pxnew, &pynew, &pznew);
258 // Normalize:
259  G4double pxre=pxinc-pxnew;
260  G4double pyre=pyinc-pynew;
261  G4double pzre=pzinc-pznew;
262  G4ThreeVector it0(pxnew*GeV, pynew*GeV, pznew*GeV);
263  if(p1>0)
264  {
265  pxnew = pxnew/p1;
266  pynew = pynew/p1;
267  pznew = pznew/p1;
268  }
269  else
270  {
271  //G4double pphi = 2*pi*G4UniformRand();
272  //G4double ccth = 2*G4UniformRand()-1;
273  pxnew = 0;//std::sin(std::acos(ccth))*std::sin(pphi);
274  pynew = 0;//std::sin(std::acos(ccth))*std::cos(phi);
275  pznew = 1;//ccth;
276  }
277  if (verboseLevel > 1) {
278  G4cout << "DoIt: returning new momentum vector" << G4endl;
279  G4cout << "DoIt: "<<pxinc << " " << pyinc << " " << pzinc <<" "<<p<< G4endl;
280  G4cout << "DoIt: "<<pxnew << " " << pynew << " " << pznew <<" "<<p<< G4endl;
281  }
282 
283  if (aSecondary) {
284  aSecondary->SetMomentumDirection(pxnew, pynew, pznew);
285  } else {
286  try
287  {
288  theParticleChange.SetMomentumChange(pxnew, pynew, pznew);
289  theParticleChange.SetEnergyChange(std::sqrt(mass1*mass1+it0.mag2())-mass1);
290  }
291  catch(G4HadronicException)
292  {
293  std::cerr << "GHADException originating from components of G4LElastic"<<std::cout;
294  throw;
295  }
297  G4ThreeVector it(pxre*GeV, pyre*GeV, pzre*GeV);
298  G4DynamicParticle * aSec =
299  new G4DynamicParticle(theDef, it.unit(), it.mag2()/(2.*mass2));
301  // G4cout << "Final check ###### "<<p<<" "<<it.mag()<<" "<<p1<<G4endl;
302  }
303  return &theParticleChange;
304 }
305 
306 
307 // The following is a "translation" of a root-finding routine
308 // from GEANT3.21/GHEISHA. Some of the labelled block structure has
309 // been retained for clarity. This routine will not be needed after
310 // the planned revisions to DoIt().
311 
312 G4int
313 G4LElastic::Rtmi(G4double* x, G4double xli, G4double xri, G4double eps,
314  G4int iend,
315  G4double aa, G4double bb, G4double cc, G4double dd,
316  G4double rr)
317 {
318  G4int ier = 0;
319  G4double xl = xli;
320  G4double xr = xri;
321  *x = xl;
322  G4double tol = *x;
323  G4double f = Fctcos(tol, aa, bb, cc, dd, rr);
324  if (f == 0.) return ier;
325  G4double fl, fr;
326  fl = f;
327  *x = xr;
328  tol = *x;
329  f = Fctcos(tol, aa, bb, cc, dd, rr);
330  if (f == 0.) return ier;
331  fr = f;
332 
333 // Error return in case of wrong input data
334  if (fl*fr >= 0.) {
335  ier = 2;
336  return ier;
337  }
338 
339 // Basic assumption fl*fr less than 0 is satisfied.
340 // Generate tolerance for function values.
341  G4int i = 0;
342  G4double tolf = 100.*eps;
343 
344 // Start iteration loop
345 label4:
346  i++;
347 
348 // Start bisection loop
349  for (G4int k = 1; k <= iend; k++) {
350  *x = 0.5*(xl + xr);
351  tol = *x;
352  f = Fctcos(tol, aa, bb, cc, dd, rr);
353  if (f == 0.) return 0;
354  if (f*fr < 0.) { // Interchange xl and xr in order to get the
355  tol = xl; // same Sign in f and fr
356  xl = xr;
357  xr = tol;
358  tol = fl;
359  fl = fr;
360  fr = tol;
361  }
362  tol = f - fl;
363  G4double a = f*tol;
364  a = a + a;
365  if (a < fr*(fr - fl) && i <= iend) goto label17;
366  xr = *x;
367  fr = f;
368 
369 // Test on satisfactory accuracy in bisection loop
370  tol = eps;
371  a = std::abs(xr);
372  if (a > 1.) tol = tol*a;
373  if (std::abs(xr - xl) <= tol && std::abs(fr - fl) <= tolf) goto label14;
374  }
375 // End of bisection loop
376 
377 // No convergence after iend iteration steps followed by iend
378 // successive steps of bisection or steadily increasing function
379 // values at right bounds. Error return.
380  ier = 1;
381 
382 label14:
383  if (std::abs(fr) > std::abs(fl)) {
384  *x = xl;
385  f = fl;
386  }
387  return ier;
388 
389 // Computation of iterated x-value by inverse parabolic interp
390 label17:
391  G4double a = fr - f;
392  G4double dx = (*x - xl)*fl*(1. + f*(a - tol)/(a*(fr - fl)))/tol;
393  G4double xm = *x;
394  G4double fm = f;
395  *x = xl - dx;
396  tol = *x;
397  f = Fctcos(tol, aa, bb, cc, dd, rr);
398  if (f == 0.) return ier;
399 
400 // Test on satisfactory accuracy in iteration loop
401  tol = eps;
402  a = std::abs(*x);
403  if (a > 1) tol = tol*a;
404  if (std::abs(dx) <= tol && std::abs(f) <= tolf) return ier;
405 
406 // Preparation of next bisection loop
407  if (f*fl < 0.) {
408  xr = *x;
409  fr = f;
410  }
411  else {
412  xl = *x;
413  fl = f;
414  xr = xm;
415  fr = fm;
416  }
417  goto label4;
418 }
419 
420 
421 // Test function for root-finder
422 
423 G4double
424 G4LElastic::Fctcos(G4double t,
425  G4double aa, G4double bb, G4double cc, G4double dd,
426  G4double rr)
427 {
428  const G4double expxl = -82.;
429  const G4double expxu = 82.;
430 
431  G4double test1 = -bb*t;
432  if (test1 > expxu) test1 = expxu;
433  if (test1 < expxl) test1 = expxl;
434 
435  G4double test2 = -dd*t;
436  if (test2 > expxu) test2 = expxu;
437  if (test2 < expxl) test2 = expxl;
438 
439  return aa*std::exp(test1) + cc*std::exp(test2) - rr;
440 }
441 
442 
443 void
444 G4LElastic::Defs1(G4double p, G4double px, G4double py, G4double pz,
445  G4double pxinc, G4double pyinc, G4double pzinc,
446  G4double* pxnew, G4double* pynew, G4double* pznew)
447 {
448 // Transform scattered particle to reflect direction of incident particle
449  G4double pt2 = pxinc*pxinc + pyinc*pyinc;
450  if (pt2 > 0.) {
451  G4double cost = pzinc/p;
452  G4double sint1 = std::sqrt(std::abs((1. - cost )*(1.+cost)));
453  G4double sint2 = std::sqrt(pt2)/p;
454  G4double sint = 0.5*(sint1 + sint2);
455  G4double ph = pi*0.5;
456  if (pyinc < 0.) ph = pi*1.5;
457  if (std::abs(pxinc) > 1.e-6) ph = std::atan2(pyinc, pxinc);
458  G4double cosp = std::cos(ph);
459  G4double sinp = std::sin(ph);
460  if (verboseLevel > 1) {
461  G4cout << "cost sint " << cost << " " << sint << G4endl;
462  G4cout << "cosp sinp " << cosp << " " << sinp << G4endl;
463  }
464  *pxnew = cost*cosp*px - sinp*py + sint*cosp*pz;
465  *pynew = cost*sinp*px + cosp*py + sint*sinp*pz;
466  *pznew = -sint*px +cost*pz;
467  }
468  else {
469  G4double cost=pzinc/p;
470  *pxnew = cost*px;
471  *pynew = py;
472  *pznew = cost*pz;
473  }
474 }
475 
476 const std::pair<G4double, G4double> G4LElastic::GetFatalEnergyCheckLevels() const
477 {
478  return std::pair<G4double, G4double>(5*perCent,250*GeV); // max energy non-conservation is mass of heavy nucleus
479 }