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