Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ComponentGGHadronNucleusXsc Class Reference

#include <G4ComponentGGHadronNucleusXsc.hh>

Inheritance diagram for G4ComponentGGHadronNucleusXsc:
Collaboration diagram for G4ComponentGGHadronNucleusXsc:

Public Member Functions

 G4ComponentGGHadronNucleusXsc ()
 
virtual ~G4ComponentGGHadronNucleusXsc ()
 
virtual G4double GetTotalIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetTotalElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetInelasticIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetProductionIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double GetInelasticElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetProductionElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetElasticElementCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
 
virtual G4double GetElasticIsotopeCrossSection (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
virtual G4double ComputeQuasiElasticRatio (const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
 
G4bool IsIsoApplicable (const G4DynamicParticle *aDP, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetRatioSD (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetRatioQE (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXscPDG (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXscPDG (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXscNS (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXscNS (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHNinelasticXsc (const G4DynamicParticle *, const G4Element *)
 
G4double GetHNinelasticXsc (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHNinelasticXscVU (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double CalculateEcmValue (const G4double, const G4double, const G4double)
 
G4double CalcMandelstamS (const G4double, const G4double, const G4double)
 
G4double GetNucleusRadius (const G4DynamicParticle *, const G4Element *)
 
G4double GetNucleusRadius (G4int At)
 
G4double GetAxsc2piR2 ()
 
G4double GetModelInLog ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
G4double GetElasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetInelasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetTotalGlauberGribovXsc ()
 
G4double GetElasticGlauberGribovXsc ()
 
G4double GetInelasticGlauberGribovXsc ()
 
G4double GetProductionGlauberGribovXsc ()
 
G4double GetDiffractionGlauberGribovXsc ()
 
G4double GetRadiusConst ()
 
G4double GetParticleBarCorTot (const G4ParticleDefinition *theParticle, G4int Z)
 
G4double GetParticleBarCorIn (const G4ParticleDefinition *theParticle, G4int Z)
 
void SetEnergyLowerLimit (G4double E)
 
- Public Member Functions inherited from G4VComponentCrossSection
 G4VComponentCrossSection (const G4String &nam="")
 
virtual ~G4VComponentCrossSection ()
 
G4double GetTotalElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
G4double GetInelasticElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
G4double GetElasticElementCrossSection (const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void Description () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Static Public Member Functions

static const char * Default_Name ()
 

Detailed Description

Definition at line 51 of file G4ComponentGGHadronNucleusXsc.hh.

Constructor & Destructor Documentation

G4ComponentGGHadronNucleusXsc::G4ComponentGGHadronNucleusXsc ( )

Definition at line 46 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4OmegaMinus * OmegaMinus()
static G4KaonZeroLong * KaonZeroLong()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
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
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4SigmaMinus * SigmaMinus()
static G4AntiLambda * AntiLambda()
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
G4VComponentCrossSection(const G4String &nam="")
static G4AntiXiZero * AntiXiZero()
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static constexpr double fermi
Definition: G4SIunits.hh:103
static G4He3 * He3()
Definition: G4He3.cc:94
static G4AntiNeutron * AntiNeutron()

Here is the call graph for this function:

G4ComponentGGHadronNucleusXsc::~G4ComponentGGHadronNucleusXsc ( )
virtual

Definition at line 93 of file G4ComponentGGHadronNucleusXsc.cc.

94 {
95  if (hnXsc) delete hnXsc;
96 }

Member Function Documentation

G4double G4ComponentGGHadronNucleusXsc::CalcMandelstamS ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 1513 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::CalculateEcmValue ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 1497 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
double G4double
Definition: G4Types.hh:76
G4double G4ComponentGGHadronNucleusXsc::ComputeQuasiElasticRatio ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Reimplemented from G4VComponentCrossSection.

Definition at line 212 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
double G4double
Definition: G4Types.hh:76
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

void G4ComponentGGHadronNucleusXsc::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Definition at line 1527 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static const char* G4ComponentGGHadronNucleusXsc::Default_Name ( )
inlinestatic

Definition at line 58 of file G4ComponentGGHadronNucleusXsc.hh.

58 { return "Glauber-Gribov"; }
G4double G4ComponentGGHadronNucleusXsc::GetAxsc2piR2 ( )
inline

Definition at line 143 of file G4ComponentGGHadronNucleusXsc.hh.

143 {return fAxsc2piR2;};
G4double G4ComponentGGHadronNucleusXsc::GetDiffractionGlauberGribovXsc ( )
inline

Definition at line 155 of file G4ComponentGGHadronNucleusXsc.hh.

155 { return fDiffractionXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 184 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 226 of file G4ComponentGGHadronNucleusXsc.hh.

228 {
229  GetIsoCrossSection(dp, Z, A);
230  return fElasticXsc;
231 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribovXsc ( )
inline

Definition at line 152 of file G4ComponentGGHadronNucleusXsc.hh.

152 { return fElasticXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetElasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 198 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 525 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 542 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
G4double GetMass() const
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 746 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 764 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double GetTotalEnergy() const
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4double GetMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4ParticleTable * GetParticleTable()
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 612 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4Element *)
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 631 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4ParticleDefinition * GetDefinition() const
double B(double temperature)
int G4int
Definition: G4Types.hh:78
G4double GetMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4ParticleTable * GetParticleTable()
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 1143 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
int G4lrint(double ad)
Definition: templates.hh:163

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1157 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetHNinelasticXscVU ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 1185 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double GetTotalEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
double mag() const
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 156 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 237 of file G4ComponentGGHadronNucleusXsc.hh.

239 {
240  GetIsoCrossSection(dp, Z, A);
241  return fInelasticXsc;
242 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribovXsc ( )
inline

Definition at line 153 of file G4ComponentGGHadronNucleusXsc.hh.

153 { return fInelasticXsc; };

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetInelasticIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 128 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)

Definition at line 294 of file G4ComponentGGHadronNucleusXsc.cc.

299 {
300  G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
301  G4double hpInXsc(0.), hnInXsc(0.);
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 }
G4double GetElasticHadronNucleonXsc()
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
const int N
Definition: mixmax.h:43
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetKaonNucleonXscGG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetParticleBarCorIn(const G4ParticleDefinition *theParticle, G4int Z)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
static constexpr double fermi
Definition: G4SIunits.hh:103
G4double GetParticleBarCorTot(const G4ParticleDefinition *theParticle, G4int Z)
G4double GetInelasticHadronNucleonXsc()

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetModelInLog ( )
inline

Definition at line 144 of file G4ComponentGGHadronNucleusXsc.hh.

144 {return fModelInLog;};
G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( const G4DynamicParticle ,
const G4Element anElement 
)

Definition at line 1397 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double GetN() const
Definition: G4Element.hh:135
int G4int
Definition: G4Types.hh:78
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76
const G4double oneThird

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetNucleusRadius ( G4int  At)

Definition at line 1454 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
const G4double oneThird

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorIn ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 270 of file G4ComponentGGHadronNucleusXsc.hh.

272 {
273  if(Z >= 2 && Z <= 92)
274  {
275  if( theParticle == theProton ) return fProtonBarCorrectionIn[Z];
276  else if( theParticle == theNeutron) return fNeutronBarCorrectionIn[Z];
277  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionIn[Z];
278  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionIn[Z];
279  else return 1.0;
280  }
281  else return 1.0;
282 }

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetParticleBarCorTot ( const G4ParticleDefinition theParticle,
G4int  Z 
)
inline

Definition at line 250 of file G4ComponentGGHadronNucleusXsc.hh.

252 {
253  if(Z >= 2 && Z <= 92)
254  {
255  if( theParticle == theProton ) return fProtonBarCorrectionTot[Z];
256  else if( theParticle == theNeutron) return fNeutronBarCorrectionTot[Z];
257  else if( theParticle == thePiPlus ) return fPionPlusBarCorrectionTot[Z];
258  else if( theParticle == thePiMinus) return fPionMinusBarCorrectionTot[Z];
259  else return 1.0;
260  }
261  else return 1.0;
262 }

Here is the caller graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetProductionElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Definition at line 170 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetProductionGlauberGribovXsc ( )
inline

Definition at line 154 of file G4ComponentGGHadronNucleusXsc.hh.

154 { return fProductionXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetProductionIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Definition at line 142 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetRadiusConst ( )
inline

Definition at line 156 of file G4ComponentGGHadronNucleusXsc.hh.

156 { return fRadiusConst; };
G4double G4ComponentGGHadronNucleusXsc::GetRatioQE ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 478 of file G4ComponentGGHadronNucleusXsc.cc.

479 {
480  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
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 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double GetHNinelasticXsc(const G4DynamicParticle *, const G4Element *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetRatioSD ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 436 of file G4ComponentGGHadronNucleusXsc.cc.

437 {
438  G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
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 }
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4ParticleDefinition * GetDefinition() const
double A(double temperature)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4Element *)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetTotalElementCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4double  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 114 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
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)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4double G4ComponentGGHadronNucleusXsc::GetTotalGlauberGribovXsc ( )
inline

Definition at line 151 of file G4ComponentGGHadronNucleusXsc.hh.

151 { return fTotalXsc; };
G4double G4ComponentGGHadronNucleusXsc::GetTotalIsotopeCrossSection ( const G4ParticleDefinition aParticle,
G4double  kinEnergy,
G4int  Z,
G4int  A 
)
virtual

Implements G4VComponentCrossSection.

Definition at line 100 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
double A(double temperature)
G4ThreeVector G4ParticleMomentum

Here is the call graph for this function:

G4bool G4ComponentGGHadronNucleusXsc::IsIsoApplicable ( const G4DynamicParticle aDP,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)

Definition at line 236 of file G4ComponentGGHadronNucleusXsc.cc.

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 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4ComponentGGHadronNucleusXsc::SetEnergyLowerLimit ( G4double  E)
inline

Definition at line 161 of file G4ComponentGGHadronNucleusXsc.hh.

161 {fLowerLimit=E;};

The documentation for this class was generated from the following files: