Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ComponentGGHadronNucleusXsc.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 // author: V. Grichine
27 //
28 // 25.04.12 V. Grichine - first implementation
29 
31 
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4ParticleTable.hh"
35 #include "G4IonTable.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4DynamicParticle.hh"
38 #include "G4HadronNucleonXsc.hh"
39 #include "G4Log.hh"
40 #include "G4Exp.hh"
41 #include "G4Pow.hh"
42 
44 //
45 
47  : G4VComponentCrossSection(Default_Name()),
48 // fUpperLimit(100000*GeV),
49  fLowerLimit(10.*MeV),// fLowerLimit(3*GeV),
50  fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
51  fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
52  fDiffractionXsc(0.0)
53 // , fHadronNucleonXsc(0.0)
54 {
55  theGamma = G4Gamma::Gamma();
56  theProton = G4Proton::Proton();
57  theNeutron = G4Neutron::Neutron();
58  theAProton = G4AntiProton::AntiProton();
59  theANeutron = G4AntiNeutron::AntiNeutron();
60  thePiPlus = G4PionPlus::PionPlus();
61  thePiMinus = G4PionMinus::PionMinus();
62  thePiZero = G4PionZero::PionZero();
63  theKPlus = G4KaonPlus::KaonPlus();
64  theKMinus = G4KaonMinus::KaonMinus();
67  theL = G4Lambda::Lambda();
68  theAntiL = G4AntiLambda::AntiLambda();
69  theSPlus = G4SigmaPlus::SigmaPlus();
70  theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
71  theSMinus = G4SigmaMinus::SigmaMinus();
72  theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
73  theS0 = G4SigmaZero::SigmaZero();
75  theXiMinus = G4XiMinus::XiMinus();
76  theXi0 = G4XiZero::XiZero();
77  theAXiMinus = G4AntiXiMinus::AntiXiMinus();
78  theAXi0 = G4AntiXiZero::AntiXiZero();
79  theOmega = G4OmegaMinus::OmegaMinus();
81  theD = G4Deuteron::Deuteron();
82  theT = G4Triton::Triton();
83  theA = G4Alpha::Alpha();
84  theHe3 = G4He3::He3();
85 
86  hnXsc = new G4HadronNucleonXsc();
87 }
88 
90 //
91 //
92 
94 {
95  if (hnXsc) delete hnXsc;
96 }
97 
99 
101  G4double kinEnergy,
102  G4int Z, G4int A)
103 {
104  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
105  kinEnergy);
106  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
107  delete aDP;
108 
109  return fTotalXsc;
110 }
111 
113 
115  G4double kinEnergy,
116  G4int Z, G4double A)
117 {
118  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
119  kinEnergy);
120  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
121  delete aDP;
122 
123  return fTotalXsc;
124 }
125 
127 
129  G4double kinEnergy,
130  G4int Z, G4int A)
131 {
132  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
133  kinEnergy);
134  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
135  delete aDP;
136 
137  return fInelasticXsc;
138 }
139 
141 
143  G4double kinEnergy,
144  G4int Z, G4int A)
145 {
146  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
147  kinEnergy);
148  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
149  delete aDP;
150 
151  return fProductionXsc;
152 }
153 
155 
157  G4double kinEnergy,
158  G4int Z, G4double A)
159 {
160  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
161  kinEnergy);
162  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
163  delete aDP;
164 
165  return fInelasticXsc;
166 }
167 
169 
171  G4double kinEnergy,
172  G4int Z, G4double A)
173 {
174  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
175  kinEnergy);
176  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
177  delete aDP;
178 
179  return fProductionXsc;
180 }
181 
183 
185  G4double kinEnergy,
186  G4int Z, G4double A)
187 {
188  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
189  kinEnergy);
190  fTotalXsc = GetIsoCrossSection(aDP, Z, G4int(A));
191  delete aDP;
192 
193  return fElasticXsc;
194 }
195 
197 
199  G4double kinEnergy,
200  G4int Z, G4int A)
201 {
202  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
203  kinEnergy);
204  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
205  delete aDP;
206 
207  return fElasticXsc;
208 }
209 
211 
213  G4double kinEnergy,
214  G4int Z, G4int A)
215 {
216  G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
217  kinEnergy);
218  fTotalXsc = GetIsoCrossSection(aDP, Z, A);
219  delete aDP;
220  G4double ratio = 0.;
221 
222  if(fInelasticXsc > 0.)
223  {
224  ratio = (fInelasticXsc - fProductionXsc)/fInelasticXsc;
225  if(ratio < 0.) ratio = 0.;
226  }
227  return ratio;
228 }
229 
230 
231 
232 
234 
235 G4bool
237  G4int Z, G4int /*A*/,
238  const G4Element*,
239  const G4Material*)
240 {
241  G4bool applicable = false;
242  // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
243  G4double kineticEnergy = aDP->GetKineticEnergy();
244 
245  const G4ParticleDefinition* theParticle = aDP->GetDefinition();
246 
247  if (
248  Z >= 1 // >= H for kaons
249  &&
250  (
251  kineticEnergy >= fLowerLimit
252  &&
253  // Z > 1 && // >= He
254  (
255  theParticle == theAProton ||
256  theParticle == theGamma ||
257  theParticle == theSMinus ||
258  theParticle == theProton ||
259  theParticle == theNeutron ||
260  theParticle == thePiPlus ||
261  theParticle == thePiMinus
262  )
263  )
264  )
265  applicable = true;
266 
267  if (
268  Z >= 1 // >= H for kaons
269  &&
270  (
271  kineticEnergy >= 0.01*fLowerLimit
272  &&
273  (
274  theParticle == theKPlus ||
275  theParticle == theKMinus ||
276  theParticle == theK0L ||
277  theParticle == theK0S
278  )
279  )
280  )
281  applicable = true;
282 
283  return applicable;
284 }
285 
287 //
288 // Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
289 // Glauber model with Gribov correction calculated in the dipole approximation on
290 // light cone. Gaussian density of point-like nucleons helps to calculate rest integrals of the model.
291 // [1] B.Z. Kopeliovich, nucl-th/0306044 + simplification above
292 
293 G4double
295  G4int Z, G4int A,
296  const G4Isotope*,
297  const G4Element*,
298  const G4Material*)
299 {
300  G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
301  G4double hpInXsc(0.), hnInXsc(0.);
302  G4double R = GetNucleusRadius(A);
303 
304  G4int N = A - Z; // number of neutrons
305  if (N < 0) N = 0;
306 
307  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
308 
309  if( theParticle == theProton ||
310  theParticle == theNeutron ||
311  theParticle == thePiPlus ||
312  theParticle == thePiMinus )
313  {
314  // sigma = GetHadronNucleonXscNS(aParticle, A, Z);
315 
316  sigma = Z*hnXsc->GetHadronNucleonXscNS(aParticle, theProton);
317 
318  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
319 
320  sigma += N*hnXsc->GetHadronNucleonXscNS(aParticle, theNeutron);
321 
322  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
323 
324  cofInelastic = 2.4;
325  cofTotal = 2.0;
326  }
327  else if( theParticle == theKPlus ||
328  theParticle == theKMinus ||
329  theParticle == theK0S ||
330  theParticle == theK0L )
331  {
332  // sigma = GetKaonNucleonXscVector(aParticle, A, Z);
333 
334  sigma = Z*hnXsc->GetKaonNucleonXscGG(aParticle, theProton);
335 
336  hpInXsc = hnXsc->GetInelasticHadronNucleonXsc();
337 
338  sigma += N*hnXsc->GetKaonNucleonXscGG(aParticle, theNeutron);
339 
340  hnInXsc = hnXsc->GetInelasticHadronNucleonXsc();
341 
342  cofInelastic = 2.2;
343  cofTotal = 2.0;
344  R = 1.3*fermi;
345  R *= G4Pow::GetInstance()->powA(G4double(A), 0.3333);
346  }
347  else
348  {
349  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
350  cofInelastic = 2.2;
351  cofTotal = 2.0;
352  }
353  // cofInelastic = 2.0;
354 
355  if( A > 1 )
356  {
357  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
358  ratio = sigma/nucleusSquare;
359 
360  xsection = nucleusSquare*G4Log( 1. + ratio );
361 
362  xsection *= GetParticleBarCorTot(theParticle, Z);
363 
364  fTotalXsc = xsection;
365 
366  // inelastic xsc
367 
368  fAxsc2piR2 = cofInelastic*ratio;
369 
370  fModelInLog = G4Log( 1. + fAxsc2piR2 );
371 
372  fInelasticXsc = nucleusSquare*fModelInLog/cofInelastic;
373 
374  fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
375 
376  fElasticXsc = fTotalXsc - fInelasticXsc;
377 
378  if( fElasticXsc < 0. ) fElasticXsc = 0.;
379 
380  G4double difratio = ratio/(1.+ratio);
381 
382  fDiffractionXsc = 0.5*nucleusSquare*( difratio - G4Log( 1. + difratio ) );
383 
384 
385  // sigma = GetHNinelasticXsc(aParticle, A, Z);
386 
387  sigma = Z*hpInXsc + N*hnInXsc;
388 
389  ratio = sigma/nucleusSquare;
390 
391  fProductionXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
392 
393  fProductionXsc *= GetParticleBarCorIn(theParticle, Z);
394 
395  if (fElasticXsc < 0.) fElasticXsc = 0.;
396  }
397  else // H
398  {
399  fTotalXsc = sigma;
400  xsection = sigma;
401 
402  fInelasticXsc = hnXsc->GetInelasticHadronNucleonXsc();
403 
404  if ( theParticle != theAProton )
405  {
406  fElasticXsc = hnXsc->GetElasticHadronNucleonXsc();
407 
408  // sigma = GetHNinelasticXsc(aParticle, A, Z);
409  // fInelasticXsc = sigma;
410  // fElasticXsc = fTotalXsc - fInelasticXsc;
411  }
412  else if( theParticle == theKPlus ||
413  theParticle == theKMinus ||
414  theParticle == theK0S ||
415  theParticle == theK0L )
416  {
417  fInelasticXsc = hpInXsc;
418  fElasticXsc = fTotalXsc - fInelasticXsc;
419  }
420  else
421  {
422  fInelasticXsc = hpInXsc;
423  fElasticXsc = fTotalXsc - fInelasticXsc;
424  }
425  if (fElasticXsc < 0.) fElasticXsc = 0.;
426 
427  }
428  return xsection;
429 }
430 
432 //
433 // Return single-diffraction/inelastic cross-section ratio
434 
437 {
438  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
439  G4double R = GetNucleusRadius(A);
440 
441  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
442 
443  if( theParticle == theProton ||
444  theParticle == theNeutron ||
445  theParticle == thePiPlus ||
446  theParticle == thePiMinus )
447  {
448  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
449  cofInelastic = 2.4;
450  cofTotal = 2.0;
451  }
452  else
453  {
454  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
455  cofInelastic = 2.2;
456  cofTotal = 2.0;
457  }
458  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
459  ratio = sigma/nucleusSquare;
460 
461  fInelasticXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
462 
463  G4double difratio = ratio/(1.+ratio);
464 
465  fDiffractionXsc = 0.5*nucleusSquare*( difratio - G4Log( 1. + difratio ) );
466 
467  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
468  else ratio = 0.;
469 
470  return ratio;
471 }
472 
474 //
475 // Return suasi-elastic/inelastic cross-section ratio
476 
479 {
480  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
481  G4double R = GetNucleusRadius(A);
482 
483  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
484 
485  if( theParticle == theProton ||
486  theParticle == theNeutron ||
487  theParticle == thePiPlus ||
488  theParticle == thePiMinus )
489  {
490  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
491  cofInelastic = 2.4;
492  cofTotal = 2.0;
493  }
494  else
495  {
496  sigma = GetHadronNucleonXscNS(aParticle, A, Z);
497  cofInelastic = 2.2;
498  cofTotal = 2.0;
499  }
500  nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
501  ratio = sigma/nucleusSquare;
502 
503  fInelasticXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
504 
505  sigma = GetHNinelasticXsc(aParticle, A, Z);
506  ratio = sigma/nucleusSquare;
507 
508  fProductionXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )/cofInelastic;
509 
510  if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
511  else ratio = 0.;
512  if ( ratio < 0. ) ratio = 0.;
513 
514  return ratio;
515 }
516 
518 //
519 // Returns hadron-nucleon Xsc according to differnt parametrisations:
520 // [2] E. Levin, hep-ph/9710546
521 // [3] U. Dersch, et al, hep-ex/9910052
522 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
523 
524 G4double
526  const G4Element* anElement)
527 {
528  G4int At = G4lrint(anElement->GetN()); // number of nucleons
529  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
530 
531  return GetHadronNucleonXsc(aParticle, At, Zt);
532 }
533 
535 //
536 // Returns hadron-nucleon Xsc according to differnt parametrisations:
537 // [2] E. Levin, hep-ph/9710546
538 // [3] U. Dersch, et al, hep-ex/9910052
539 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
540 
541 G4double
543  G4int At, G4int /*Zt*/)
544 {
545  G4double xsection;
546 
547  //G4double targ_mass = G4NucleiProperties::GetNuclearMass(At, Zt);
548 
549  G4double targ_mass = 0.939*GeV; // ~mean neutron and proton ???
550 
551  G4double proj_mass = aParticle->GetMass();
552  G4double proj_momentum = aParticle->GetMomentum().mag();
553  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
554 
555  sMand /= GeV*GeV; // in GeV for parametrisation
556  proj_momentum /= GeV;
557 
558  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
559 
560  G4double aa = At;
561 
562  if(theParticle == theGamma)
563  {
564  xsection = aa*(0.0677*G4Pow::GetInstance()->powA(sMand,0.0808) + 0.129*G4Pow::GetInstance()->powA(sMand,-0.4525));
565  }
566  else if(theParticle == theNeutron) // as proton ???
567  {
568  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
569  }
570  else if(theParticle == theProton)
571  {
572  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
573  // xsection = At*( 49.51*G4Pow::GetInstance()->powA(sMand,-0.097) + 0.314*G4Log(sMand)*G4Log(sMand) );
574  // xsection = At*( 38.4 + 0.85*std::abs(G4Pow::GetInstance()->powA(log(sMand),1.47)) );
575  }
576  else if(theParticle == theAProton)
577  {
578  xsection = aa*( 21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 98.39*G4Pow::GetInstance()->powA(sMand,-0.4525));
579  }
580  else if(theParticle == thePiPlus)
581  {
582  xsection = aa*(13.63*G4Pow::GetInstance()->powA(sMand,0.0808) + 27.56*G4Pow::GetInstance()->powA(sMand,-0.4525));
583  }
584  else if(theParticle == thePiMinus)
585  {
586  // xsection = At*( 55.2*G4Pow::GetInstance()->powA(sMand,-0.255) + 0.346*G4Log(sMand)*G4Log(sMand) );
587  xsection = aa*(13.63*G4Pow::GetInstance()->powA(sMand,0.0808) + 36.02*G4Pow::GetInstance()->powA(sMand,-0.4525));
588  }
589  else if(theParticle == theKPlus)
590  {
591  xsection = aa*(11.82*G4Pow::GetInstance()->powA(sMand,0.0808) + 8.15*G4Pow::GetInstance()->powA(sMand,-0.4525));
592  }
593  else if(theParticle == theKMinus)
594  {
595  xsection = aa*(11.82*G4Pow::GetInstance()->powA(sMand,0.0808) + 26.36*G4Pow::GetInstance()->powA(sMand,-0.4525));
596  }
597  else // as proton ???
598  {
599  xsection = aa*(21.70*G4Pow::GetInstance()->powA(sMand,0.0808) + 56.08*G4Pow::GetInstance()->powA(sMand,-0.4525));
600  }
601  xsection *= millibarn;
602  return xsection;
603 }
604 
605 
607 //
608 // Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
609 // http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
610 
611 G4double
613  const G4Element* anElement)
614 {
615  G4int At = G4lrint(anElement->GetN()); // number of nucleons
616  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
617 
618  return GetHadronNucleonXscPDG(aParticle, At, Zt);
619 }
620 
621 
622 
623 
625 //
626 // Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
627 // http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
628 // At = number of nucleons, Zt = number of protons
629 
630 G4double
632  G4int At, G4int Zt)
633 {
634  G4double xsection;
635 
636  G4int Nt = At-Zt; // number of neutrons
637  if (Nt < 0) Nt = 0;
638 
639  G4double zz = Zt;
640  G4double aa = At;
641  G4double nn = Nt;
642 
644  GetIonTable()->GetIonMass(Zt, At);
645 
646  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
647 
648  G4double proj_mass = aParticle->GetMass();
649  G4double proj_momentum = aParticle->GetMomentum().mag();
650 
651  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
652 
653  sMand /= GeV*GeV; // in GeV for parametrisation
654 
655  // General PDG fit constants
656 
657  G4double s0 = 5.38*5.38; // in Gev^2
658  G4double eta1 = 0.458;
659  G4double eta2 = 0.458;
660  G4double B = 0.308;
661 
662 
663  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
664 
665 
666  if(theParticle == theNeutron) // proton-neutron fit
667  {
668  xsection = zz*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
669  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
670  xsection += nn*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
671  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2)); // pp for nn
672  }
673  else if(theParticle == theProton)
674  {
675 
676  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
677  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
678 
679  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
680  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
681  }
682  else if(theParticle == theAProton)
683  {
684  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
685  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) + 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
686 
687  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
688  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) + 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
689  }
690  else if(theParticle == thePiPlus)
691  {
692  xsection = aa*( 20.86 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
693  + 19.24*G4Pow::GetInstance()->powA(sMand,-eta1) - 6.03*G4Pow::GetInstance()->powA(sMand,-eta2));
694  }
695  else if(theParticle == thePiMinus)
696  {
697  xsection = aa*( 20.86 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
698  + 19.24*G4Pow::GetInstance()->powA(sMand,-eta1) + 6.03*G4Pow::GetInstance()->powA(sMand,-eta2));
699  }
700  else if(theParticle == theKPlus || theParticle == theK0L )
701  {
702  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
703  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) - 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
704 
705  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
706  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) - 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
707  }
708  else if(theParticle == theKMinus || theParticle == theK0S )
709  {
710  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
711  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) + 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
712 
713  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
714  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) + 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
715  }
716  else if(theParticle == theSMinus)
717  {
718  xsection = aa*( 35.20 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
719  - 199.*G4Pow::GetInstance()->powA(sMand,-eta1) + 264.*G4Pow::GetInstance()->powA(sMand,-eta2));
720  }
721  else if(theParticle == theGamma) // modify later on
722  {
723  xsection = aa*( 0.0 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
724  + 0.032*G4Pow::GetInstance()->powA(sMand,-eta1) - 0.0*G4Pow::GetInstance()->powA(sMand,-eta2));
725 
726  }
727  else // as proton ???
728  {
729  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
730  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
731 
732  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
733  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
734  }
735  xsection *= millibarn; // parametrised in mb
736  return xsection;
737 }
738 
739 
741 //
742 // Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
743 // data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
744 
745 G4double
747  const G4Element* anElement)
748 {
749  G4int At = G4lrint(anElement->GetN()); // number of nucleons
750  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
751 
752  return GetHadronNucleonXscNS(aParticle, At, Zt);
753 }
754 
755 
756 
757 
759 //
760 // Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
761 // data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
762 
763 G4double
765  G4int At, G4int Zt)
766 {
767  G4double xsection(0);
768  // G4double Delta; DHW 19 May 2011: variable set but not used
769  G4double A0, B0;
770  G4double hpXscv(0);
771  G4double hnXscv(0);
772 
773  G4int Nt = At-Zt; // number of neutrons
774  if (Nt < 0) Nt = 0;
775 
776  G4double aa = At;
777  G4double zz = Zt;
778  G4double nn = Nt;
779 
781  GetIonTable()->GetIonMass(Zt, At);
782 
783  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
784 
785  G4double proj_mass = aParticle->GetMass();
786  G4double proj_energy = aParticle->GetTotalEnergy();
787  G4double proj_momentum = aParticle->GetMomentum().mag();
788 
789  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
790 
791  sMand /= GeV*GeV; // in GeV for parametrisation
792  proj_momentum /= GeV;
793  proj_energy /= GeV;
794  proj_mass /= GeV;
795 
796  // General PDG fit constants
797 
798  G4double s0 = 5.38*5.38; // in Gev^2
799  G4double eta1 = 0.458;
800  G4double eta2 = 0.458;
801  G4double B = 0.308;
802 
803 
804  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
805 
806 
807  if(theParticle == theNeutron)
808  {
809  if( proj_momentum >= 373.)
810  {
811  return GetHadronNucleonXscPDG(aParticle,At,Zt);
812  }
813  else if( proj_momentum >= 10.)
814  // if( proj_momentum >= 2.)
815  {
816  // Delta = 1.; // DHW 19 May 2011: variable set but not used
817  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
818 
819  //AR-12Aug2016 if(proj_momentum >= 10.)
820  {
821  B0 = 7.5;
822  A0 = 100. - B0*G4Log(3.0e7);
823 
824  xsection = A0 + B0*G4Log(proj_energy) - 11
825  + 103*G4Pow::GetInstance()->powA(2*0.93827*proj_energy + proj_mass*proj_mass+
826  0.93827*0.93827,-0.165); // mb
827  }
828  xsection *= zz + nn;
829  }
830  else
831  {
832  // nn to be pp
833 
834  if( proj_momentum < 0.73 )
835  {
836  hnXscv = 23 + 50*( G4Pow::GetInstance()->powA( G4Log(0.73/proj_momentum), 3.5 ) );
837  }
838  else if( proj_momentum < 1.05 )
839  {
840  hnXscv = 23 + 40*(G4Log(proj_momentum/0.73))*
841  (G4Log(proj_momentum/0.73));
842  }
843  else // if( proj_momentum < 10. )
844  {
845  hnXscv = 39.0+
846  75*(proj_momentum - 1.2)/(G4Pow::GetInstance()->powA(proj_momentum,3.0) + 0.15);
847  }
848  // pn to be np
849 
850  if( proj_momentum < 0.8 )
851  {
852  hpXscv = 33+30*G4Pow::GetInstance()->powA(G4Log(proj_momentum/1.3),4.0);
853  }
854  else if( proj_momentum < 1.4 )
855  {
856  hpXscv = 33+30*G4Pow::GetInstance()->powA(G4Log(proj_momentum/0.95),2.0);
857  }
858  else // if( proj_momentum < 10. )
859  {
860  hpXscv = 33.3+
861  20.8*(G4Pow::GetInstance()->powA(proj_momentum,2.0)-1.35)/
862  (G4Pow::GetInstance()->powA(proj_momentum,2.50)+0.95);
863  }
864  xsection = hpXscv*zz + hnXscv*nn;
865  }
866  }
867  else if(theParticle == theProton)
868  {
869  if( proj_momentum >= 373.)
870  {
871  return GetHadronNucleonXscPDG(aParticle,At,Zt);
872  }
873  else if( proj_momentum >= 10.)
874  // if( proj_momentum >= 2.)
875  {
876  // Delta = 1.; DHW 19 May 2011: variable set but not used
877  // if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
878 
879  //AR-12Aug2016 if(proj_momentum >= 10.)
880  {
881  B0 = 7.5;
882  A0 = 100. - B0*G4Log(3.0e7);
883 
884  xsection = A0 + B0*G4Log(proj_energy) - 11
885  + 103*G4Pow::GetInstance()->powA(2*0.93827*proj_energy + proj_mass*proj_mass+
886  0.93827*0.93827,-0.165); // mb
887  }
888  xsection *= zz + nn;
889  }
890  else
891  {
892  // pp
893 
894  if( proj_momentum < 0.73 )
895  {
896  hpXscv = 23 + 50*( G4Pow::GetInstance()->powA( G4Log(0.73/proj_momentum), 3.5 ) );
897  }
898  else if( proj_momentum < 1.05 )
899  {
900  hpXscv = 23 + 40*(G4Log(proj_momentum/0.73))*
901  (G4Log(proj_momentum/0.73));
902  }
903  else // if( proj_momentum < 10. )
904  {
905  hpXscv = 39.0+
906  75*(proj_momentum - 1.2)/(G4Pow::GetInstance()->powA(proj_momentum,3.0) + 0.15);
907  }
908  // pn to be np
909 
910  if( proj_momentum < 0.8 )
911  {
912  hnXscv = 33+30*G4Pow::GetInstance()->powA(G4Log(proj_momentum/1.3),4.0);
913  }
914  else if( proj_momentum < 1.4 )
915  {
916  hnXscv = 33+30*G4Pow::GetInstance()->powA(G4Log(proj_momentum/0.95),2.0);
917  }
918  else // if( proj_momentum < 10. )
919  {
920  hnXscv = 33.3+
921  20.8*(G4Pow::GetInstance()->powA(proj_momentum,2.0)-1.35)/
922  (G4Pow::GetInstance()->powA(proj_momentum,2.50)+0.95);
923  }
924  xsection = hpXscv*zz + hnXscv*nn;
925  // xsection = hpXscv*(Zt + Nt);
926  // xsection = hnXscv*(Zt + Nt);
927  }
928  // xsection *= 0.95;
929  }
930  else if( theParticle == theAProton )
931  {
932  // xsection = Zt*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
933  // + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) + 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
934 
935  // xsection += Nt*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
936  // + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) + 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
937 
938  G4double logP = G4Log(proj_momentum);
939 
940  if( proj_momentum <= 1.0 )
941  {
942  xsection = zz*(65.55 + 53.84/(proj_momentum+1.e-6) );
943  }
944  else
945  {
946  xsection = zz*( 41.1 + 77.2*G4Pow::GetInstance()->powA( proj_momentum, -0.68)
947  + 0.293*logP*logP - 1.82*logP );
948  }
949  if ( nn > 0.)
950  {
951  xsection += nn*( 41.9 + 96.2*G4Pow::GetInstance()->powA( proj_momentum, -0.99) - 0.154*logP);
952  }
953  else // H
954  {
955  fInelasticXsc = 38.0 + 38.0*G4Pow::GetInstance()->powA( proj_momentum, -0.96)
956  - 0.169*logP*logP;
957  fInelasticXsc *= millibarn;
958  }
959  }
960  else if( theParticle == thePiPlus )
961  {
962  if(proj_momentum < 0.4)
963  {
964  G4double Ex3 = 180*G4Exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
965  hpXscv = Ex3+20.0;
966  }
967  else if( proj_momentum < 1.15 )
968  {
969  G4double Ex4 = 88*(G4Log(proj_momentum/0.75))*(G4Log(proj_momentum/0.75));
970  hpXscv = Ex4+14.0;
971  }
972  else if(proj_momentum < 3.5)
973  {
974  G4double Ex1 = 3.2*G4Exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
975  G4double Ex2 = 12*G4Exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
976  hpXscv = Ex1+Ex2+27.5;
977  }
978  else // if(proj_momentum > 3.5) // mb
979  {
980  hpXscv = 10.6+2.*G4Log(proj_energy)+25*G4Pow::GetInstance()->powA(proj_energy,-0.43);
981  }
982  // pi+n = pi-p??
983 
984  if(proj_momentum < 0.37)
985  {
986  hnXscv = 28.0 + 40*G4Exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
987  }
988  else if(proj_momentum<0.65)
989  {
990  hnXscv = 26+110*(G4Log(proj_momentum/0.48))*(G4Log(proj_momentum/0.48));
991  }
992  else if(proj_momentum<1.3)
993  {
994  hnXscv = 36.1+
995  10*G4Exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
996  24*G4Exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
997  }
998  else if(proj_momentum<3.0)
999  {
1000  hnXscv = 36.1+0.079-4.313*G4Log(proj_momentum)+
1001  3*G4Exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1002  1.5*G4Exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1003  }
1004  else // mb
1005  {
1006  hnXscv = 10.6+2*G4Log(proj_energy)+30*G4Pow::GetInstance()->powA(proj_energy,-0.43);
1007  }
1008  xsection = hpXscv*zz + hnXscv*nn;
1009  }
1010  else if(theParticle == thePiMinus)
1011  {
1012  // pi-n = pi+p??
1013 
1014  if(proj_momentum < 0.4)
1015  {
1016  G4double Ex3 = 180*G4Exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
1017  hnXscv = Ex3+20.0;
1018  }
1019  else if(proj_momentum < 1.15)
1020  {
1021  G4double Ex4 = 88*(G4Log(proj_momentum/0.75))*(G4Log(proj_momentum/0.75));
1022  hnXscv = Ex4+14.0;
1023  }
1024  else if(proj_momentum < 3.5)
1025  {
1026  G4double Ex1 = 3.2*G4Exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
1027  G4double Ex2 = 12*G4Exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
1028  hnXscv = Ex1+Ex2+27.5;
1029  }
1030  else // if(proj_momentum > 3.5) // mb
1031  {
1032  hnXscv = 10.6+2.*G4Log(proj_energy)+25*G4Pow::GetInstance()->powA(proj_energy,-0.43);
1033  }
1034  // pi-p
1035 
1036  if(proj_momentum < 0.37)
1037  {
1038  hpXscv = 28.0 + 40*G4Exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
1039  }
1040  else if(proj_momentum<0.65)
1041  {
1042  hpXscv = 26+110*(G4Log(proj_momentum/0.48))*(G4Log(proj_momentum/0.48));
1043  }
1044  else if(proj_momentum<1.3)
1045  {
1046  hpXscv = 36.1+
1047  10*G4Exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
1048  24*G4Exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
1049  }
1050  else if(proj_momentum<3.0)
1051  {
1052  hpXscv = 36.1+0.079-4.313*G4Log(proj_momentum)+
1053  3*G4Exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1054  1.5*G4Exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1055  }
1056  else // mb
1057  {
1058  hpXscv = 10.6+2*G4Log(proj_energy)+30*G4Pow::GetInstance()->powA(proj_energy,-0.43);
1059  }
1060  xsection = hpXscv*zz + hnXscv*nn;
1061  }
1062  else if(theParticle == theKPlus)
1063  {
1064  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1065  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) - 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
1066 
1067  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1068  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) - 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
1069  }
1070  else if(theParticle == theKMinus)
1071  {
1072  xsection = zz*( 17.91 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1073  + 7.14*G4Pow::GetInstance()->powA(sMand,-eta1) + 13.45*G4Pow::GetInstance()->powA(sMand,-eta2));
1074 
1075  xsection += nn*( 17.87 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1076  + 5.17*G4Pow::GetInstance()->powA(sMand,-eta1) + 7.23*G4Pow::GetInstance()->powA(sMand,-eta2));
1077  }
1078  else if(theParticle == theSMinus)
1079  {
1080  xsection = aa*( 35.20 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1081  - 199.*G4Pow::GetInstance()->powA(sMand,-eta1) + 264.*G4Pow::GetInstance()->powA(sMand,-eta2));
1082  }
1083  else if(theParticle == theGamma) // modify later on
1084  {
1085  xsection = aa*( 0.0 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1086  + 0.032*G4Pow::GetInstance()->powA(sMand,-eta1) - 0.0*G4Pow::GetInstance()->powA(sMand,-eta2));
1087 
1088  }
1089  else // as proton ???
1090  {
1091  xsection = zz*( 35.45 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1092  + 42.53*G4Pow::GetInstance()->powA(sMand,-eta1) - 33.34*G4Pow::GetInstance()->powA(sMand,-eta2));
1093 
1094  xsection += nn*( 35.80 + B*G4Pow::GetInstance()->powA(G4Log(sMand/s0),2.)
1095  + 40.15*G4Pow::GetInstance()->powA(sMand,-eta1) - 30.*G4Pow::GetInstance()->powA(sMand,-eta2));
1096  }
1097  xsection *= millibarn; // parametrised in mb
1098  return xsection;
1099 }
1100 
1101 /*
1102 G4double
1103 G4ComponentGGHadronNucleusXsc::GetKaonNucleonXscVector(const G4DynamicParticle* aParticle,
1104  G4int At, G4int Zt)
1105 {
1106  G4double Tkin, logTkin, xsc, xscP, xscN;
1107  const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
1108 
1109  G4int Nt = At-Zt; // number of neutrons
1110  if (Nt < 0) Nt = 0;
1111 
1112  Tkin = aParticle->GetKineticEnergy(); // Tkin in MeV
1113 
1114  if( Tkin > 70*GeV ) return GetHadronNucleonXscPDG(aParticle,At,Zt);
1115 
1116  logTkin = G4Log(Tkin); // Tkin in MeV!!!
1117 
1118  if( theParticle == theKPlus )
1119  {
1120  xscP = hnXsc->GetKpProtonTotXscVector(logTkin);
1121  xscN = hnXsc->GetKpNeutronTotXscVector(logTkin);
1122  }
1123  else if( theParticle == theKMinus )
1124  {
1125  xscP = hnXsc->GetKmProtonTotXscVector(logTkin);
1126  xscN = hnXsc->GetKmNeutronTotXscVector(logTkin);
1127  }
1128  else // K-zero as half of K+ and K-
1129  {
1130  xscP = (hnXsc->GetKpProtonTotXscVector(logTkin)+hnXsc->GetKmProtonTotXscVector(logTkin))*0.5;
1131  xscN = (hnXsc->GetKpNeutronTotXscVector(logTkin)+hnXsc->GetKmNeutronTotXscVector(logTkin))*0.5;
1132  }
1133  xsc = xscP*Zt + xscN*Nt;
1134  return xsc;
1135 }
1136 */
1137 
1139 //
1140 // Returns hadron-nucleon inelastic cross-section based on proper parametrisation
1141 
1142 G4double
1144  const G4Element* anElement)
1145 {
1146  G4int At = G4lrint(anElement->GetN()); // number of nucleons
1147  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
1148 
1149  return GetHNinelasticXsc(aParticle, At, Zt);
1150 }
1151 
1153 //
1154 // Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1155 
1156 G4double
1158  G4int At, G4int Zt)
1159 {
1160  const G4ParticleDefinition* hadron = aParticle->GetDefinition();
1161  G4double sumInelastic;
1162  G4int Nt = At - Zt;
1163  if(Nt < 0) Nt = 0;
1164 
1165  if( hadron == theKPlus )
1166  {
1167  sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
1168  }
1169  else
1170  {
1171  //sumInelastic = Zt*GetHadronNucleonXscMK(aParticle, theProton);
1172  // sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);
1173  sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
1174  sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
1175  }
1176  return sumInelastic;
1177 }
1178 
1179 
1181 //
1182 // Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1183 
1184 G4double
1186  G4int At, G4int Zt)
1187 {
1188  G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
1189  G4int absPDGcode = std::abs(PDGcode);
1190 
1191  G4double Elab = aParticle->GetTotalEnergy();
1192  // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
1193  G4double Plab = aParticle->GetMomentum().mag();
1194  // std::sqrt(Elab * Elab - 0.88);
1195 
1196  Elab /= GeV;
1197  Plab /= GeV;
1198 
1199  G4double LogPlab = G4Log( Plab );
1200  G4double sqrLogPlab = LogPlab * LogPlab;
1201 
1202  //G4cout<<"Plab = "<<Plab<<G4endl;
1203 
1204  G4double NumberOfTargetProtons = G4double(Zt);
1205  G4double NumberOfTargetNucleons = G4double(At);
1206  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
1207 
1208  if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
1209 
1210  G4double Xtotal, Xelastic, Xinelastic;
1211 
1212  if( absPDGcode > 1000 ) //------Projectile is baryon --------
1213  {
1214  G4double XtotPP = 48.0 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1215  0.522*sqrLogPlab - 4.51*LogPlab;
1216 
1217  G4double XtotPN = 47.3 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1218  0.513*sqrLogPlab - 4.27*LogPlab;
1219 
1220  G4double XelPP = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1221  0.169*sqrLogPlab - 1.85*LogPlab;
1222 
1223  G4double XelPN = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1224  0.169*sqrLogPlab - 1.85*LogPlab;
1225 
1226  Xtotal = (NumberOfTargetProtons * XtotPP +
1227  NumberOfTargetNeutrons * XtotPN);
1228 
1229  Xelastic = (NumberOfTargetProtons * XelPP +
1230  NumberOfTargetNeutrons * XelPN);
1231  }
1232  else if( PDGcode == 211 ) //------Projectile is PionPlus -------
1233  {
1234  G4double XtotPiP = 16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1235  0.19 *sqrLogPlab - 0.0 *LogPlab;
1236 
1237  G4double XtotPiN = 33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1238  0.456*sqrLogPlab - 4.03*LogPlab;
1239 
1240  G4double XelPiP = 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1241  0.079*sqrLogPlab - 0.0 *LogPlab;
1242 
1243  G4double XelPiN = 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1244  0.043*sqrLogPlab - 0.0 *LogPlab;
1245 
1246  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1247  NumberOfTargetNeutrons * XtotPiN );
1248 
1249  Xelastic = ( NumberOfTargetProtons * XelPiP +
1250  NumberOfTargetNeutrons * XelPiN );
1251  }
1252  else if( PDGcode == -211 ) //------Projectile is PionMinus -------
1253  {
1254  G4double XtotPiP = 33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1255  0.456*sqrLogPlab - 4.03*LogPlab;
1256 
1257  G4double XtotPiN = 16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1258  0.19 *sqrLogPlab - 0.0 *LogPlab;
1259 
1260  G4double XelPiP = 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1261  0.043*sqrLogPlab - 0.0 *LogPlab;
1262 
1263  G4double XelPiN = 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1264  0.079*sqrLogPlab - 0.0 *LogPlab;
1265 
1266  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1267  NumberOfTargetNeutrons * XtotPiN );
1268 
1269  Xelastic = ( NumberOfTargetProtons * XelPiP +
1270  NumberOfTargetNeutrons * XelPiN );
1271  }
1272  else if( PDGcode == 111 ) //------Projectile is PionZero -------
1273  {
1274  G4double XtotPiP =(16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1275  0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1276  33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1277  0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1278 
1279  G4double XtotPiN =(33.0 + 14.0 *G4Pow::GetInstance()->powA(Plab,-1.36) +
1280  0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1281  16.4 + 19.3 *G4Pow::GetInstance()->powA(Plab,-0.42) +
1282  0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1283 
1284  G4double XelPiP =( 0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1285  0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1286  1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1287  0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1288 
1289  G4double XelPiN =( 1.76 + 11.2*G4Pow::GetInstance()->powA(Plab,-0.64) +
1290  0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1291  0.0 + 11.4*G4Pow::GetInstance()->powA(Plab,-0.40) +
1292  0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1293 
1294  Xtotal = ( NumberOfTargetProtons * XtotPiP +
1295  NumberOfTargetNeutrons * XtotPiN );
1296 
1297  Xelastic = ( NumberOfTargetProtons * XelPiP +
1298  NumberOfTargetNeutrons * XelPiN );
1299  }
1300  else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
1301  {
1302  G4double XtotKP = 18.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1303  0.26 *sqrLogPlab - 1.0 *LogPlab;
1304  G4double XtotKN = 18.7 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1305  0.21 *sqrLogPlab - 0.89*LogPlab;
1306 
1307  G4double XelKP = 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1308  0.16 *sqrLogPlab - 1.3 *LogPlab;
1309 
1310  G4double XelKN = 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1311  0.29 *sqrLogPlab - 2.4 *LogPlab;
1312 
1313  Xtotal = ( NumberOfTargetProtons * XtotKP +
1314  NumberOfTargetNeutrons * XtotKN );
1315 
1316  Xelastic = ( NumberOfTargetProtons * XelKP +
1317  NumberOfTargetNeutrons * XelKN );
1318  }
1319  else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
1320  {
1321  G4double XtotKP = 32.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1322  0.66 *sqrLogPlab - 5.6 *LogPlab;
1323  G4double XtotKN = 25.2 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1324  0.38 *sqrLogPlab - 2.9 *LogPlab;
1325 
1326  G4double XelKP = 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1327  0.29 *sqrLogPlab - 2.4 *LogPlab;
1328 
1329  G4double XelKN = 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1330  0.16 *sqrLogPlab - 1.3 *LogPlab;
1331 
1332  Xtotal = ( NumberOfTargetProtons * XtotKP +
1333  NumberOfTargetNeutrons * XtotKN );
1334 
1335  Xelastic = ( NumberOfTargetProtons * XelKP +
1336  NumberOfTargetNeutrons * XelKN );
1337  }
1338  else if( PDGcode == 311 ) //------Projectile is KaonZero ------
1339  {
1340  G4double XtotKP = ( 18.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1341  0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1342  32.1 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1343  0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1344 
1345  G4double XtotKN = ( 18.7 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1346  0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1347  25.2 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1348  0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1349 
1350  G4double XelKP = ( 5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 )
1351  + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1352  7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1353  0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1354 
1355  G4double XelKN = ( 7.3 + 0. *G4Pow::GetInstance()->powA(Plab,-0. ) +
1356  0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1357  5.0 + 8.1*G4Pow::GetInstance()->powA(Plab,-1.8 ) +
1358  0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1359 
1360  Xtotal = ( NumberOfTargetProtons * XtotKP +
1361  NumberOfTargetNeutrons * XtotKN );
1362 
1363  Xelastic = ( NumberOfTargetProtons * XelKP +
1364  NumberOfTargetNeutrons * XelKN );
1365  }
1366  else //------Projectile is undefined, Nucleon assumed
1367  {
1368  G4double XtotPP = 48.0 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1369  0.522*sqrLogPlab - 4.51*LogPlab;
1370 
1371  G4double XtotPN = 47.3 + 0. *G4Pow::GetInstance()->powA(Plab, 0. ) +
1372  0.513*sqrLogPlab - 4.27*LogPlab;
1373 
1374  G4double XelPP = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1375  0.169*sqrLogPlab - 1.85*LogPlab;
1376  G4double XelPN = 11.9 + 26.9*G4Pow::GetInstance()->powA(Plab,-1.21) +
1377  0.169*sqrLogPlab - 1.85*LogPlab;
1378 
1379  Xtotal = ( NumberOfTargetProtons * XtotPP +
1380  NumberOfTargetNeutrons * XtotPN );
1381 
1382  Xelastic = ( NumberOfTargetProtons * XelPP +
1383  NumberOfTargetNeutrons * XelPN );
1384  }
1385  Xinelastic = Xtotal - Xelastic;
1386 
1387  if( Xinelastic < 0.) Xinelastic = 0.;
1388 
1389  return Xinelastic*= millibarn;
1390 }
1391 
1393 //
1394 //
1395 
1396 G4double
1398  const G4Element* anElement)
1399 {
1400  G4int At = G4lrint(anElement->GetN());
1401  G4double oneThird = 1.0/3.0;
1402  G4double cubicrAt = G4Pow::GetInstance()->powA(G4double(At), oneThird);
1403 
1404  G4double R; // = fRadiusConst*cubicrAt;
1405  /*
1406  G4double tmp = G4Pow::GetInstance()->powA( cubicrAt-1., 3.);
1407  tmp += At;
1408  tmp *= 0.5;
1409 
1410  if (At > 20.) // 20.
1411  {
1412  R = fRadiusConst*G4Pow::GetInstance()->powA (tmp, oneThird);
1413  }
1414  else
1415  {
1416  R = fRadiusConst*cubicrAt;
1417  }
1418  */
1419 
1420  R = fRadiusConst*cubicrAt;
1421 
1422  G4double meanA = 21.;
1423 
1424  G4double tauA1 = 40.;
1425  G4double tauA2 = 10.;
1426  G4double tauA3 = 5.;
1427 
1428  G4double a1 = 0.85;
1429  G4double b1 = 1. - a1;
1430 
1431  G4double b2 = 0.3;
1432  G4double b3 = 4.;
1433 
1434  if (At > 20) // 20.
1435  {
1436  R *= ( a1 + b1*G4Exp( -(At - meanA)/tauA1) );
1437  }
1438  else if (At > 3)
1439  {
1440  R *= ( 1.0 + b2*( 1. - G4Exp( (At - meanA)/tauA2) ) );
1441  }
1442  else
1443  {
1444  R *= ( 1.0 + b3*( 1. - G4Exp( (At - meanA)/tauA3) ) );
1445  }
1446  return R;
1447 
1448 }
1450 //
1451 //
1452 
1453 G4double
1455 {
1456  G4double oneThird = 1.0/3.0;
1457  G4double cubicrAt = G4Pow::GetInstance()->powA(G4double(At), oneThird);
1458 
1459  G4double R; // = fRadiusConst*cubicrAt;
1460 
1461  /*
1462  G4double tmp = G4Pow::GetInstance()->powA( cubicrAt-1., 3.);
1463  tmp += At;
1464  tmp *= 0.5;
1465 
1466  if (At > 20.)
1467  {
1468  R = fRadiusConst*G4Pow::GetInstance()->powA (tmp, oneThird);
1469  }
1470  else
1471  {
1472  R = fRadiusConst*cubicrAt;
1473  }
1474  */
1475 
1476  R = fRadiusConst*cubicrAt;
1477 
1478  G4double meanA = 20.;
1479  G4double tauA = 20.;
1480 
1481  if (At > 20) // 20.
1482  {
1483  R *= ( 0.8 + 0.2*G4Exp( -(G4double(At) - meanA)/tauA) );
1484  }
1485  else
1486  {
1487  R *= ( 1.0 + 0.1*( 1. - G4Exp( (G4double(At) - meanA)/tauA) ) );
1488  }
1489 
1490  return R;
1491 }
1492 
1494 //
1495 //
1496 
1498  const G4double mt ,
1499  const G4double Plab )
1500 {
1501  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1502  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1503  // G4double Pcm = Plab * mt / Ecm;
1504  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1505 
1506  return Ecm ; // KEcm;
1507 }
1508 
1510 //
1511 //
1512 
1514  const G4double mt ,
1515  const G4double Plab )
1516 {
1517  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1518  G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1519 
1520  return sMand;
1521 }
1522 
1524 //
1525 //
1526 
1528 {
1529  outFile << "G4ComponentGGHadronNucleusXsc calculates total, inelastic and\n"
1530  << "elastic cross sections for hadron-nucleus cross sections using\n"
1531  << "the Glauber model with Gribov corrections. It is valid for all\n"
1532  << "targets except hydrogen, and for incident p, pbar, n, sigma-,\n"
1533  << "pi+, pi-, K+, K- and gammas with energies above 3 GeV. This is\n"
1534  << "a cross section component which is to be used to build a cross\n"
1535  << "data set.\n";
1536 }
1537 
1538 
1540 //
1541 // Correction arrays for GG <-> Bar changea at ~ 90 GeV
1542 
1543 const G4double G4ComponentGGHadronNucleusXsc::fNeutronBarCorrectionTot[93] = {
1544 
1545  1.0, 1.0, 1.42517e+00, // 1.118517e+00,
1546 1.082002e+00, 1.116171e+00, 1.078747e+00, 1.061315e+00,
1547 1.058205e+00, 1.082663e+00, 1.068500e+00, 1.076912e+00, 1.083475e+00, 1.079117e+00,
1548 1.071856e+00, 1.071990e+00, 1.073774e+00, 1.079356e+00, 1.081314e+00, 1.082056e+00,
1549 1.090772e+00, 1.096776e+00, 1.095828e+00, 1.097678e+00, 1.099157e+00, 1.103677e+00,
1550 1.105132e+00, 1.109806e+00, 1.110816e+00, 1.117378e+00, 1.115165e+00, 1.115710e+00,
1551 1.111855e+00, 1.110482e+00, 1.110112e+00, 1.106676e+00, 1.108706e+00, 1.105549e+00,
1552 1.106318e+00, 1.106242e+00, 1.107672e+00, 1.107342e+00, 1.108119e+00, 1.106655e+00,
1553 1.102588e+00, 1.096657e+00, 1.092920e+00, 1.086629e+00, 1.083592e+00, 1.076030e+00,
1554 1.083777e+00, 1.089460e+00, 1.086545e+00, 1.079924e+00, 1.082218e+00, 1.077798e+00,
1555 1.077062e+00, 1.072825e+00, 1.072241e+00, 1.072104e+00, 1.072490e+00, 1.069829e+00,
1556 1.070398e+00, 1.065458e+00, 1.064968e+00, 1.060524e+00, 1.060048e+00, 1.057620e+00,
1557 1.056428e+00, 1.055366e+00, 1.055017e+00, 1.052304e+00, 1.051767e+00, 1.049728e+00,
1558 1.048745e+00, 1.047399e+00, 1.045876e+00, 1.042972e+00, 1.041824e+00, 1.039993e+00,
1559 1.039021e+00, 1.036627e+00, 1.034176e+00, 1.032526e+00, 1.033633e+00, 1.036107e+00,
1560 1.037803e+00, 1.031266e+00, 1.032991e+00, 1.033284e+00, 1.035015e+00, 1.033945e+00,
1561 1.037075e+00, 1.034721e+00
1562 
1563 };
1564 
1565 const G4double G4ComponentGGHadronNucleusXsc::fNeutronBarCorrectionIn[93] = {
1566 
1567 1.0, 1.0, 1.167421e+00, 1.156250e+00, 1.205364e+00, 1.154225e+00, 1.120391e+00, // 6
1568 1.124632e+00, 1.129460e+00, 1.107863e+00, 1.102152e+00, 1.104593e+00, 1.100285e+00, // 12
1569 1.098450e+00, 1.092677e+00, 1.101124e+00, 1.106461e+00, 1.115049e+00, 1.123903e+00, // 18
1570 1.126661e+00, 1.131259e+00, 1.133949e+00, 1.134185e+00, 1.133767e+00, 1.132813e+00, // 24
1571 1.131515e+00, 1.144338e+00, // 1.130338e+00,
1572 1.134171e+00, 1.139206e+00, 1.148474e+00, // 1.141474e+00,
1573 1.142189e+00,
1574 1.140725e+00, 1.140100e+00, 1.139848e+00, 1.137674e+00, 1.138645e+00, 1.136339e+00,
1575 1.136439e+00, 1.135946e+00, 1.136431e+00, 1.135702e+00, 1.135703e+00, 1.134113e+00,
1576 1.131935e+00, 1.128381e+00, 1.126373e+00, 1.122453e+00, 1.120908e+00, 1.115953e+00,
1577 1.115947e+00, 1.114426e+00, 1.111749e+00, 1.106207e+00, 1.107494e+00, 1.103622e+00,
1578 1.102576e+00, 1.098816e+00, 1.097889e+00, 1.097306e+00, 1.097130e+00, 1.094578e+00,
1579 1.094552e+00, 1.090222e+00, 1.089358e+00, 1.085409e+00, 1.084560e+00, 1.082182e+00,
1580 1.080773e+00, 1.079464e+00, 1.078724e+00, 1.076121e+00, 1.075235e+00, 1.073159e+00,
1581 1.071920e+00, 1.070395e+00, 1.069503e+00, 1.067525e+00, 1.066919e+00, 1.065779e+00,
1582 1.065319e+00, 1.063730e+00, 1.062092e+00, 1.061085e+00, 1.059908e+00, 1.059815e+00,
1583 1.059109e+00, 1.051920e+00, 1.051258e+00, 1.049473e+00, 1.048823e+00, 1.045984e+00,
1584 1.046435e+00, 1.042614e+00
1585 
1586 };
1587 
1588 const G4double G4ComponentGGHadronNucleusXsc::fProtonBarCorrectionTot[93] = {
1589 
1590 1.0, 1.0,
1591 1.118515e+00, 1.082000e+00, 1.116169e+00, 1.078745e+00, 1.061313e+00, 1.058203e+00,
1592 1.082661e+00, 1.068498e+00, 1.076910e+00, 1.083474e+00, 1.079115e+00, 1.071854e+00,
1593 1.071988e+00, 1.073772e+00, 1.079355e+00, 1.081312e+00, 1.082054e+00, 1.090770e+00,
1594 1.096774e+00, 1.095827e+00, 1.097677e+00, 1.099156e+00, 1.103676e+00, 1.105130e+00,
1595 1.109805e+00, 1.110814e+00, 1.117377e+00, 1.115163e+00, 1.115708e+00, 1.111853e+00,
1596 1.110480e+00, 1.110111e+00, 1.106674e+00, 1.108705e+00, 1.105548e+00, 1.106317e+00,
1597 1.106241e+00, 1.107671e+00, 1.107341e+00, 1.108118e+00, 1.106654e+00, 1.102586e+00,
1598 1.096655e+00, 1.092918e+00, 1.086628e+00, 1.083590e+00, 1.076028e+00, 1.083776e+00,
1599 1.089458e+00, 1.086543e+00, 1.079923e+00, 1.082216e+00, 1.077797e+00, 1.077061e+00,
1600 1.072824e+00, 1.072239e+00, 1.072103e+00, 1.072488e+00, 1.069828e+00, 1.070396e+00,
1601 1.065456e+00, 1.064966e+00, 1.060523e+00, 1.060047e+00, 1.057618e+00, 1.056427e+00,
1602 1.055365e+00, 1.055016e+00, 1.052303e+00, 1.051766e+00, 1.049727e+00, 1.048743e+00,
1603 1.047397e+00, 1.045875e+00, 1.042971e+00, 1.041823e+00, 1.039992e+00, 1.039019e+00,
1604 1.036626e+00, 1.034175e+00, 1.032525e+00, 1.033632e+00, 1.036106e+00, 1.037802e+00,
1605 1.031265e+00, 1.032990e+00, 1.033283e+00, 1.035014e+00, 1.033944e+00, 1.037074e+00,
1606 1.034720e+00
1607 
1608 };
1609 
1610 const G4double G4ComponentGGHadronNucleusXsc::fProtonBarCorrectionIn[93] = {
1611 
1612 1.0, 1.0,
1613 1.147419e+00, // 1.167419e+00,
1614 1.156248e+00, 1.205362e+00, 1.154224e+00, 1.120390e+00, 1.124630e+00, // 7
1615 1.129459e+00, 1.107861e+00, 1.102151e+00, 1.104591e+00, 1.100284e+00, 1.098449e+00, // 13
1616 1.092675e+00, 1.101122e+00, 1.106460e+00, 1.115048e+00, 1.123902e+00, 1.126659e+00, // 19
1617 1.131258e+00, 1.133948e+00, 1.134183e+00, 1.133766e+00, 1.132812e+00, 1.131514e+00, // 25
1618 1.140337e+00, // 1.130337e+00,
1619 
1620 1.134170e+00, 1.139205e+00, 1.151472e+00, // 1.141472e+00,
1621 1.142188e+00, 1.140724e+00,
1622 1.140099e+00, 1.139847e+00, 1.137672e+00, 1.138644e+00, 1.136338e+00, 1.136438e+00,
1623 1.135945e+00, 1.136429e+00, 1.135701e+00, 1.135702e+00, 1.134112e+00, 1.131934e+00,
1624 1.128380e+00, 1.126371e+00, 1.122452e+00, 1.120907e+00, 1.115952e+00, 1.115946e+00,
1625 1.114425e+00, 1.111748e+00, 1.106205e+00, 1.107493e+00, 1.103621e+00, 1.102575e+00,
1626 1.098815e+00, 1.097888e+00, 1.097305e+00, 1.097129e+00, 1.094577e+00, 1.094551e+00,
1627 1.090221e+00, 1.089357e+00, 1.085408e+00, 1.084559e+00, 1.082181e+00, 1.080772e+00,
1628 1.079463e+00, 1.078723e+00, 1.076120e+00, 1.075234e+00, 1.073158e+00, 1.071919e+00,
1629 1.070394e+00, 1.069502e+00, 1.067524e+00, 1.066918e+00, 1.065778e+00, 1.065318e+00,
1630 1.063729e+00, 1.062091e+00, 1.061084e+00, 1.059907e+00, 1.059814e+00, 1.059108e+00,
1631 1.051919e+00, 1.051257e+00, 1.049472e+00, 1.048822e+00, 1.045983e+00, 1.046434e+00,
1632 1.042613e+00
1633 
1634 };
1635 
1636 
1637 const G4double G4ComponentGGHadronNucleusXsc::fPionPlusBarCorrectionTot[93] = {
1638 
1639 1.0, 1.0,
1640 1.075927e+00, 1.074407e+00, 1.126098e+00, 1.100127e+00, 1.089742e+00, 1.083536e+00,
1641 1.089988e+00, 1.103566e+00, 1.096922e+00, 1.126573e+00, 1.132734e+00, 1.136512e+00,
1642 1.136629e+00, 1.133086e+00, 1.132428e+00, 1.129299e+00, 1.125622e+00, 1.126992e+00,
1643 1.127840e+00, 1.162670e+00, 1.160392e+00, 1.157864e+00, 1.157227e+00, 1.154627e+00,
1644 1.192555e+00, 1.197243e+00, 1.197911e+00, 1.200326e+00, 1.220053e+00, 1.215019e+00,
1645 1.211703e+00, 1.209080e+00, 1.204248e+00, 1.203328e+00, 1.198671e+00, 1.196840e+00,
1646 1.194392e+00, 1.193037e+00, 1.190408e+00, 1.188583e+00, 1.206127e+00, 1.210028e+00,
1647 1.206434e+00, 1.204456e+00, 1.200547e+00, 1.199058e+00, 1.200174e+00, 1.200276e+00,
1648 1.198912e+00, 1.213048e+00, 1.207160e+00, 1.208020e+00, 1.203814e+00, 1.202380e+00,
1649 1.198306e+00, 1.197002e+00, 1.196027e+00, 1.195449e+00, 1.192563e+00, 1.192135e+00,
1650 1.187556e+00, 1.186308e+00, 1.182124e+00, 1.180900e+00, 1.178224e+00, 1.176471e+00,
1651 1.174811e+00, 1.173702e+00, 1.170827e+00, 1.169581e+00, 1.167205e+00, 1.165626e+00,
1652 1.180244e+00, 1.177626e+00, 1.175121e+00, 1.173903e+00, 1.172192e+00, 1.171128e+00,
1653 1.168997e+00, 1.166826e+00, 1.164130e+00, 1.165412e+00, 1.165504e+00, 1.165020e+00,
1654 1.158462e+00, 1.158014e+00, 1.156519e+00, 1.156081e+00, 1.153602e+00, 1.154190e+00,
1655 1.152974e+00
1656 
1657 };
1658 
1659 const G4double G4ComponentGGHadronNucleusXsc::fPionPlusBarCorrectionIn[93] = {
1660 
1661 1.0, 1.0,
1662 1.140246e+00, 1.097872e+00, 1.104301e+00, 1.068722e+00, 1.056495e+00, 1.062622e+00, // 7
1663 1.047987e+00, 1.037032e+00, 1.035686e+00, 1.042870e+00, 1.052222e+00, 1.075100e+00, // 13
1664 1.084480e+00, 1.078286e+00, 1.081488e+00, 1.089713e+00, 1.099105e+00, 1.098003e+00, // 19
1665 1.102175e+00, 1.117707e+00, 1.121734e+00, 1.125229e+00, 1.126457e+00, 1.128905e+00, // 25
1666 1.163312e+00, 1.126263e+00, 1.126459e+00, 1.135191e+00, 1.116986e+00, 1.117184e+00, // 31
1667 1.117037e+00, 1.116777e+00, 1.115858e+00, 1.115745e+00, 1.114489e+00, 1.113993e+00, // 37
1668 1.113226e+00, 1.112818e+00, 1.111890e+00, 1.111238e+00, 1.111209e+00, 1.111775e+00, // 43
1669 1.110256e+00, 1.109414e+00, 1.107647e+00, 1.106980e+00, 1.106096e+00, 1.107331e+00, // 49
1670 1.107849e+00, 1.106407e+00, 1.103426e+00, 1.103896e+00, 1.101756e+00, 1.101031e+00, // 55
1671 1.098915e+00, 1.098260e+00, 1.097768e+00, 1.097487e+00, 1.095964e+00, 1.095773e+00, // 61
1672 1.093348e+00, 1.092687e+00, 1.090465e+00, 1.089821e+00, 1.088394e+00, 1.087462e+00, // 67
1673 1.086571e+00, 1.085997e+00, 1.084451e+00, 1.083798e+00, 1.082513e+00, 1.081670e+00, // 73
1674 1.080735e+00, 1.075659e+00, 1.074341e+00, 1.073689e+00, 1.072787e+00, 1.072237e+00, // 79
1675 1.071107e+00, 1.069955e+00, 1.074856e+00, 1.065873e+00, 1.065938e+00, 1.065694e+00,
1676 1.062192e+00, 1.061967e+00, 1.061180e+00, 1.060960e+00, 1.059646e+00, 1.059975e+00,
1677 1.059658e+00
1678 
1679 };
1680 
1681 
1682 const G4double G4ComponentGGHadronNucleusXsc::fPionMinusBarCorrectionTot[93] = {
1683 
1684 1.0, 1.0,
1685 1.3956e+00, 1.077959e+00, 1.129145e+00, 1.102088e+00, 1.089765e+00, 1.083542e+00, // 7
1686 1.089995e+00, 1.104895e+00, 1.097154e+00, 1.127663e+00, 1.133063e+00, 1.137425e+00, // 13
1687 1.136724e+00, 1.133859e+00, 1.132498e+00, 1.130276e+00, 1.127896e+00, 1.127656e+00, // 19
1688 1.127905e+00, 1.164210e+00, 1.162259e+00, 1.160075e+00, 1.158978e+00, 1.156649e+00, // 25
1689 1.194157e+00, 1.199177e+00, 1.198983e+00, 1.202325e+00, 1.221967e+00, 1.217548e+00,
1690 1.214389e+00, 1.211760e+00, 1.207335e+00, 1.206081e+00, 1.201766e+00, 1.199779e+00,
1691 1.197283e+00, 1.195706e+00, 1.193071e+00, 1.191115e+00, 1.208838e+00, 1.212681e+00,
1692 1.209235e+00, 1.207163e+00, 1.203451e+00, 1.201807e+00, 1.203283e+00, 1.203388e+00,
1693 1.202244e+00, 1.216509e+00, 1.211066e+00, 1.211504e+00, 1.207539e+00, 1.205991e+00,
1694 1.202143e+00, 1.200724e+00, 1.199595e+00, 1.198815e+00, 1.196025e+00, 1.195390e+00,
1695 1.191137e+00, 1.189791e+00, 1.185888e+00, 1.184575e+00, 1.181996e+00, 1.180229e+00,
1696 1.178545e+00, 1.177355e+00, 1.174616e+00, 1.173312e+00, 1.171016e+00, 1.169424e+00,
1697 1.184120e+00, 1.181478e+00, 1.179085e+00, 1.177817e+00, 1.176124e+00, 1.175003e+00,
1698 1.172947e+00, 1.170858e+00, 1.168170e+00, 1.169397e+00, 1.169304e+00, 1.168706e+00,
1699 1.162774e+00, 1.162217e+00, 1.160740e+00, 1.160196e+00, 1.157857e+00, 1.158220e+00,
1700 1.157267e+00
1701 };
1702 
1703 
1704 const G4double G4ComponentGGHadronNucleusXsc::fPionMinusBarCorrectionIn[93] = {
1705 
1706 1.0, 1.0,
1707 1.463e+00, 1.100898e+00, 1.106773e+00, 1.070289e+00, 1.040514e+00, 1.062628e+00, // 7
1708 1.047992e+00, 1.038041e+00, 1.035862e+00, 1.043679e+00, 1.052466e+00, 1.065780e+00, // 13
1709 1.070551e+00, 1.078869e+00, 1.081541e+00, 1.090455e+00, 1.100847e+00, 1.098511e+00, // 19
1710 1.102226e+00, 1.118865e+00, 1.123143e+00, 1.126904e+00, 1.127785e+00, 1.130444e+00, // 25
1711 1.148502e+00, 1.127678e+00, 1.127244e+00, 1.123634e+00, 1.118347e+00, 1.118988e+00,
1712 1.118957e+00, 1.118696e+00, 1.118074e+00, 1.117722e+00, 1.116717e+00, 1.116111e+00,
1713 1.115311e+00, 1.114745e+00, 1.113814e+00, 1.113069e+00, 1.113141e+00, 1.113660e+00,
1714 1.112249e+00, 1.111343e+00, 1.109718e+00, 1.108942e+00, 1.108310e+00, 1.109549e+00,
1715 1.110227e+00, 1.108846e+00, 1.106183e+00, 1.106354e+00, 1.104388e+00, 1.103583e+00,
1716 1.101632e+00, 1.100896e+00, 1.100296e+00, 1.099873e+00, 1.098420e+00, 1.098082e+00,
1717 1.095892e+00, 1.095162e+00, 1.093144e+00, 1.092438e+00, 1.091083e+00, 1.090142e+00,
1718 1.089236e+00, 1.088604e+00, 1.087159e+00, 1.086465e+00, 1.085239e+00, 1.084388e+00,
1719 1.083473e+00, 1.078373e+00, 1.077136e+00, 1.076450e+00, 1.075561e+00, 1.074973e+00,
1720 1.073898e+00, 1.072806e+00, 1.067706e+00, 1.068684e+00, 1.068618e+00, 1.068294e+00,
1721 1.065241e+00, 1.064939e+00, 1.064166e+00, 1.063872e+00, 1.062659e+00, 1.062828e+00,
1722 1.062699e+00
1723 
1724 };
1725 
1726 
1727 //
1728 //
G4double GetElasticHadronNucleonXsc()
G4double GetRatioQE(const G4DynamicParticle *, G4int At, G4int Zt)
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
const int N
Definition: mixmax.h:43
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
virtual G4double GetProductionElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
static G4AntiOmegaMinus * AntiOmegaMinus()
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetN() const
Definition: G4Element.hh:135
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
static G4OmegaMinus * OmegaMinus()
G4double GetZ() const
Definition: G4Element.hh:131
static G4KaonZeroLong * KaonZeroLong()
G4double GetRatioSD(const G4DynamicParticle *, G4int At, G4int Zt)
G4ParticleDefinition * GetDefinition() const
double B(double temperature)
static G4AntiSigmaPlus * AntiSigmaPlus()
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetKaonNucleonXscGG(const G4DynamicParticle *, const G4ParticleDefinition *)
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4AntiXiMinus * AntiXiMinus()
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual G4double GetInelasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
static G4SigmaMinus * SigmaMinus()
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual void CrossSectionDescription(std::ostream &) const
static G4ParticleTable * GetParticleTable()
static G4AntiLambda * AntiLambda()
int G4lrint(double ad)
Definition: templates.hh:163
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
virtual G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
static G4AntiXiZero * AntiXiZero()
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
G4bool IsIsoApplicable(const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
static constexpr double pi
Definition: G4SIunits.hh:75
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
double G4double
Definition: G4Types.hh:76
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
static constexpr double fermi
Definition: G4SIunits.hh:103
G4ThreeVector G4ParticleMomentum
static G4He3 * He3()
Definition: G4He3.cc:94
double mag() const
const G4double oneThird
virtual G4double ComputeQuasiElasticRatio(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
static G4AntiNeutron * AntiNeutron()
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4double GetInelasticHadronNucleonXsc()
G4ThreeVector GetMomentum() const
virtual G4double GetProductionIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)